| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | /* | ||
| 2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
| 4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
| 5 | |||
| 6 | This file is part of SLEPc. | ||
| 7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
| 8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 9 | */ | ||
| 10 | /* | ||
| 11 | SVD routines for setting up the solver | ||
| 12 | */ | ||
| 13 | |||
| 14 | #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/ | ||
| 15 | |||
| 16 | /*@ | ||
| 17 | SVDSetOperators - Set the matrices associated with the singular value problem. | ||
| 18 | |||
| 19 | Collective | ||
| 20 | |||
| 21 | Input Parameters: | ||
| 22 | + svd - the singular value solver context | ||
| 23 | . A - the matrix associated with the singular value problem | ||
| 24 | - B - the second matrix in the case of GSVD | ||
| 25 | |||
| 26 | Level: beginner | ||
| 27 | |||
| 28 | .seealso: [](ch:svd), `SVDSolve()`, `SVDGetOperators()` | ||
| 29 | @*/ | ||
| 30 | 2617 | PetscErrorCode SVDSetOperators(SVD svd,Mat A,Mat B) | |
| 31 | { | ||
| 32 | 2617 | PetscInt Ma,Na,Mb,Nb,ma,na,mb,nb,M0,N0,m0,n0; | |
| 33 | 2617 | PetscBool samesize=PETSC_TRUE; | |
| 34 | |||
| 35 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2617 | PetscFunctionBegin; |
| 36 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2617 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
| 37 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2617 | PetscValidHeaderSpecific(A,MAT_CLASSID,2); |
| 38 |
4/14✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
2617 | if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); |
| 39 |
13/32✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
2617 | PetscCheckSameComm(svd,1,A,2); |
| 40 |
15/34✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✗ Branch 32 not taken.
✗ Branch 33 not taken.
|
2617 | if (B) PetscCheckSameComm(svd,1,B,3); |
| 41 | |||
| 42 | /* Check matrix sizes */ | ||
| 43 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2617 | PetscCall(MatGetSize(A,&Ma,&Na)); |
| 44 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2617 | PetscCall(MatGetLocalSize(A,&ma,&na)); |
| 45 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2617 | if (svd->OP) { |
| 46 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
274 | PetscCall(MatGetSize(svd->OP,&M0,&N0)); |
| 47 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
274 | PetscCall(MatGetLocalSize(svd->OP,&m0,&n0)); |
| 48 |
4/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
|
274 | if (M0!=Ma || N0!=Na || m0!=ma || n0!=na) samesize = PETSC_FALSE; |
| 49 | } | ||
| 50 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2617 | if (B) { |
| 51 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
985 | PetscCall(MatGetSize(B,&Mb,&Nb)); |
| 52 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
985 | PetscCall(MatGetLocalSize(B,&mb,&nb)); |
| 53 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
985 | PetscCheck(Na==Nb,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Different number of columns in A (%" PetscInt_FMT ") and B (%" PetscInt_FMT ")",Na,Nb); |
| 54 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
985 | PetscCheck(na==nb,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Different local column size in A (%" PetscInt_FMT ") and B (%" PetscInt_FMT ")",na,nb); |
| 55 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
985 | if (svd->OPb) { |
| 56 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
134 | PetscCall(MatGetSize(svd->OPb,&M0,&N0)); |
| 57 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
134 | PetscCall(MatGetLocalSize(svd->OPb,&m0,&n0)); |
| 58 |
4/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
|
134 | if (M0!=Mb || N0!=Nb || m0!=mb || n0!=nb) samesize = PETSC_FALSE; |
| 59 | } | ||
| 60 | } | ||
| 61 | |||
| 62 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2617 | PetscCall(PetscObjectReference((PetscObject)A)); |
| 63 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
2617 | if (B) PetscCall(PetscObjectReference((PetscObject)B)); |
| 64 |
3/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
2617 | if (svd->state && !samesize) PetscCall(SVDReset(svd)); |
| 65 | else { | ||
| 66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2617 | PetscCall(MatDestroy(&svd->OP)); |
| 67 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2617 | PetscCall(MatDestroy(&svd->OPb)); |
| 68 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2617 | PetscCall(MatDestroy(&svd->A)); |
| 69 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2617 | PetscCall(MatDestroy(&svd->B)); |
| 70 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2617 | PetscCall(MatDestroy(&svd->AT)); |
| 71 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2617 | PetscCall(MatDestroy(&svd->BT)); |
| 72 | } | ||
| 73 | 2617 | svd->nrma = 0.0; | |
| 74 | 2617 | svd->nrmb = 0.0; | |
| 75 | 2617 | svd->OP = A; | |
| 76 | 2617 | svd->OPb = B; | |
| 77 | 2617 | svd->state = SVD_STATE_INITIAL; | |
| 78 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2617 | PetscFunctionReturn(PETSC_SUCCESS); |
| 79 | } | ||
| 80 | |||
| 81 | /*@ | ||
| 82 | SVDGetOperators - Get the matrices associated with the singular value problem. | ||
| 83 | |||
| 84 | Collective | ||
| 85 | |||
| 86 | Input Parameter: | ||
| 87 | . svd - the singular value solver context | ||
| 88 | |||
| 89 | Output Parameters: | ||
| 90 | + A - the matrix associated with the singular value problem | ||
| 91 | - B - the second matrix in the case of GSVD | ||
| 92 | |||
| 93 | Level: intermediate | ||
| 94 | |||
| 95 | .seealso: [](ch:svd), `SVDSolve()`, `SVDSetOperators()` | ||
| 96 | @*/ | ||
| 97 | 161 | PetscErrorCode SVDGetOperators(SVD svd,Mat *A,Mat *B) | |
| 98 | { | ||
| 99 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
161 | PetscFunctionBegin; |
| 100 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
161 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
| 101 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
161 | if (A) *A = svd->OP; |
| 102 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 10 times.
|
161 | if (B) *B = svd->OPb; |
| 103 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
161 | PetscFunctionReturn(PETSC_SUCCESS); |
| 104 | } | ||
| 105 | |||
| 106 | /*@ | ||
| 107 | SVDSetSignature - Set the signature matrix defining a hyperbolic singular value problem. | ||
| 108 | |||
| 109 | Collective | ||
| 110 | |||
| 111 | Input Parameters: | ||
| 112 | + svd - the singular value solver context | ||
| 113 | - omega - a vector containing the diagonal elements of the signature matrix (or `NULL`) | ||
| 114 | |||
| 115 | Notes: | ||
| 116 | The signature matrix is relevant only for hyperbolic problems (HSVD). | ||
| 117 | Use `NULL` to reset a previously set signature. | ||
| 118 | |||
| 119 | Level: intermediate | ||
| 120 | |||
| 121 | .seealso: [](ch:svd), `SVDSetProblemType()`, `SVDSetOperators()`, `SVDGetSignature()`, `SVD_HYPERBOLIC` | ||
| 122 | @*/ | ||
| 123 | 262 | PetscErrorCode SVDSetSignature(SVD svd,Vec omega) | |
| 124 | { | ||
| 125 | 262 | PetscInt N,Ma,n,ma; | |
| 126 | |||
| 127 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
262 | PetscFunctionBegin; |
| 128 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
262 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
| 129 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
262 | if (omega) { |
| 130 |
2/12✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
262 | PetscValidHeaderSpecific(omega,VEC_CLASSID,2); |
| 131 |
13/32✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
262 | PetscCheckSameComm(svd,1,omega,2); |
| 132 | } | ||
| 133 | |||
| 134 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 9 times.
✓ Branch 3 taken 9 times.
|
262 | if (omega && svd->OP) { /* Check sizes */ |
| 135 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
252 | PetscCall(VecGetSize(omega,&N)); |
| 136 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
252 | PetscCall(VecGetLocalSize(omega,&n)); |
| 137 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
252 | PetscCall(MatGetSize(svd->OP,&Ma,NULL)); |
| 138 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
252 | PetscCall(MatGetLocalSize(svd->OP,&ma,NULL)); |
| 139 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
252 | PetscCheck(N==Ma,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_SIZ,"Global size of signature (%" PetscInt_FMT ") does not match the row size of A (%" PetscInt_FMT ")",N,Ma); |
| 140 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
252 | PetscCheck(n==ma,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_SIZ,"Local size of signature (%" PetscInt_FMT ") does not match the local row size of A (%" PetscInt_FMT ")",n,ma); |
| 141 | } | ||
| 142 | |||
| 143 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
262 | PetscCall(VecDestroy(&svd->omega)); |
| 144 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
262 | if (omega) { |
| 145 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
262 | PetscCall(VecDuplicate(omega,&svd->omega)); |
| 146 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
262 | PetscCall(VecCopy(omega,svd->omega)); |
| 147 | } | ||
| 148 | 262 | svd->state = SVD_STATE_INITIAL; | |
| 149 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
262 | PetscFunctionReturn(PETSC_SUCCESS); |
| 150 | } | ||
| 151 | |||
| 152 | /*@ | ||
| 153 | SVDGetSignature - Get the signature matrix defining a hyperbolic singular value problem. | ||
| 154 | |||
| 155 | Not Collective | ||
| 156 | |||
| 157 | Input Parameter: | ||
| 158 | . svd - the singular value solver context | ||
| 159 | |||
| 160 | Output Parameter: | ||
| 161 | . omega - a vector containing the diagonal elements of the signature matrix | ||
| 162 | |||
| 163 | Notes: | ||
| 164 | The signature matrix is relevant only for hyperbolic problems (HSVD). | ||
| 165 | If no signature has been set, this function will return a vector of all ones. | ||
| 166 | |||
| 167 | The user should pass a previously created `Vec` with the appropriate dimension. | ||
| 168 | |||
| 169 | Level: intermediate | ||
| 170 | |||
| 171 | .seealso: [](ch:svd), `SVDSetSignature()`, `SVD_HYPERBOLIC` | ||
| 172 | @*/ | ||
| 173 | 10 | PetscErrorCode SVDGetSignature(SVD svd,Vec omega) | |
| 174 | { | ||
| 175 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
| 176 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
| 177 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(omega,VEC_CLASSID,2); |
| 178 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
10 | if (svd->omega) PetscCall(VecCopy(svd->omega,omega)); |
| 179 | ✗ | else PetscCall(VecSet(omega,1.0)); | |
| 180 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
| 181 | } | ||
| 182 | |||
| 183 | /*@ | ||
| 184 | SVDSetDSType - Sets the type of the internal `DS` object based on the current | ||
| 185 | settings of the singular value solver. | ||
| 186 | |||
| 187 | Collective | ||
| 188 | |||
| 189 | Input Parameter: | ||
| 190 | . svd - the singular value solver context | ||
| 191 | |||
| 192 | Note: | ||
| 193 | This function need not be called explicitly, since it will be called at | ||
| 194 | both `SVDSetFromOptions()` and `SVDSetUp()`. | ||
| 195 | |||
| 196 | Level: developer | ||
| 197 | |||
| 198 | .seealso: [](ch:svd), `SVDSetFromOptions()`, `SVDSetUp()` | ||
| 199 | @*/ | ||
| 200 | 4994 | PetscErrorCode SVDSetDSType(SVD svd) | |
| 201 | { | ||
| 202 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4994 | PetscFunctionBegin; |
| 203 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
4994 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
| 204 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
4994 | PetscTryTypeMethod(svd,setdstype); |
| 205 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
984 | PetscFunctionReturn(PETSC_SUCCESS); |
| 206 | } | ||
| 207 | |||
| 208 | /*@ | ||
| 209 | SVDSetUp - Sets up all the internal data structures necessary for the | ||
| 210 | execution of the singular value solver. | ||
| 211 | |||
| 212 | Collective | ||
| 213 | |||
| 214 | Input Parameter: | ||
| 215 | . svd - the singular value solver context | ||
| 216 | |||
| 217 | Notes: | ||
| 218 | This function need not be called explicitly in most cases, since `SVDSolve()` | ||
| 219 | calls it. It can be useful when one wants to measure the set-up time | ||
| 220 | separately from the solve time. | ||
| 221 | |||
| 222 | Level: developer | ||
| 223 | |||
| 224 | .seealso: [](ch:svd), `SVDCreate()`, `SVDSolve()`, `SVDDestroy()` | ||
| 225 | @*/ | ||
| 226 | 2773 | PetscErrorCode SVDSetUp(SVD svd) | |
| 227 | { | ||
| 228 | 2773 | PetscBool flg; | |
| 229 | 2773 | PetscInt M,N,P=0,k,maxnsol,m,Nom,nom; | |
| 230 | 2773 | SlepcSC sc; | |
| 231 | 2773 | Vec *T; | |
| 232 | 2773 | BV bv; | |
| 233 | 2773 | SVDStoppingCtx ctx; | |
| 234 | |||
| 235 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2773 | PetscFunctionBegin; |
| 236 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2773 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
| 237 |
2/14✓ Branch 0 taken 8 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
2773 | if (svd->state) PetscFunctionReturn(PETSC_SUCCESS); |
| 238 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(PetscLogEventBegin(SVD_SetUp,svd,0,0,0)); |
| 239 | |||
| 240 | /* reset the convergence flag from the previous solves */ | ||
| 241 | 2773 | svd->reason = SVD_CONVERGED_ITERATING; | |
| 242 | |||
| 243 | /* set default solver type (SVDSetFromOptions was not called) */ | ||
| 244 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
2773 | if (!((PetscObject)svd)->type_name) PetscCall(SVDSetType(svd,SVDCROSS)); |
| 245 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
2773 | if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds)); |
| 246 | |||
| 247 | /* check matrices */ | ||
| 248 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
2773 | PetscCheck(svd->OP,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONGSTATE,"SVDSetOperators() must be called first"); |
| 249 | |||
| 250 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | if (!svd->problem_type) { /* set default problem type */ |
| 251 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1336 | if (svd->OPb) { |
| 252 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
204 | PetscCheck(!svd->omega,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"There is no support yet for generalized hyperbolic problems"); |
| 253 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
204 | PetscCall(SVDSetProblemType(svd,SVD_GENERALIZED)); |
| 254 | } else { | ||
| 255 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1132 | if (svd->omega) PetscCall(SVDSetProblemType(svd,SVD_HYPERBOLIC)); |
| 256 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1122 | else PetscCall(SVDSetProblemType(svd,SVD_STANDARD)); |
| 257 | } | ||
| 258 | } else { /* check consistency of problem type set by user */ | ||
| 259 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1437 | if (svd->OPb) { |
| 260 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
781 | PetscCheck(svd->isgeneralized,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_INCOMP,"Inconsistent SVD state: the problem type does not match the number of matrices"); |
| 261 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
781 | PetscCheck(!svd->omega,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"There is no support yet for generalized hyperbolic problems"); |
| 262 | } else { | ||
| 263 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
656 | PetscCheck(!svd->isgeneralized,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_INCOMP,"Inconsistent SVD state: the problem type does not match the number of matrices"); |
| 264 |
3/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
656 | if (svd->omega) PetscCheck(svd->ishyperbolic,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_INCOMP,"Inconsistent SVD state: the problem type must be set to hyperbolic when passing a signature with SVDSetSignature()"); |
| 265 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
404 | else PetscCheck(!svd->ishyperbolic,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_INCOMP,"Inconsistent SVD state: a hyperbolic problem requires passing a signature with SVDSetSignature()"); |
| 266 | } | ||
| 267 | } | ||
| 268 | |||
| 269 | /* set DS type once the default problem type has been determined */ | ||
| 270 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(SVDSetDSType(svd)); |
| 271 | |||
| 272 | /* determine how to handle the transpose */ | ||
| 273 | 2773 | svd->expltrans = PETSC_TRUE; | |
| 274 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | if (svd->impltrans) svd->expltrans = PETSC_FALSE; |
| 275 | else { | ||
| 276 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2635 | PetscCall(MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg)); |
| 277 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2635 | if (!flg) svd->expltrans = PETSC_FALSE; |
| 278 | else { | ||
| 279 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2625 | PetscCall(PetscObjectTypeCompareAny((PetscObject)svd,&flg,SVDLAPACK,SVDSCALAPACK,SVDKSVD,SVDELEMENTAL,"")); |
| 280 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2625 | if (flg) svd->expltrans = PETSC_FALSE; |
| 281 | } | ||
| 282 | } | ||
| 283 | |||
| 284 | /* get matrix dimensions */ | ||
| 285 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(MatGetSize(svd->OP,&M,&N)); |
| 286 |
3/16✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2773 | if (M == 0 || N == 0) PetscFunctionReturn(PETSC_SUCCESS); |
| 287 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | if (svd->isgeneralized) { |
| 288 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
985 | PetscCall(MatGetSize(svd->OPb,&P,NULL)); |
| 289 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
985 | PetscCheck(M+P>=N,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"The case when [A;B] has less rows than columns is not supported"); |
| 290 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1788 | } else if (svd->ishyperbolic) { |
| 291 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
262 | PetscCall(VecGetSize(svd->omega,&Nom)); |
| 292 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
262 | PetscCall(VecGetLocalSize(svd->omega,&nom)); |
| 293 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
262 | PetscCall(MatGetLocalSize(svd->OP,&m,NULL)); |
| 294 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
262 | PetscCheck(Nom==M,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_SIZ,"Global size of signature (%" PetscInt_FMT ") does not match the row size of A (%" PetscInt_FMT ")",Nom,M); |
| 295 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
262 | PetscCheck(nom==m,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_SIZ,"Local size of signature (%" PetscInt_FMT ") does not match the local row size of A (%" PetscInt_FMT ")",nom,m); |
| 296 | } | ||
| 297 | |||
| 298 | /* build transpose matrix */ | ||
| 299 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(MatDestroy(&svd->A)); |
| 300 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(MatDestroy(&svd->AT)); |
| 301 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(PetscObjectReference((PetscObject)svd->OP)); |
| 302 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | if (svd->expltrans) { |
| 303 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
2346 | if (svd->isgeneralized || M>=N) { |
| 304 | 1978 | svd->A = svd->OP; | |
| 305 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1978 | PetscCall(MatHermitianTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->AT)); |
| 306 | } else { | ||
| 307 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
368 | PetscCall(MatHermitianTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->A)); |
| 308 | 368 | svd->AT = svd->OP; | |
| 309 | } | ||
| 310 | } else { | ||
| 311 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
427 | if (svd->isgeneralized || M>=N) { |
| 312 | 366 | svd->A = svd->OP; | |
| 313 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
366 | PetscCall(MatCreateHermitianTranspose(svd->OP,&svd->AT)); |
| 314 | } else { | ||
| 315 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
61 | PetscCall(MatCreateHermitianTranspose(svd->OP,&svd->A)); |
| 316 | 61 | svd->AT = svd->OP; | |
| 317 | } | ||
| 318 | } | ||
| 319 | |||
| 320 | /* build transpose matrix B for GSVD */ | ||
| 321 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | if (svd->isgeneralized) { |
| 322 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
985 | PetscCall(MatDestroy(&svd->B)); |
| 323 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
985 | PetscCall(MatDestroy(&svd->BT)); |
| 324 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
985 | PetscCall(PetscObjectReference((PetscObject)svd->OPb)); |
| 325 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
985 | if (svd->expltrans) { |
| 326 | 915 | svd->B = svd->OPb; | |
| 327 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
915 | PetscCall(MatHermitianTranspose(svd->OPb,MAT_INITIAL_MATRIX,&svd->BT)); |
| 328 | } else { | ||
| 329 | 70 | svd->B = svd->OPb; | |
| 330 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(MatCreateHermitianTranspose(svd->OPb,&svd->BT)); |
| 331 | } | ||
| 332 | } | ||
| 333 | |||
| 334 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
2773 | if (!svd->isgeneralized && M<N) { |
| 335 | /* swap initial vectors */ | ||
| 336 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
429 | if (svd->nini || svd->ninil) { |
| 337 | ✗ | T=svd->ISL; svd->ISL=svd->IS; svd->IS=T; | |
| 338 | ✗ | k=svd->ninil; svd->ninil=svd->nini; svd->nini=k; | |
| 339 | } | ||
| 340 | /* swap basis vectors */ | ||
| 341 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
429 | if (!svd->swapped) { /* only the first time in case of multiple calls */ |
| 342 | 359 | bv=svd->V; svd->V=svd->U; svd->U=bv; | |
| 343 | 359 | svd->swapped = PETSC_TRUE; | |
| 344 | } | ||
| 345 | } | ||
| 346 | |||
| 347 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | maxnsol = svd->isgeneralized? PetscMin(PetscMin(M,N),P): PetscMin(M,N); |
| 348 | 2773 | svd->ncv = PetscMin(svd->ncv,maxnsol); | |
| 349 | 2773 | svd->nsv = PetscMin(svd->nsv,maxnsol); | |
| 350 |
3/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
2773 | PetscCheck(svd->ncv==PETSC_DETERMINE || svd->nsv<=svd->ncv,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv"); |
| 351 | |||
| 352 | /* relative convergence criterion is not allowed in GSVD */ | ||
| 353 |
9/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
4016 | if (svd->conv==(SVDConv)-1) PetscCall(SVDSetConvergenceTest(svd,svd->isgeneralized?SVD_CONV_NORM:SVD_CONV_REL)); |
| 354 |
3/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
2773 | PetscCheck(!svd->isgeneralized || svd->conv!=SVD_CONV_REL,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Relative convergence criterion is not allowed in GSVD"); |
| 355 | |||
| 356 | /* initialization of matrix norm (standard case only, for GSVD it is done inside setup()) */ | ||
| 357 |
9/12✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 10 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 8 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2773 | if (!svd->isgeneralized && svd->conv==SVD_CONV_NORM && !svd->nrma) PetscCall(MatNorm(svd->OP,NORM_INFINITY,&svd->nrma)); |
| 358 | |||
| 359 | /* threshold stopping test */ | ||
| 360 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 10 times.
|
2773 | if (svd->stop==SVD_STOP_THRESHOLD) { |
| 361 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
96 | PetscCheck(svd->thres>0.0,PetscObjectComm((PetscObject)svd),PETSC_ERR_USER_INPUT,"Must give a threshold value with SVDSetThreshold()"); |
| 362 |
1/6✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
96 | PetscCheck(svd->which==SVD_LARGEST || !svd->threlative,PetscObjectComm((PetscObject)svd),PETSC_ERR_USER_INPUT,"Can only use a relative threshold when which=SVD_LARGEST"); |
| 363 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
96 | PetscCall(PetscNew(&ctx)); |
| 364 | 96 | ctx->thres = svd->thres; | |
| 365 | 96 | ctx->threlative = svd->threlative; | |
| 366 | 96 | ctx->which = svd->which; | |
| 367 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
96 | PetscCall(SVDSetStoppingTestFunction(svd,SVDStoppingThreshold,ctx,PetscCtxDestroyDefault)); |
| 368 | } | ||
| 369 | |||
| 370 | /* call specific solver setup */ | ||
| 371 |
5/10✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
2773 | PetscUseTypeMethod(svd,setup); |
| 372 | |||
| 373 | /* set tolerance if not yet set */ | ||
| 374 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | if (svd->tol==(PetscReal)PETSC_DETERMINE) svd->tol = SLEPC_DEFAULT_TOL; |
| 375 | |||
| 376 | /* fill sorting criterion context */ | ||
| 377 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(DSGetSlepcSC(svd->ds,&sc)); |
| 378 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | sc->comparison = (svd->which==SVD_LARGEST)? SlepcCompareLargestReal: SlepcCompareSmallestReal; |
| 379 | 2773 | sc->comparisonctx = NULL; | |
| 380 | 2773 | sc->map = NULL; | |
| 381 | 2773 | sc->mapobj = NULL; | |
| 382 | |||
| 383 | /* process initial vectors */ | ||
| 384 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | if (svd->nini<0) { |
| 385 | 274 | k = -svd->nini; | |
| 386 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
274 | PetscCheck(k<=svd->ncv,PetscObjectComm((PetscObject)svd),PETSC_ERR_USER_INPUT,"The number of initial vectors is larger than ncv"); |
| 387 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
274 | PetscCall(BVInsertVecs(svd->V,0,&k,svd->IS,PETSC_TRUE)); |
| 388 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
274 | PetscCall(SlepcBasisDestroy_Private(&svd->nini,&svd->IS)); |
| 389 | 274 | svd->nini = k; | |
| 390 | } | ||
| 391 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | if (svd->ninil<0) { |
| 392 | 331 | k = 0; | |
| 393 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
331 | if (svd->leftbasis) { |
| 394 | 250 | k = -svd->ninil; | |
| 395 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
250 | PetscCheck(k<=svd->ncv,PetscObjectComm((PetscObject)svd),PETSC_ERR_USER_INPUT,"The number of left initial vectors is larger than ncv"); |
| 396 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
250 | PetscCall(BVInsertVecs(svd->U,0,&k,svd->ISL,PETSC_TRUE)); |
| 397 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
81 | } else PetscCall(PetscInfo(svd,"Ignoring initial left vectors\n")); |
| 398 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
331 | PetscCall(SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL)); |
| 399 | 331 | svd->ninil = k; | |
| 400 | } | ||
| 401 | |||
| 402 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(PetscLogEventEnd(SVD_SetUp,svd,0,0,0)); |
| 403 | 2773 | svd->state = SVD_STATE_SETUP; | |
| 404 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2773 | PetscFunctionReturn(PETSC_SUCCESS); |
| 405 | } | ||
| 406 | |||
| 407 | /*@ | ||
| 408 | SVDSetInitialSpaces - Specify two bases of vectors that constitute the initial | ||
| 409 | right and/or left spaces. | ||
| 410 | |||
| 411 | Collective | ||
| 412 | |||
| 413 | Input Parameters: | ||
| 414 | + svd - the singular value solver context | ||
| 415 | . nr - number of right vectors | ||
| 416 | . isr - set of basis vectors of the right initial space | ||
| 417 | . nl - number of left vectors | ||
| 418 | - isl - set of basis vectors of the left initial space | ||
| 419 | |||
| 420 | Notes: | ||
| 421 | The initial right and left spaces are rough approximations to the right and/or | ||
| 422 | left singular subspaces from which the solver starts to iterate. | ||
| 423 | It is not necessary to provide both sets of vectors. | ||
| 424 | |||
| 425 | Some solvers start to iterate on a single vector (initial vector). In that case, | ||
| 426 | the other vectors are ignored. | ||
| 427 | |||
| 428 | These vectors do not persist from one `SVDSolve()` call to the other, so the | ||
| 429 | initial spaces should be set every time. | ||
| 430 | |||
| 431 | The vectors do not need to be mutually orthonormal, since they are explicitly | ||
| 432 | orthonormalized internally. | ||
| 433 | |||
| 434 | Common usage of this function is when the user can provide a rough approximation | ||
| 435 | of the wanted singular spaces. Then, convergence may be faster. | ||
| 436 | |||
| 437 | Level: intermediate | ||
| 438 | |||
| 439 | .seealso: [](ch:svd), `SVDSetUp()` | ||
| 440 | @*/ | ||
| 441 | 379 | PetscErrorCode SVDSetInitialSpaces(SVD svd,PetscInt nr,Vec isr[],PetscInt nl,Vec isl[]) | |
| 442 | { | ||
| 443 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
379 | PetscFunctionBegin; |
| 444 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
379 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
| 445 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
379 | PetscValidLogicalCollectiveInt(svd,nr,2); |
| 446 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
379 | PetscValidLogicalCollectiveInt(svd,nl,4); |
| 447 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
379 | PetscCheck(nr>=0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument nr cannot be negative"); |
| 448 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
379 | if (nr>0) { |
| 449 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
379 | PetscAssertPointer(isr,3); |
| 450 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
379 | PetscValidHeaderSpecific(*isr,VEC_CLASSID,3); |
| 451 | } | ||
| 452 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
379 | PetscCheck(nl>=0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument nl cannot be negative"); |
| 453 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
379 | if (nl>0) { |
| 454 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
379 | PetscAssertPointer(isl,5); |
| 455 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
379 | PetscValidHeaderSpecific(*isl,VEC_CLASSID,5); |
| 456 | } | ||
| 457 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
379 | PetscCall(SlepcBasisReference_Private(nr,isr,&svd->nini,&svd->IS)); |
| 458 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
379 | PetscCall(SlepcBasisReference_Private(nl,isl,&svd->ninil,&svd->ISL)); |
| 459 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
379 | if (nr>0 || nl>0) svd->state = SVD_STATE_INITIAL; |
| 460 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
68 | PetscFunctionReturn(PETSC_SUCCESS); |
| 461 | } | ||
| 462 | |||
| 463 | /* | ||
| 464 | SVDSetDimensions_Default - Set reasonable values for ncv, mpd if not set | ||
| 465 | by the user. This is called at setup. | ||
| 466 | */ | ||
| 467 | 1444 | PetscErrorCode SVDSetDimensions_Default(SVD svd) | |
| 468 | { | ||
| 469 | 1444 | PetscInt N,M,P,maxnsol; | |
| 470 | |||
| 471 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1444 | PetscFunctionBegin; |
| 472 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1444 | PetscCall(MatGetSize(svd->OP,&M,&N)); |
| 473 | 1444 | maxnsol = PetscMin(M,N); | |
| 474 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1444 | if (svd->isgeneralized) { |
| 475 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
624 | PetscCall(MatGetSize(svd->OPb,&P,NULL)); |
| 476 | 624 | maxnsol = PetscMin(maxnsol,P); | |
| 477 | } | ||
| 478 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
|
1444 | if (svd->nsv==0 && svd->stop!=SVD_STOP_THRESHOLD) svd->nsv = 1; |
| 479 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1444 | if (svd->ncv!=PETSC_DETERMINE) { /* ncv set */ |
| 480 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
600 | PetscCheck(svd->ncv>=svd->nsv,PetscObjectComm((PetscObject)svd),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nsv"); |
| 481 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 10 times.
|
844 | } else if (svd->mpd!=PETSC_DETERMINE) { /* mpd set */ |
| 482 | 12 | svd->ncv = PetscMin(maxnsol,svd->nsv+svd->mpd); | |
| 483 | } else { /* neither set: defaults depend on nsv being small or large */ | ||
| 484 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
832 | if (svd->nsv<500) svd->ncv = PetscMin(maxnsol,PetscMax(2*svd->nsv,10)); |
| 485 | else { | ||
| 486 | ✗ | svd->mpd = 500; | |
| 487 | ✗ | svd->ncv = PetscMin(maxnsol,svd->nsv+svd->mpd); | |
| 488 | } | ||
| 489 | } | ||
| 490 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1444 | if (svd->mpd==PETSC_DETERMINE) svd->mpd = svd->ncv; |
| 491 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
292 | PetscFunctionReturn(PETSC_SUCCESS); |
| 492 | } | ||
| 493 | |||
| 494 | /*@ | ||
| 495 | SVDAllocateSolution - Allocate memory storage for common variables such | ||
| 496 | as the singular values and the basis vectors. | ||
| 497 | |||
| 498 | Collective | ||
| 499 | |||
| 500 | Input Parameters: | ||
| 501 | + svd - the singular value solver context | ||
| 502 | - extra - number of additional positions, used for methods that require a | ||
| 503 | working basis slightly larger than `ncv` | ||
| 504 | |||
| 505 | Developer Notes: | ||
| 506 | This is `SLEPC_EXTERN` because it may be required by user plugin `SVD` | ||
| 507 | implementations. | ||
| 508 | |||
| 509 | This is called at setup after setting the value of `ncv` and the flag `leftbasis`. | ||
| 510 | |||
| 511 | Level: developer | ||
| 512 | |||
| 513 | .seealso: [](ch:svd), `SVDSetUp()`, `SVDSetDimensions()` | ||
| 514 | @*/ | ||
| 515 | 2773 | PetscErrorCode SVDAllocateSolution(SVD svd,PetscInt extra) | |
| 516 | { | ||
| 517 | 2773 | PetscInt oldsize,requested; | |
| 518 | 2773 | Vec tr,tl; | |
| 519 | |||
| 520 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2773 | PetscFunctionBegin; |
| 521 | 2773 | requested = svd->ncv + extra; | |
| 522 | |||
| 523 | /* oldsize is zero if this is the first time setup is called */ | ||
| 524 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2773 | PetscCall(BVGetSizes(svd->V,NULL,NULL,&oldsize)); |
| 525 | |||
| 526 | /* allocate sigma */ | ||
| 527 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
2773 | if (requested != oldsize || !svd->sigma) { |
| 528 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2411 | PetscCall(PetscFree3(svd->sigma,svd->perm,svd->errest)); |
| 529 |
1/10✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
2411 | if (svd->sign) PetscCall(PetscFree(svd->sign)); |
| 530 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2411 | PetscCall(PetscMalloc3(requested,&svd->sigma,requested,&svd->perm,requested,&svd->errest)); |
| 531 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
2411 | if (svd->ishyperbolic) PetscCall(PetscMalloc1(requested,&svd->sign)); |
| 532 | } | ||
| 533 | /* allocate V */ | ||
| 534 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
2773 | if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,NULL)); |
| 535 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2773 | if (!oldsize) { |
| 536 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
2313 | if (!((PetscObject)svd->V)->type_name) PetscCall(BVSetType(svd->V,BVMAT)); |
| 537 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2313 | PetscCall(MatCreateVecsEmpty(svd->A,&tr,NULL)); |
| 538 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2313 | PetscCall(BVSetSizesFromVec(svd->V,tr,requested)); |
| 539 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2313 | PetscCall(VecDestroy(&tr)); |
| 540 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
460 | } else PetscCall(BVResize(svd->V,requested,PETSC_FALSE)); |
| 541 | /* allocate U */ | ||
| 542 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
2773 | if (svd->leftbasis && !svd->isgeneralized) { |
| 543 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1336 | if (!svd->U) PetscCall(SVDGetBV(svd,NULL,&svd->U)); |
| 544 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1336 | if (!oldsize) { |
| 545 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1090 | if (!((PetscObject)svd->U)->type_name) PetscCall(BVSetType(svd->U,((PetscObject)svd->V)->type_name)); |
| 546 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1090 | PetscCall(MatCreateVecsEmpty(svd->A,NULL,&tl)); |
| 547 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1090 | PetscCall(BVSetSizesFromVec(svd->U,tl,requested)); |
| 548 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1090 | PetscCall(VecDestroy(&tl)); |
| 549 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
246 | } else PetscCall(BVResize(svd->U,requested,PETSC_FALSE)); |
| 550 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1437 | } else if (svd->isgeneralized) { /* left basis for the GSVD */ |
| 551 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
985 | if (!svd->U) PetscCall(SVDGetBV(svd,NULL,&svd->U)); |
| 552 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
985 | if (!oldsize) { |
| 553 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
851 | if (!((PetscObject)svd->U)->type_name) PetscCall(BVSetType(svd->U,((PetscObject)svd->V)->type_name)); |
| 554 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
851 | PetscCall(SVDCreateLeftTemplate(svd,&tl)); |
| 555 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
851 | PetscCall(BVSetSizesFromVec(svd->U,tl,requested)); |
| 556 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
851 | PetscCall(VecDestroy(&tl)); |
| 557 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
134 | } else PetscCall(BVResize(svd->U,requested,PETSC_FALSE)); |
| 558 | } | ||
| 559 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
547 | PetscFunctionReturn(PETSC_SUCCESS); |
| 560 | } | ||
| 561 | |||
| 562 | /*@ | ||
| 563 | SVDReallocateSolution - Reallocate memory storage for common variables such | ||
| 564 | as the singular values and the basis vectors. | ||
| 565 | |||
| 566 | Collective | ||
| 567 | |||
| 568 | Input Parameters: | ||
| 569 | + svd - the singular value solver context | ||
| 570 | - newsize - new size | ||
| 571 | |||
| 572 | Developer Notes: | ||
| 573 | This is `SLEPC_EXTERN` because it may be required by user plugin SVD | ||
| 574 | implementations. | ||
| 575 | |||
| 576 | This is called during the iteration in case the threshold stopping test has | ||
| 577 | been selected. | ||
| 578 | |||
| 579 | Level: developer | ||
| 580 | |||
| 581 | .seealso: [](ch:svd), `SVDAllocateSolution()`, `SVDSetThreshold()` | ||
| 582 | @*/ | ||
| 583 | 47 | PetscErrorCode SVDReallocateSolution(SVD svd,PetscInt newsize) | |
| 584 | { | ||
| 585 | 47 | PetscInt oldsize,*nperm; | |
| 586 | 47 | PetscReal *nsigma,*nerrest,*nsign; | |
| 587 | |||
| 588 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
47 | PetscFunctionBegin; |
| 589 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
47 | PetscCall(BVGetSizes(svd->V,NULL,NULL,&oldsize)); |
| 590 |
2/14✓ Branch 0 taken 6 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
47 | if (oldsize>=newsize) PetscFunctionReturn(PETSC_SUCCESS); |
| 591 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
47 | PetscCall(PetscInfo(svd,"Reallocating basis vectors to %" PetscInt_FMT " columns\n",newsize)); |
| 592 | |||
| 593 | /* reallocate sigma */ | ||
| 594 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
47 | PetscCall(PetscMalloc3(newsize,&nsigma,newsize,&nperm,newsize,&nerrest)); |
| 595 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
47 | PetscCall(PetscArraycpy(nsigma,svd->sigma,oldsize)); |
| 596 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
47 | PetscCall(PetscArraycpy(nperm,svd->perm,oldsize)); |
| 597 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
47 | PetscCall(PetscArraycpy(nerrest,svd->errest,oldsize)); |
| 598 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
47 | PetscCall(PetscFree3(svd->sigma,svd->perm,svd->errest)); |
| 599 | 47 | svd->sigma = nsigma; | |
| 600 | 47 | svd->perm = nperm; | |
| 601 | 47 | svd->errest = nerrest; | |
| 602 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
|
47 | if (svd->ishyperbolic) { |
| 603 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18 | PetscCall(PetscMalloc1(newsize,&nsign)); |
| 604 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18 | PetscCall(PetscArraycpy(nsign,svd->sign,oldsize)); |
| 605 |
5/8✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
18 | PetscCall(PetscFree(svd->sign)); |
| 606 | 18 | svd->sign = nsign; | |
| 607 | } | ||
| 608 | /* reallocate V,U */ | ||
| 609 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
47 | PetscCall(BVResize(svd->V,newsize,PETSC_TRUE)); |
| 610 |
7/10✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 6 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
47 | if (svd->leftbasis || svd->isgeneralized) PetscCall(BVResize(svd->U,newsize,PETSC_TRUE)); |
| 611 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
11 | PetscFunctionReturn(PETSC_SUCCESS); |
| 612 | } | ||
| 613 |