Line data Source code
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 package
12 : */
13 :
14 : #include <slepc/private/epsimpl.h> /*I "slepceps.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
21 : #else
22 : #define PRIMME_DRIVER zprimme
23 : #endif
24 : #else
25 : #if defined(PETSC_USE_REAL_SINGLE)
26 : #define PRIMME_DRIVER sprimme
27 : #else
28 : #define PRIMME_DRIVER dprimme
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_params primme; /* param struct */
38 : PetscInt bs; /* block size */
39 : primme_preset_method method; /* primme method */
40 : Mat A,B; /* problem matrices */
41 : KSP ksp; /* linear solver and preconditioner */
42 : Vec x,y; /* auxiliary vectors */
43 : double target; /* a copy of eps's target */
44 : } EPS_PRIMME;
45 :
46 886 : static void par_GlobalSumReal(void *sendBuf,void *recvBuf,int *count,primme_params *primme,int *ierr)
47 : {
48 886 : if (sendBuf == recvBuf) {
49 1772 : *ierr = MPI_Allreduce(MPI_IN_PLACE,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));
50 : } else {
51 0 : *ierr = MPI_Allreduce(sendBuf,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));
52 : }
53 886 : }
54 :
55 : #if defined(SLEPC_HAVE_PRIMME3)
56 2332 : static void par_broadcastReal(void *buf,int *count,primme_params *primme,int *ierr)
57 : {
58 2332 : *ierr = MPI_Bcast(buf,*count,MPIU_REAL,0/*root*/,PetscObjectComm((PetscObject)primme->commInfo));
59 2332 : }
60 : #endif
61 :
62 : #if defined(SLEPC_HAVE_PRIMME2p2)
63 15752 : static void convTestFun(double *eval,void *evec,double *resNorm,int *isconv,primme_params *primme,int *err)
64 : {
65 15752 : PetscErrorCode ierr;
66 15752 : EPS eps = (EPS)primme->commInfo;
67 15752 : PetscScalar eigvr = eval?*eval:0.0;
68 15752 : PetscReal r = resNorm?*resNorm:HUGE_VAL,errest;
69 :
70 15752 : ierr = (*eps->converged)(eps,eigvr,0.0,r,&errest,eps->convergedctx);
71 15752 : if (ierr) *err = 1;
72 : else {
73 15752 : *isconv = (errest<=eps->tol?1:0);
74 15752 : *err = 0;
75 : }
76 15752 : }
77 :
78 14873 : static void monitorFun(void *basisEvals,int *basisSize,int *basisFlags,int *iblock,int *blockSize,void *basisNorms,int *numConverged,void *lockedEvals,int *numLocked,int *lockedFlags,void *lockedNorms,int *inner_its,void *LSRes,
79 : #if defined(SLEPC_HAVE_PRIMME3)
80 : const char *msg,double *time,
81 : #endif
82 : primme_event *event,struct primme_params *primme,int *err)
83 : {
84 14873 : PetscErrorCode ierr = PETSC_SUCCESS;
85 14873 : EPS eps = (EPS)primme->commInfo;
86 14873 : PetscInt i,k,nerrest;
87 :
88 14873 : switch (*event) {
89 1041 : case primme_event_outer_iteration:
90 : /* Update EPS */
91 1041 : eps->its = primme->stats.numOuterIterations;
92 1041 : eps->nconv = primme->initSize;
93 1041 : k=0;
94 2798 : if (lockedEvals && numLocked) for (i=0; i<*numLocked && k<eps->ncv; i++) eps->eigr[k++] = ((PetscReal*)lockedEvals)[i];
95 1041 : nerrest = k;
96 1041 : if (iblock && blockSize) {
97 2635 : for (i=0; i<*blockSize && k+iblock[i]<eps->ncv; i++) eps->errest[k+iblock[i]] = ((PetscReal*)basisNorms)[i];
98 1041 : nerrest = k+(*blockSize>0?1+iblock[*blockSize-1]:0);
99 : }
100 25127 : if (basisEvals && basisSize) for (i=0; i<*basisSize && k<eps->ncv; i++) eps->eigr[k++] = ((PetscReal*)basisEvals)[i];
101 : /* Show progress */
102 1041 : ierr = EPSMonitor(eps,eps->its,numConverged?*numConverged:0,eps->eigr,eps->eigi,eps->errest,nerrest);
103 1041 : break;
104 : #if defined(SLEPC_HAVE_PRIMME3)
105 0 : case primme_event_message:
106 : /* Print PRIMME information messages */
107 0 : ierr = PetscInfo(eps,"%s\n",msg);
108 0 : break;
109 : #endif
110 : default:
111 : break;
112 : }
113 14873 : *err = (ierr!=0)? 1: 0;
114 14873 : }
115 : #endif /* SLEPC_HAVE_PRIMME2p2 */
116 :
117 14287 : static void matrixMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
118 : {
119 14287 : PetscInt i;
120 14287 : EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
121 14287 : Vec x = ops->x,y = ops->y;
122 14287 : Mat A = ops->A;
123 :
124 14287 : PetscFunctionBegin;
125 30672 : for (i=0;i<*blockSize;i++) {
126 16385 : PetscCallAbort(PetscObjectComm((PetscObject)A),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i));
127 16385 : PetscCallAbort(PetscObjectComm((PetscObject)A),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i));
128 16385 : PetscCallAbort(PetscObjectComm((PetscObject)A),MatMult(A,x,y));
129 16385 : PetscCallAbort(PetscObjectComm((PetscObject)A),VecResetArray(x));
130 16385 : PetscCallAbort(PetscObjectComm((PetscObject)A),VecResetArray(y));
131 : }
132 14287 : PetscFunctionReturnVoid();
133 : }
134 :
135 : #if defined(SLEPC_HAVE_PRIMME3)
136 22653 : static void massMatrixMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
137 : {
138 22653 : PetscInt i;
139 22653 : EPS_PRIMME *ops = (EPS_PRIMME*)primme->massMatrix;
140 22653 : Vec x = ops->x,y = ops->y;
141 22653 : Mat B = ops->B;
142 :
143 22653 : PetscFunctionBegin;
144 45306 : for (i=0;i<*blockSize;i++) {
145 22653 : PetscCallAbort(PetscObjectComm((PetscObject)B),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i));
146 22653 : PetscCallAbort(PetscObjectComm((PetscObject)B),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i));
147 22653 : PetscCallAbort(PetscObjectComm((PetscObject)B),MatMult(B,x,y));
148 22653 : PetscCallAbort(PetscObjectComm((PetscObject)B),VecResetArray(x));
149 22653 : PetscCallAbort(PetscObjectComm((PetscObject)B),VecResetArray(y));
150 : }
151 22653 : PetscFunctionReturnVoid();
152 : }
153 : #endif
154 :
155 13436 : static void applyPreconditioner_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
156 : {
157 13436 : PetscInt i;
158 13436 : EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
159 13436 : Vec x = ops->x,y = ops->y;
160 :
161 13436 : PetscFunctionBegin;
162 28786 : for (i=0;i<*blockSize;i++) {
163 15350 : PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i));
164 15350 : PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i));
165 15350 : PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),KSPSolve(ops->ksp,x,y));
166 15350 : PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecResetArray(x));
167 15350 : PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecResetArray(y));
168 : }
169 13436 : PetscFunctionReturnVoid();
170 : }
171 :
172 10 : static PetscErrorCode EPSSetUp_PRIMME(EPS eps)
173 : {
174 10 : PetscMPIInt numProcs,procID;
175 10 : EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
176 10 : primme_params *primme = &ops->primme;
177 10 : PetscBool flg;
178 :
179 10 : PetscFunctionBegin;
180 10 : EPSCheckHermitianDefinite(eps);
181 10 : PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&numProcs));
182 10 : PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&procID));
183 :
184 : /* Check some constraints and set some default values */
185 10 : if (eps->max_it==PETSC_DEFAULT) eps->max_it = PETSC_MAX_INT;
186 10 : PetscCall(STGetMatrix(eps->st,0,&ops->A));
187 10 : if (eps->isgeneralized) {
188 : #if defined(SLEPC_HAVE_PRIMME3)
189 1 : PetscCall(STGetMatrix(eps->st,1,&ops->B));
190 : #else
191 : SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This version of PRIMME is not available for generalized problems");
192 : #endif
193 : }
194 10 : EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
195 10 : EPSCheckIgnored(eps,EPS_FEATURE_BALANCE);
196 10 : if (!eps->which) eps->which = EPS_LARGEST_REAL;
197 : #if !defined(SLEPC_HAVE_PRIMME2p2)
198 : if (eps->converged != EPSConvergedAbsolute) PetscCall(PetscInfo(eps,"Warning: using absolute convergence test\n"));
199 : #else
200 10 : EPSCheckIgnored(eps,EPS_FEATURE_CONVERGENCE);
201 : #endif
202 :
203 : /* Transfer SLEPc options to PRIMME options */
204 10 : primme_free(primme);
205 10 : primme_initialize(primme);
206 10 : primme->n = eps->n;
207 10 : primme->nLocal = eps->nloc;
208 10 : primme->numEvals = eps->nev;
209 10 : primme->matrix = ops;
210 10 : primme->matrixMatvec = matrixMatvec_PRIMME;
211 : #if defined(SLEPC_HAVE_PRIMME3)
212 10 : if (eps->isgeneralized) {
213 1 : primme->massMatrix = ops;
214 1 : primme->massMatrixMatvec = massMatrixMatvec_PRIMME;
215 : }
216 : #endif
217 10 : primme->commInfo = eps;
218 10 : primme->maxOuterIterations = eps->max_it;
219 : #if !defined(SLEPC_HAVE_PRIMME2p2)
220 : primme->eps = SlepcDefaultTol(eps->tol);
221 : #endif
222 10 : primme->numProcs = numProcs;
223 10 : primme->procID = procID;
224 10 : primme->printLevel = 1;
225 10 : primme->correctionParams.precondition = 1;
226 10 : primme->globalSumReal = par_GlobalSumReal;
227 : #if defined(SLEPC_HAVE_PRIMME3)
228 10 : primme->broadcastReal = par_broadcastReal;
229 : #endif
230 : #if defined(SLEPC_HAVE_PRIMME2p2)
231 10 : primme->convTestFun = convTestFun;
232 10 : primme->monitorFun = monitorFun;
233 : #endif
234 10 : if (ops->bs > 0) primme->maxBlockSize = ops->bs;
235 :
236 10 : switch (eps->which) {
237 3 : case EPS_LARGEST_REAL:
238 3 : primme->target = primme_largest;
239 3 : break;
240 2 : case EPS_SMALLEST_REAL:
241 2 : primme->target = primme_smallest;
242 2 : break;
243 1 : case EPS_LARGEST_MAGNITUDE:
244 1 : primme->target = primme_largest_abs;
245 1 : ops->target = 0.0;
246 1 : primme->numTargetShifts = 1;
247 1 : primme->targetShifts = &ops->target;
248 1 : break;
249 1 : case EPS_SMALLEST_MAGNITUDE:
250 1 : primme->target = primme_closest_abs;
251 1 : ops->target = 0.0;
252 1 : primme->numTargetShifts = 1;
253 1 : primme->targetShifts = &ops->target;
254 1 : break;
255 3 : case EPS_TARGET_MAGNITUDE:
256 : case EPS_TARGET_REAL:
257 3 : primme->target = primme_closest_abs;
258 3 : primme->numTargetShifts = 1;
259 3 : ops->target = PetscRealPart(eps->target);
260 3 : primme->targetShifts = &ops->target;
261 3 : break;
262 0 : default:
263 0 : SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'which' value not supported by PRIMME");
264 : }
265 :
266 10 : switch (eps->extraction) {
267 7 : case EPS_RITZ:
268 7 : primme->projectionParams.projection = primme_proj_RR;
269 7 : break;
270 1 : case EPS_HARMONIC:
271 1 : primme->projectionParams.projection = primme_proj_harmonic;
272 1 : break;
273 2 : case EPS_REFINED:
274 2 : primme->projectionParams.projection = primme_proj_refined;
275 2 : break;
276 0 : default:
277 0 : SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'extraction' value not supported by PRIMME");
278 : }
279 :
280 : /* If user sets mpd or ncv, maxBasisSize is modified */
281 10 : if (eps->mpd!=PETSC_DEFAULT) {
282 1 : primme->maxBasisSize = eps->mpd;
283 1 : if (eps->ncv!=PETSC_DEFAULT) PetscCall(PetscInfo(eps,"Warning: 'ncv' is ignored by PRIMME\n"));
284 9 : } else if (eps->ncv!=PETSC_DEFAULT) primme->maxBasisSize = eps->ncv;
285 :
286 10 : PetscCheck(primme_set_method(ops->method,primme)>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME method not valid");
287 :
288 10 : eps->mpd = primme->maxBasisSize;
289 10 : eps->ncv = (primme->locking?eps->nev:0)+primme->maxBasisSize;
290 10 : ops->bs = primme->maxBlockSize;
291 :
292 : /* Set workspace */
293 10 : PetscCall(EPSAllocateSolution(eps,0));
294 :
295 : /* Setup the preconditioner */
296 10 : if (primme->correctionParams.precondition) {
297 10 : PetscCall(STGetKSP(eps->st,&ops->ksp));
298 10 : PetscCall(PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&flg));
299 10 : if (!flg) PetscCall(PetscInfo(eps,"Warning: ignoring KSP, should use KSPPREONLY\n"));
300 10 : primme->preconditioner = NULL;
301 10 : primme->applyPreconditioner = applyPreconditioner_PRIMME;
302 : }
303 :
304 : /* Prepare auxiliary vectors */
305 10 : if (!ops->x) PetscCall(MatCreateVecsEmpty(ops->A,&ops->x,&ops->y));
306 10 : PetscFunctionReturn(PETSC_SUCCESS);
307 : }
308 :
309 10 : static PetscErrorCode EPSSolve_PRIMME(EPS eps)
310 : {
311 10 : EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
312 10 : PetscScalar *a;
313 10 : PetscInt i,ld,ierrprimme;
314 10 : PetscReal *evals,*rnorms;
315 :
316 10 : PetscFunctionBegin;
317 : /* Reset some parameters left from previous runs */
318 : #if defined(SLEPC_HAVE_PRIMME2p2)
319 10 : ops->primme.aNorm = 0.0;
320 : #else
321 : /* Force PRIMME to stop by absolute error */
322 : ops->primme.aNorm = 1.0;
323 : #endif
324 10 : ops->primme.initSize = eps->nini;
325 10 : ops->primme.iseed[0] = -1;
326 10 : ops->primme.iseed[1] = -1;
327 10 : ops->primme.iseed[2] = -1;
328 10 : ops->primme.iseed[3] = -1;
329 10 : PetscCall(BVGetLeadingDimension(eps->V,&ld));
330 10 : ops->primme.ldevecs = ld;
331 :
332 : /* Call PRIMME solver */
333 10 : PetscCall(BVGetArray(eps->V,&a));
334 10 : PetscCall(PetscMalloc2(eps->ncv,&evals,eps->ncv,&rnorms));
335 10 : ierrprimme = PRIMME_DRIVER(evals,a,rnorms,&ops->primme);
336 267 : for (i=0;i<eps->ncv;i++) eps->eigr[i] = evals[i];
337 257 : for (i=0;i<eps->ncv;i++) eps->errest[i] = rnorms[i];
338 10 : PetscCall(PetscFree2(evals,rnorms));
339 10 : PetscCall(BVRestoreArray(eps->V,&a));
340 :
341 10 : eps->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
342 10 : eps->reason = eps->nconv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
343 10 : eps->its = ops->primme.stats.numOuterIterations;
344 :
345 : /* Process PRIMME error code */
346 10 : switch (ierrprimme) {
347 : case 0: /* no error */
348 : break;
349 0 : case -1:
350 0 : SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": unexpected error",ierrprimme);
351 0 : case -2:
352 0 : SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": allocation error",ierrprimme);
353 : case -3: /* stop due to maximum number of iterations or matvecs */
354 : break;
355 0 : default:
356 0 : PetscCheck(ierrprimme<-39,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": configuration error; check PRIMME's manual",ierrprimme);
357 0 : PetscCheck(ierrprimme>=-39,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": runtime error; check PRIMME's manual",ierrprimme);
358 : }
359 10 : PetscFunctionReturn(PETSC_SUCCESS);
360 : }
361 :
362 9 : static PetscErrorCode EPSReset_PRIMME(EPS eps)
363 : {
364 9 : EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
365 :
366 9 : PetscFunctionBegin;
367 9 : primme_free(&ops->primme);
368 9 : PetscCall(VecDestroy(&ops->x));
369 9 : PetscCall(VecDestroy(&ops->y));
370 9 : PetscFunctionReturn(PETSC_SUCCESS);
371 : }
372 :
373 9 : static PetscErrorCode EPSDestroy_PRIMME(EPS eps)
374 : {
375 9 : PetscFunctionBegin;
376 9 : PetscCall(PetscFree(eps->data));
377 9 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL));
378 9 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL));
379 9 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL));
380 9 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL));
381 9 : PetscFunctionReturn(PETSC_SUCCESS);
382 : }
383 :
384 0 : static PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
385 : {
386 0 : PetscBool isascii;
387 0 : EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
388 0 : PetscMPIInt rank;
389 :
390 0 : PetscFunctionBegin;
391 0 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
392 0 : if (isascii) {
393 0 : PetscCall(PetscViewerASCIIPrintf(viewer," block size=%" PetscInt_FMT "\n",ctx->bs));
394 0 : PetscCall(PetscViewerASCIIPrintf(viewer," solver method: %s\n",EPSPRIMMEMethods[(EPSPRIMMEMethod)ctx->method]));
395 :
396 : /* Display PRIMME params */
397 0 : PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank));
398 0 : if (!rank) primme_display_params(ctx->primme);
399 : }
400 0 : PetscFunctionReturn(PETSC_SUCCESS);
401 : }
402 :
403 8 : static PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps,PetscOptionItems *PetscOptionsObject)
404 : {
405 8 : EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
406 8 : PetscInt bs;
407 8 : EPSPRIMMEMethod meth;
408 8 : PetscBool flg;
409 :
410 8 : PetscFunctionBegin;
411 8 : PetscOptionsHeadBegin(PetscOptionsObject,"EPS PRIMME Options");
412 :
413 8 : PetscCall(PetscOptionsInt("-eps_primme_blocksize","Maximum block size","EPSPRIMMESetBlockSize",ctx->bs,&bs,&flg));
414 8 : if (flg) PetscCall(EPSPRIMMESetBlockSize(eps,bs));
415 :
416 8 : PetscCall(PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg));
417 8 : if (flg) PetscCall(EPSPRIMMESetMethod(eps,meth));
418 :
419 8 : PetscOptionsHeadEnd();
420 8 : PetscFunctionReturn(PETSC_SUCCESS);
421 : }
422 :
423 5 : static PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
424 : {
425 5 : EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
426 :
427 5 : PetscFunctionBegin;
428 5 : if (bs == PETSC_DEFAULT) ops->bs = 0;
429 : else {
430 5 : PetscCheck(bs>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
431 5 : ops->bs = bs;
432 : }
433 5 : PetscFunctionReturn(PETSC_SUCCESS);
434 : }
435 :
436 : /*@
437 : EPSPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.
438 :
439 : Logically Collective
440 :
441 : Input Parameters:
442 : + eps - the eigenproblem solver context
443 : - bs - block size
444 :
445 : Options Database Key:
446 : . -eps_primme_blocksize - Sets the max allowed block size value
447 :
448 : Notes:
449 : If the block size is not set, the value established by primme_initialize
450 : is used.
451 :
452 : The user should set the block size based on the architecture specifics
453 : of the target computer, as well as any a priori knowledge of multiplicities.
454 : The code does NOT require bs > 1 to find multiple eigenvalues. For some
455 : methods, keeping bs = 1 yields the best overall performance.
456 :
457 : Level: advanced
458 :
459 : .seealso: EPSPRIMMEGetBlockSize()
460 : @*/
461 5 : PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
462 : {
463 5 : PetscFunctionBegin;
464 5 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
465 20 : PetscValidLogicalCollectiveInt(eps,bs,2);
466 5 : PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
467 5 : PetscFunctionReturn(PETSC_SUCCESS);
468 : }
469 :
470 4 : static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
471 : {
472 4 : EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
473 :
474 4 : PetscFunctionBegin;
475 4 : *bs = ops->bs;
476 4 : PetscFunctionReturn(PETSC_SUCCESS);
477 : }
478 :
479 : /*@
480 : EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
481 :
482 : Not Collective
483 :
484 : Input Parameter:
485 : . eps - the eigenproblem solver context
486 :
487 : Output Parameter:
488 : . bs - returned block size
489 :
490 : Level: advanced
491 :
492 : .seealso: EPSPRIMMESetBlockSize()
493 : @*/
494 4 : PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
495 : {
496 4 : PetscFunctionBegin;
497 4 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
498 4 : PetscAssertPointer(bs,2);
499 4 : PetscUseMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
500 4 : PetscFunctionReturn(PETSC_SUCCESS);
501 : }
502 :
503 5 : static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
504 : {
505 5 : EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
506 :
507 5 : PetscFunctionBegin;
508 5 : ops->method = (primme_preset_method)method;
509 5 : PetscFunctionReturn(PETSC_SUCCESS);
510 : }
511 :
512 : /*@
513 : EPSPRIMMESetMethod - Sets the method for the PRIMME library.
514 :
515 : Logically Collective
516 :
517 : Input Parameters:
518 : + eps - the eigenproblem solver context
519 : - method - method that will be used by PRIMME
520 :
521 : Options Database Key:
522 : . -eps_primme_method - Sets the method for the PRIMME library
523 :
524 : Note:
525 : If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.
526 :
527 : Level: advanced
528 :
529 : .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
530 : @*/
531 5 : PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
532 : {
533 5 : PetscFunctionBegin;
534 5 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
535 20 : PetscValidLogicalCollectiveEnum(eps,method,2);
536 5 : PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
537 5 : PetscFunctionReturn(PETSC_SUCCESS);
538 : }
539 :
540 4 : static PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
541 : {
542 4 : EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
543 :
544 4 : PetscFunctionBegin;
545 4 : *method = (EPSPRIMMEMethod)ops->method;
546 4 : PetscFunctionReturn(PETSC_SUCCESS);
547 : }
548 :
549 : /*@
550 : EPSPRIMMEGetMethod - Gets the method for the PRIMME library.
551 :
552 : Not Collective
553 :
554 : Input Parameter:
555 : . eps - the eigenproblem solver context
556 :
557 : Output Parameter:
558 : . method - method that will be used by PRIMME
559 :
560 : Level: advanced
561 :
562 : .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
563 : @*/
564 4 : PetscErrorCode EPSPRIMMEGetMethod(EPS eps,EPSPRIMMEMethod *method)
565 : {
566 4 : PetscFunctionBegin;
567 4 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
568 4 : PetscAssertPointer(method,2);
569 4 : PetscUseMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
570 4 : PetscFunctionReturn(PETSC_SUCCESS);
571 : }
572 :
573 9 : SLEPC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
574 : {
575 9 : EPS_PRIMME *primme;
576 :
577 9 : PetscFunctionBegin;
578 9 : PetscCall(PetscNew(&primme));
579 9 : eps->data = (void*)primme;
580 :
581 9 : primme_initialize(&primme->primme);
582 9 : primme->primme.globalSumReal = par_GlobalSumReal;
583 : #if defined(SLEPC_HAVE_PRIMME3)
584 9 : primme->primme.broadcastReal = par_broadcastReal;
585 : #endif
586 : #if defined(SLEPC_HAVE_PRIMME2p2)
587 9 : primme->primme.convTestFun = convTestFun;
588 9 : primme->primme.monitorFun = monitorFun;
589 : #endif
590 9 : primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;
591 :
592 9 : eps->categ = EPS_CATEGORY_PRECOND;
593 :
594 9 : eps->ops->solve = EPSSolve_PRIMME;
595 9 : eps->ops->setup = EPSSetUp_PRIMME;
596 9 : eps->ops->setupsort = EPSSetUpSort_Basic;
597 9 : eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
598 9 : eps->ops->destroy = EPSDestroy_PRIMME;
599 9 : eps->ops->reset = EPSReset_PRIMME;
600 9 : eps->ops->view = EPSView_PRIMME;
601 9 : eps->ops->backtransform = EPSBackTransform_Default;
602 9 : eps->ops->setdefaultst = EPSSetDefaultST_GMRES;
603 :
604 9 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME));
605 9 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME));
606 9 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME));
607 9 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME));
608 9 : PetscFunctionReturn(PETSC_SUCCESS);
609 : }
|