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 | This file implements a wrapper to the PRIMME SVD solver | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/ | ||
15 | |||
16 | #include <primme.h> | ||
17 | |||
18 | #if defined(PETSC_USE_COMPLEX) | ||
19 | #if defined(PETSC_USE_REAL_SINGLE) | ||
20 | #define PRIMME_DRIVER cprimme_svds | ||
21 | #else | ||
22 | #define PRIMME_DRIVER zprimme_svds | ||
23 | #endif | ||
24 | #else | ||
25 | #if defined(PETSC_USE_REAL_SINGLE) | ||
26 | #define PRIMME_DRIVER sprimme_svds | ||
27 | #else | ||
28 | #define PRIMME_DRIVER dprimme_svds | ||
29 | #endif | ||
30 | #endif | ||
31 | |||
32 | #if defined(PRIMME_VERSION_MAJOR) && PRIMME_VERSION_MAJOR*100+PRIMME_VERSION_MINOR >= 202 | ||
33 | #define SLEPC_HAVE_PRIMME2p2 | ||
34 | #endif | ||
35 | |||
36 | typedef struct { | ||
37 | primme_svds_params primme; /* param struct */ | ||
38 | PetscInt bs; /* block size */ | ||
39 | primme_svds_preset_method method; /* primme method */ | ||
40 | SVD svd; /* reference to the solver */ | ||
41 | Vec x,y; /* auxiliary vectors */ | ||
42 | } SVD_PRIMME; | ||
43 | |||
44 | static void multMatvec_PRIMME(void*,PRIMME_INT*,void*,PRIMME_INT*,int*,int*,struct primme_svds_params*,int*); | ||
45 | |||
46 | 13786 | static void par_GlobalSumReal(void *sendBuf,void *recvBuf,int *count,primme_svds_params *primme,int *ierr) | |
47 | { | ||
48 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
13786 | if (sendBuf == recvBuf) { |
49 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
27572 | *ierr = MPI_Allreduce(MPI_IN_PLACE,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo)); |
50 | } else { | ||
51 | ✗ | *ierr = MPI_Allreduce(sendBuf,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo)); | |
52 | } | ||
53 | 13786 | } | |
54 | |||
55 | #if defined(SLEPC_HAVE_PRIMME3) | ||
56 | 7376 | static void par_broadcastReal(void *buf,int *count,primme_svds_params *primme,int *ierr) | |
57 | { | ||
58 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
7376 | *ierr = MPI_Bcast(buf,*count,MPIU_REAL,0/*root*/,PetscObjectComm((PetscObject)primme->commInfo)); |
59 | 7376 | } | |
60 | #endif | ||
61 | |||
62 | #if defined(SLEPC_HAVE_PRIMME2p2) | ||
63 | 5736 | static void convTestFun(double *sval,void *leftsvec,void *rightsvec,double *resNorm, | |
64 | #if defined(SLEPC_HAVE_PRIMME3) | ||
65 | int *method, | ||
66 | #endif | ||
67 | int *isconv,struct primme_svds_params *primme,int *err) | ||
68 | { | ||
69 | 5736 | PetscErrorCode ierr; | |
70 | 5736 | SVD svd = (SVD)primme->commInfo; | |
71 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
5736 | PetscReal sigma = sval?*sval:0.0; |
72 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
5736 | PetscReal r = resNorm?*resNorm:HUGE_VAL,errest; |
73 | |||
74 | 5736 | ierr = (*svd->converged)(svd,sigma,r,&errest,svd->convergedctx); | |
75 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
5736 | if (ierr) *err = 1; |
76 | else { | ||
77 | 5736 | *isconv = (errest<=svd->tol?1:0); | |
78 | 5736 | *err = 0; | |
79 | } | ||
80 | 5736 | } | |
81 | |||
82 | 4639 | static void monitorFun(void *basisSvals,int *basisSize,int *basisFlags,int *iblock,int *blockSize,void *basisNorms,int *numConverged,void *lockedSvals,int *numLocked,int *lockedFlags,void *lockedNorms,int *inner_its,void *LSRes, | |
83 | #if defined(SLEPC_HAVE_PRIMME3) | ||
84 | const char *msg,double *time, | ||
85 | #endif | ||
86 | primme_event *event,int *stage,struct primme_svds_params *primme,int *err) | ||
87 | { | ||
88 | |||
89 | 4639 | PetscErrorCode ierr = PETSC_SUCCESS; | |
90 | 4639 | SVD svd = (SVD)primme->commInfo; | |
91 | 4639 | PetscInt i,k,nerrest; | |
92 | |||
93 | 4639 | *err = 1; | |
94 |
2/3✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 6 times.
|
4639 | switch (*event) { |
95 | 845 | case primme_event_outer_iteration: | |
96 | /* Update SVD */ | ||
97 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
845 | PetscCallVoid(PetscIntCast(primme->stats.numOuterIterations,&svd->its)); |
98 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
845 | if (numConverged) svd->nconv = *numConverged; |
99 | 845 | k = 0; | |
100 |
1/6✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
845 | if (lockedSvals && numLocked) for (i=0; i<*numLocked && k<svd->ncv; i++) svd->sigma[k++] = ((PetscReal*)lockedSvals)[i]; |
101 | 845 | nerrest = k; | |
102 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
845 | if (iblock && blockSize) { |
103 |
3/4✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
|
2308 | for (i=0; i<*blockSize && k+iblock[i]<svd->ncv; i++) svd->errest[k+iblock[i]] = ((PetscReal*)basisNorms)[i]; |
104 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
845 | nerrest = k+(*blockSize>0?1+iblock[*blockSize-1]:0); |
105 | } | ||
106 |
4/6✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 6 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 6 times.
✗ Branch 5 not taken.
|
6296 | if (basisSvals && basisSize) for (i=0; i<*basisSize && k<svd->ncv; i++) svd->sigma[k++] = ((PetscReal*)basisSvals)[i]; |
107 | /* Show progress */ | ||
108 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
845 | ierr = SVDMonitor(svd,svd->its,numConverged?*numConverged:0,svd->sigma,svd->errest,nerrest); |
109 | 845 | break; | |
110 | #if defined(SLEPC_HAVE_PRIMME3) | ||
111 | ✗ | case primme_event_message: | |
112 | /* Print PRIMME information messages */ | ||
113 | ✗ | ierr = PetscInfo(svd,"%s\n",msg); | |
114 | ✗ | break; | |
115 | #endif | ||
116 | default: | ||
117 | break; | ||
118 | } | ||
119 | 4639 | *err = (ierr!=0)? 1: 0; | |
120 | } | ||
121 | #endif /* SLEPC_HAVE_PRIMME2p2 */ | ||
122 | |||
123 | 8466 | static void multMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,int *transpose,struct primme_svds_params *primme,int *ierr) | |
124 | { | ||
125 | 8466 | PetscInt i; | |
126 | 8466 | SVD_PRIMME *ops = (SVD_PRIMME*)primme->matrix; | |
127 | 8466 | Vec x = ops->x,y = ops->y; | |
128 | 8466 | SVD svd = ops->svd; | |
129 | |||
130 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
8466 | PetscFunctionBegin; |
131 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
21634 | for (i=0;i<*blockSize;i++) { |
132 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
13168 | if (*transpose) { |
133 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6500 | PetscCallAbort(PetscObjectComm((PetscObject)svd),VecPlaceArray(y,(PetscScalar*)xa+(*ldx)*i)); |
134 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6500 | PetscCallAbort(PetscObjectComm((PetscObject)svd),VecPlaceArray(x,(PetscScalar*)ya+(*ldy)*i)); |
135 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6500 | PetscCallAbort(PetscObjectComm((PetscObject)svd),MatMult(svd->AT,y,x)); |
136 | } else { | ||
137 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6668 | PetscCallAbort(PetscObjectComm((PetscObject)svd),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i)); |
138 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6668 | PetscCallAbort(PetscObjectComm((PetscObject)svd),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i)); |
139 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6668 | PetscCallAbort(PetscObjectComm((PetscObject)svd),MatMult(svd->A,x,y)); |
140 | } | ||
141 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13168 | PetscCallAbort(PetscObjectComm((PetscObject)svd),VecResetArray(x)); |
142 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13168 | PetscCallAbort(PetscObjectComm((PetscObject)svd),VecResetArray(y)); |
143 | } | ||
144 |
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.
|
8466 | PetscFunctionReturnVoid(); |
145 | } | ||
146 | |||
147 | 42 | static PetscErrorCode SVDSetUp_PRIMME(SVD svd) | |
148 | { | ||
149 | 42 | PetscMPIInt numProcs,procID; | |
150 | 42 | PetscInt n,m,nloc,mloc; | |
151 | 42 | SVD_PRIMME *ops = (SVD_PRIMME*)svd->data; | |
152 | 42 | primme_svds_params *primme = &ops->primme; | |
153 | |||
154 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
42 | PetscFunctionBegin; |
155 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
42 | SVDCheckStandard(svd); |
156 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
42 | SVDCheckDefinite(svd); |
157 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 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.
✓ 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.
|
42 | PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)svd),&numProcs)); |
158 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 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.
✓ 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.
|
42 | PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&procID)); |
159 | |||
160 | /* Check some constraints and set some default values */ | ||
161 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(MatGetSize(svd->A,&m,&n)); |
162 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(MatGetLocalSize(svd->A,&mloc,&nloc)); |
163 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(SVDSetDimensions_Default(svd)); |
164 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
42 | if (svd->max_it==PETSC_DETERMINE) svd->max_it = PETSC_INT_MAX; |
165 | 42 | svd->leftbasis = PETSC_TRUE; | |
166 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
42 | SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING); |
167 | #if !defined(SLEPC_HAVE_PRIMME2p2) | ||
168 | if (svd->converged != SVDConvergedAbsolute) PetscCall(PetscInfo(svd,"Warning: using absolute convergence test\n")); | ||
169 | #endif | ||
170 | |||
171 | /* Transfer SLEPc options to PRIMME options */ | ||
172 | 42 | primme_svds_free(primme); | |
173 | 42 | primme_svds_initialize(primme); | |
174 | 42 | primme->m = (PRIMME_INT)m; | |
175 | 42 | primme->n = (PRIMME_INT)n; | |
176 | 42 | primme->mLocal = (PRIMME_INT)mloc; | |
177 | 42 | primme->nLocal = (PRIMME_INT)nloc; | |
178 | 42 | primme->numSvals = (int)svd->nsv; | |
179 | 42 | primme->matrix = ops; | |
180 | 42 | primme->commInfo = svd; | |
181 | 42 | primme->maxMatvecs = (PRIMME_INT)svd->max_it; | |
182 | #if !defined(SLEPC_HAVE_PRIMME2p2) | ||
183 | primme->eps = SlepcDefaultTol(svd->tol); | ||
184 | #endif | ||
185 | 42 | primme->numProcs = numProcs; | |
186 | 42 | primme->procID = procID; | |
187 | 42 | primme->printLevel = 1; | |
188 | 42 | primme->matrixMatvec = multMatvec_PRIMME; | |
189 | 42 | primme->globalSumReal = par_GlobalSumReal; | |
190 | #if defined(SLEPC_HAVE_PRIMME3) | ||
191 | 42 | primme->broadcastReal = par_broadcastReal; | |
192 | #endif | ||
193 | #if defined(SLEPC_HAVE_PRIMME2p2) | ||
194 | 42 | primme->convTestFun = convTestFun; | |
195 | 42 | primme->monitorFun = monitorFun; | |
196 | #endif | ||
197 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
42 | if (ops->bs > 0) primme->maxBlockSize = (int)ops->bs; |
198 | |||
199 |
2/3✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
|
42 | switch (svd->which) { |
200 | 36 | case SVD_LARGEST: | |
201 | 36 | primme->target = primme_svds_largest; | |
202 | 36 | break; | |
203 | 6 | case SVD_SMALLEST: | |
204 | 6 | primme->target = primme_svds_smallest; | |
205 | 6 | break; | |
206 | } | ||
207 | |||
208 | /* If user sets mpd or ncv, maxBasisSize is modified */ | ||
209 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
42 | if (svd->mpd!=PETSC_DETERMINE) { |
210 | 42 | primme->maxBasisSize = (int)svd->mpd; | |
211 |
5/8✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
42 | if (svd->ncv!=PETSC_DETERMINE) PetscCall(PetscInfo(svd,"Warning: 'ncv' is ignored by PRIMME\n")); |
212 | ✗ | } else if (svd->ncv!=PETSC_DETERMINE) primme->maxBasisSize = (int)svd->ncv; | |
213 | |||
214 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
42 | PetscCheck(primme_svds_set_method(ops->method,(primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME,PRIMME_DEFAULT_METHOD,primme)>=0,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"PRIMME method not valid"); |
215 | |||
216 | 42 | svd->mpd = (PetscInt)primme->maxBasisSize; | |
217 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
42 | svd->ncv = (PetscInt)(primme->locking?svd->nsv:0)+primme->maxBasisSize; |
218 | 42 | ops->bs = (PetscInt)primme->maxBlockSize; | |
219 | |||
220 | /* Set workspace */ | ||
221 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(SVDAllocateSolution(svd,0)); |
222 | |||
223 | /* Prepare auxiliary vectors */ | ||
224 |
6/8✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
42 | if (!ops->x) PetscCall(MatCreateVecsEmpty(svd->A,&ops->x,&ops->y)); |
225 |
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.
|
14 | PetscFunctionReturn(PETSC_SUCCESS); |
226 | } | ||
227 | |||
228 | 42 | static PetscErrorCode SVDSolve_PRIMME(SVD svd) | |
229 | { | ||
230 | 42 | SVD_PRIMME *ops = (SVD_PRIMME*)svd->data; | |
231 | 42 | PetscScalar *svecs, *a; | |
232 | 42 | PetscInt i,ierrprimme,ld; | |
233 | 42 | PetscReal *svals,*rnorms; | |
234 | |||
235 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
42 | PetscFunctionBegin; |
236 | /* Reset some parameters left from previous runs */ | ||
237 | 42 | ops->primme.aNorm = 0.0; | |
238 | 42 | ops->primme.initSize = (int)svd->nini; | |
239 | 42 | ops->primme.iseed[0] = -1; | |
240 | 42 | ops->primme.iseed[1] = -1; | |
241 | 42 | ops->primme.iseed[2] = -1; | |
242 | 42 | ops->primme.iseed[3] = -1; | |
243 | |||
244 | /* Allocating left and right singular vectors contiguously */ | ||
245 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(PetscCalloc1(ops->primme.numSvals*(ops->primme.mLocal+ops->primme.nLocal),&svecs)); |
246 | |||
247 | /* Call PRIMME solver */ | ||
248 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(PetscMalloc2(svd->ncv,&svals,svd->ncv,&rnorms)); |
249 | 42 | ierrprimme = PRIMME_DRIVER(svals,svecs,rnorms,&ops->primme); | |
250 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
684 | for (i=0;i<svd->ncv;i++) svd->sigma[i] = svals[i]; |
251 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
642 | for (i=0;i<svd->ncv;i++) svd->errest[i] = rnorms[i]; |
252 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(PetscFree2(svals,rnorms)); |
253 | |||
254 | /* Copy left and right singular vectors into svd */ | ||
255 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(BVGetLeadingDimension(svd->U,&ld)); |
256 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(BVGetArray(svd->U,&a)); |
257 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 6 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
210 | for (i=0;i<ops->primme.initSize;i++) PetscCall(PetscArraycpy(a+i*ld,svecs+i*ops->primme.mLocal,ops->primme.mLocal)); |
258 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(BVRestoreArray(svd->U,&a)); |
259 | |||
260 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(BVGetLeadingDimension(svd->V,&ld)); |
261 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(BVGetArray(svd->V,&a)); |
262 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 6 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
210 | for (i=0;i<ops->primme.initSize;i++) PetscCall(PetscArraycpy(a+i*ld,svecs+ops->primme.mLocal*ops->primme.initSize+i*ops->primme.nLocal,ops->primme.nLocal)); |
263 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(BVRestoreArray(svd->V,&a)); |
264 | |||
265 |
5/8✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
42 | PetscCall(PetscFree(svecs)); |
266 | |||
267 | 42 | svd->nconv = ops->primme.initSize >= 0 ? (PetscInt)ops->primme.initSize : 0; | |
268 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
42 | svd->reason = svd->nconv >= svd->nsv ? SVD_CONVERGED_TOL: SVD_DIVERGED_ITS; |
269 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(PetscIntCast(ops->primme.stats.numOuterIterations,&svd->its)); |
270 | |||
271 | /* Process PRIMME error code */ | ||
272 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
42 | if (ierrprimme != 0) { |
273 | ✗ | switch (ierrprimme%100) { | |
274 | ✗ | case -1: | |
275 | ✗ | SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": unexpected error",ierrprimme); | |
276 | ✗ | case -2: | |
277 | ✗ | SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": allocation error",ierrprimme); | |
278 | case -3: /* stop due to maximum number of iterations or matvecs */ | ||
279 | break; | ||
280 | ✗ | default: | |
281 | ✗ | PetscCheck(ierrprimme<-39,PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": configuration error; check PRIMME's manual",ierrprimme); | |
282 | ✗ | PetscCheck(ierrprimme>=-39,PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": runtime error; check PRIMME's manual",ierrprimme); | |
283 | } | ||
284 | } | ||
285 |
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.
|
14 | PetscFunctionReturn(PETSC_SUCCESS); |
286 | } | ||
287 | |||
288 | 36 | static PetscErrorCode SVDReset_PRIMME(SVD svd) | |
289 | { | ||
290 | 36 | SVD_PRIMME *ops = (SVD_PRIMME*)svd->data; | |
291 | |||
292 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | PetscFunctionBegin; |
293 | 36 | primme_svds_free(&ops->primme); | |
294 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(VecDestroy(&ops->x)); |
295 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(VecDestroy(&ops->y)); |
296 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
297 | } | ||
298 | |||
299 | 36 | static PetscErrorCode SVDDestroy_PRIMME(SVD svd) | |
300 | { | ||
301 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | PetscFunctionBegin; |
302 |
5/8✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
36 | PetscCall(PetscFree(svd->data)); |
303 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",NULL)); |
304 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",NULL)); |
305 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",NULL)); |
306 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",NULL)); |
307 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
308 | } | ||
309 | |||
310 | ✗ | static PetscErrorCode SVDView_PRIMME(SVD svd,PetscViewer viewer) | |
311 | { | ||
312 | ✗ | PetscBool isascii; | |
313 | ✗ | SVD_PRIMME *ctx = (SVD_PRIMME*)svd->data; | |
314 | ✗ | PetscMPIInt rank; | |
315 | |||
316 | ✗ | PetscFunctionBegin; | |
317 | ✗ | PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); | |
318 | ✗ | if (isascii) { | |
319 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer," block size=%" PetscInt_FMT "\n",ctx->bs)); | |
320 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer," solver method: %s\n",SVDPRIMMEMethods[(SVDPRIMMEMethod)ctx->method])); | |
321 | |||
322 | /* Display PRIMME params */ | ||
323 | ✗ | PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&rank)); | |
324 | ✗ | if (!rank) primme_svds_display_params(ctx->primme); | |
325 | } | ||
326 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
327 | } | ||
328 | |||
329 | 30 | static PetscErrorCode SVDSetFromOptions_PRIMME(SVD svd,PetscOptionItems PetscOptionsObject) | |
330 | { | ||
331 | 30 | SVD_PRIMME *ctx = (SVD_PRIMME*)svd->data; | |
332 | 30 | PetscInt bs; | |
333 | 30 | SVDPRIMMEMethod meth; | |
334 | 30 | PetscBool flg; | |
335 | |||
336 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
337 |
1/12✗ Branch 0 not taken.
✓ Branch 1 taken 6 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.
|
30 | PetscOptionsHeadBegin(PetscOptionsObject,"SVD PRIMME Options"); |
338 | |||
339 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(PetscOptionsInt("-svd_primme_blocksize","Maximum block size","SVDPRIMMESetBlockSize",ctx->bs,&bs,&flg)); |
340 |
6/8✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
30 | if (flg) PetscCall(SVDPRIMMESetBlockSize(svd,bs)); |
341 | |||
342 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(PetscOptionsEnum("-svd_primme_method","Method for solving the singular value problem","SVDPRIMMESetMethod",SVDPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg)); |
343 |
6/8✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
30 | if (flg) PetscCall(SVDPRIMMESetMethod(svd,meth)); |
344 | |||
345 |
1/14✗ Branch 0 not taken.
✓ 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.
|
30 | PetscOptionsHeadEnd(); |
346 |
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.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
347 | } | ||
348 | |||
349 | 18 | static PetscErrorCode SVDPRIMMESetBlockSize_PRIMME(SVD svd,PetscInt bs) | |
350 | { | ||
351 | 18 | SVD_PRIMME *ops = (SVD_PRIMME*)svd->data; | |
352 | |||
353 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
18 | PetscFunctionBegin; |
354 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | if (bs == PETSC_DEFAULT || bs == PETSC_DECIDE) ops->bs = 0; |
355 | else { | ||
356 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
18 | PetscCheck(bs>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive"); |
357 | 18 | ops->bs = bs; | |
358 | } | ||
359 |
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.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
360 | } | ||
361 | |||
362 | /*@ | ||
363 | SVDPRIMMESetBlockSize - The maximum block size that PRIMME will try to use. | ||
364 | |||
365 | Logically Collective | ||
366 | |||
367 | Input Parameters: | ||
368 | + svd - the singular value solver context | ||
369 | - bs - block size | ||
370 | |||
371 | Options Database Key: | ||
372 | . -svd_primme_blocksize - Sets the max allowed block size value | ||
373 | |||
374 | Notes: | ||
375 | If the block size is not set, the value established by primme_svds_initialize | ||
376 | is used. | ||
377 | |||
378 | The user should set the block size based on the architecture specifics | ||
379 | of the target computer, as well as any a priori knowledge of multiplicities. | ||
380 | The code does NOT require bs > 1 to find multiple eigenvalues. For some | ||
381 | methods, keeping bs = 1 yields the best overall performance. | ||
382 | |||
383 | Level: advanced | ||
384 | |||
385 | .seealso: SVDPRIMMEGetBlockSize() | ||
386 | @*/ | ||
387 | 18 | PetscErrorCode SVDPRIMMESetBlockSize(SVD svd,PetscInt bs) | |
388 | { | ||
389 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
18 | PetscFunctionBegin; |
390 |
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.
|
18 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
391 |
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.
|
18 | PetscValidLogicalCollectiveInt(svd,bs,2); |
392 |
8/14✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
✓ 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.
|
18 | PetscTryMethod(svd,"SVDPRIMMESetBlockSize_C",(SVD,PetscInt),(svd,bs)); |
393 |
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.
|
18 | PetscFunctionReturn(PETSC_SUCCESS); |
394 | } | ||
395 | |||
396 | 12 | static PetscErrorCode SVDPRIMMEGetBlockSize_PRIMME(SVD svd,PetscInt *bs) | |
397 | { | ||
398 | 12 | SVD_PRIMME *ops = (SVD_PRIMME*)svd->data; | |
399 | |||
400 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
12 | PetscFunctionBegin; |
401 | 12 | *bs = ops->bs; | |
402 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
403 | } | ||
404 | |||
405 | /*@ | ||
406 | SVDPRIMMEGetBlockSize - Get the maximum block size the code will try to use. | ||
407 | |||
408 | Not Collective | ||
409 | |||
410 | Input Parameter: | ||
411 | . svd - the singular value solver context | ||
412 | |||
413 | Output Parameter: | ||
414 | . bs - returned block size | ||
415 | |||
416 | Level: advanced | ||
417 | |||
418 | .seealso: SVDPRIMMESetBlockSize() | ||
419 | @*/ | ||
420 | 12 | PetscErrorCode SVDPRIMMEGetBlockSize(SVD svd,PetscInt *bs) | |
421 | { | ||
422 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
12 | PetscFunctionBegin; |
423 |
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.
|
12 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
424 |
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.
|
12 | PetscAssertPointer(bs,2); |
425 |
9/16✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 6 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
|
12 | PetscUseMethod(svd,"SVDPRIMMEGetBlockSize_C",(SVD,PetscInt*),(svd,bs)); |
426 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
427 | } | ||
428 | |||
429 | 18 | static PetscErrorCode SVDPRIMMESetMethod_PRIMME(SVD svd,SVDPRIMMEMethod method) | |
430 | { | ||
431 | 18 | SVD_PRIMME *ops = (SVD_PRIMME*)svd->data; | |
432 | |||
433 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
18 | PetscFunctionBegin; |
434 | 18 | ops->method = (primme_svds_preset_method)method; | |
435 |
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.
|
18 | PetscFunctionReturn(PETSC_SUCCESS); |
436 | } | ||
437 | |||
438 | /*@ | ||
439 | SVDPRIMMESetMethod - Sets the method for the PRIMME SVD solver. | ||
440 | |||
441 | Logically Collective | ||
442 | |||
443 | Input Parameters: | ||
444 | + svd - the singular value solver context | ||
445 | - method - method that will be used by PRIMME | ||
446 | |||
447 | Options Database Key: | ||
448 | . -svd_primme_method - Sets the method for the PRIMME SVD solver | ||
449 | |||
450 | Note: | ||
451 | If not set, the method defaults to SVD_PRIMME_HYBRID. | ||
452 | |||
453 | Level: advanced | ||
454 | |||
455 | .seealso: SVDPRIMMEGetMethod(), SVDPRIMMEMethod | ||
456 | @*/ | ||
457 | 18 | PetscErrorCode SVDPRIMMESetMethod(SVD svd,SVDPRIMMEMethod method) | |
458 | { | ||
459 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
18 | PetscFunctionBegin; |
460 |
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.
|
18 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
461 |
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.
|
18 | PetscValidLogicalCollectiveEnum(svd,method,2); |
462 |
8/14✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
✓ 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.
|
18 | PetscTryMethod(svd,"SVDPRIMMESetMethod_C",(SVD,SVDPRIMMEMethod),(svd,method)); |
463 |
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.
|
18 | PetscFunctionReturn(PETSC_SUCCESS); |
464 | } | ||
465 | |||
466 | 12 | static PetscErrorCode SVDPRIMMEGetMethod_PRIMME(SVD svd,SVDPRIMMEMethod *method) | |
467 | { | ||
468 | 12 | SVD_PRIMME *ops = (SVD_PRIMME*)svd->data; | |
469 | |||
470 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
12 | PetscFunctionBegin; |
471 | 12 | *method = (SVDPRIMMEMethod)ops->method; | |
472 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
473 | } | ||
474 | |||
475 | /*@ | ||
476 | SVDPRIMMEGetMethod - Gets the method for the PRIMME SVD solver. | ||
477 | |||
478 | Not Collective | ||
479 | |||
480 | Input Parameter: | ||
481 | . svd - the singular value solver context | ||
482 | |||
483 | Output Parameter: | ||
484 | . method - method that will be used by PRIMME | ||
485 | |||
486 | Level: advanced | ||
487 | |||
488 | .seealso: SVDPRIMMESetMethod(), SVDPRIMMEMethod | ||
489 | @*/ | ||
490 | 12 | PetscErrorCode SVDPRIMMEGetMethod(SVD svd,SVDPRIMMEMethod *method) | |
491 | { | ||
492 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
12 | PetscFunctionBegin; |
493 |
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.
|
12 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
494 |
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.
|
12 | PetscAssertPointer(method,2); |
495 |
9/16✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 6 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
|
12 | PetscUseMethod(svd,"SVDPRIMMEGetMethod_C",(SVD,SVDPRIMMEMethod*),(svd,method)); |
496 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
497 | } | ||
498 | |||
499 | 36 | SLEPC_EXTERN PetscErrorCode SVDCreate_PRIMME(SVD svd) | |
500 | { | ||
501 | 36 | SVD_PRIMME *primme; | |
502 | |||
503 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | PetscFunctionBegin; |
504 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscNew(&primme)); |
505 | 36 | svd->data = (void*)primme; | |
506 | |||
507 | 36 | primme_svds_initialize(&primme->primme); | |
508 | 36 | primme->bs = 0; | |
509 | 36 | primme->method = (primme_svds_preset_method)SVD_PRIMME_HYBRID; | |
510 | 36 | primme->svd = svd; | |
511 | |||
512 | 36 | svd->ops->solve = SVDSolve_PRIMME; | |
513 | 36 | svd->ops->setup = SVDSetUp_PRIMME; | |
514 | 36 | svd->ops->setfromoptions = SVDSetFromOptions_PRIMME; | |
515 | 36 | svd->ops->destroy = SVDDestroy_PRIMME; | |
516 | 36 | svd->ops->reset = SVDReset_PRIMME; | |
517 | 36 | svd->ops->view = SVDView_PRIMME; | |
518 | |||
519 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",SVDPRIMMESetBlockSize_PRIMME)); |
520 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",SVDPRIMMEGetBlockSize_PRIMME)); |
521 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",SVDPRIMMESetMethod_PRIMME)); |
522 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",SVDPRIMMEGetMethod_PRIMME)); |
523 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
524 | } | ||
525 |