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