Actual source code: primme.c
slepc-3.22.2 2024-12-02
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: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&numProcs));
183: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&procID));
185: /* Check some constraints and set some default values */
186: if (eps->max_it==PETSC_DETERMINE) eps->max_it = PETSC_INT_MAX;
187: PetscCall(STGetMatrix(eps->st,0,&ops->A));
188: if (eps->isgeneralized) {
189: #if defined(SLEPC_HAVE_PRIMME3)
190: PetscCall(STGetMatrix(eps->st,1,&ops->B));
191: #else
192: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This version of PRIMME is not available for generalized problems");
193: #endif
194: }
195: EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
196: EPSCheckIgnored(eps,EPS_FEATURE_BALANCE);
197: if (!eps->which) eps->which = EPS_LARGEST_REAL;
198: #if !defined(SLEPC_HAVE_PRIMME2p2)
199: if (eps->converged != EPSConvergedAbsolute) PetscCall(PetscInfo(eps,"Warning: using absolute convergence test\n"));
200: #else
201: EPSCheckIgnored(eps,EPS_FEATURE_CONVERGENCE);
202: #endif
204: /* Transfer SLEPc options to PRIMME options */
205: primme_free(primme);
206: primme_initialize(primme);
207: primme->n = (PRIMME_INT)eps->n;
208: primme->nLocal = (PRIMME_INT)eps->nloc;
209: primme->numEvals = (int)eps->nev;
210: primme->matrix = ops;
211: primme->matrixMatvec = matrixMatvec_PRIMME;
212: #if defined(SLEPC_HAVE_PRIMME3)
213: if (eps->isgeneralized) {
214: primme->massMatrix = ops;
215: primme->massMatrixMatvec = massMatrixMatvec_PRIMME;
216: }
217: #endif
218: primme->commInfo = eps;
219: primme->maxOuterIterations = (PRIMME_INT)eps->max_it;
220: #if !defined(SLEPC_HAVE_PRIMME2p2)
221: primme->eps = SlepcDefaultTol(eps->tol);
222: #endif
223: primme->numProcs = numProcs;
224: primme->procID = procID;
225: primme->printLevel = 1;
226: primme->correctionParams.precondition = 1;
227: primme->globalSumReal = par_GlobalSumReal;
228: #if defined(SLEPC_HAVE_PRIMME3)
229: primme->broadcastReal = par_broadcastReal;
230: #endif
231: #if defined(SLEPC_HAVE_PRIMME2p2)
232: primme->convTestFun = convTestFun;
233: primme->monitorFun = monitorFun;
234: #endif
235: if (ops->bs > 0) primme->maxBlockSize = (int)ops->bs;
237: switch (eps->which) {
238: case EPS_LARGEST_REAL:
239: primme->target = primme_largest;
240: break;
241: case EPS_SMALLEST_REAL:
242: primme->target = primme_smallest;
243: break;
244: case EPS_LARGEST_MAGNITUDE:
245: primme->target = primme_largest_abs;
246: ops->target = 0.0;
247: primme->numTargetShifts = 1;
248: primme->targetShifts = &ops->target;
249: break;
250: case EPS_SMALLEST_MAGNITUDE:
251: primme->target = primme_closest_abs;
252: ops->target = 0.0;
253: primme->numTargetShifts = 1;
254: primme->targetShifts = &ops->target;
255: break;
256: case EPS_TARGET_MAGNITUDE:
257: case EPS_TARGET_REAL:
258: primme->target = primme_closest_abs;
259: primme->numTargetShifts = 1;
260: ops->target = PetscRealPart(eps->target);
261: primme->targetShifts = &ops->target;
262: break;
263: default:
264: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'which' value not supported by PRIMME");
265: }
267: switch (eps->extraction) {
268: case EPS_RITZ:
269: primme->projectionParams.projection = primme_proj_RR;
270: break;
271: case EPS_HARMONIC:
272: primme->projectionParams.projection = primme_proj_harmonic;
273: break;
274: case EPS_REFINED:
275: primme->projectionParams.projection = primme_proj_refined;
276: break;
277: default:
278: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'extraction' value not supported by PRIMME");
279: }
281: /* If user sets mpd or ncv, maxBasisSize is modified */
282: if (eps->mpd!=PETSC_DETERMINE) {
283: primme->maxBasisSize = (int)eps->mpd;
284: if (eps->ncv!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: 'ncv' is ignored by PRIMME\n"));
285: } else if (eps->ncv!=PETSC_DETERMINE) primme->maxBasisSize = (int)eps->ncv;
287: PetscCheck(primme_set_method(ops->method,primme)>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME method not valid");
289: eps->mpd = (PetscInt)primme->maxBasisSize;
290: eps->ncv = (PetscInt)(primme->locking?eps->nev:0)+primme->maxBasisSize;
291: ops->bs = (PetscInt)primme->maxBlockSize;
293: /* Set workspace */
294: PetscCall(EPSAllocateSolution(eps,0));
296: /* Setup the preconditioner */
297: if (primme->correctionParams.precondition) {
298: PetscCall(STGetKSP(eps->st,&ops->ksp));
299: PetscCall(PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&flg));
300: if (!flg) PetscCall(PetscInfo(eps,"Warning: ignoring KSP, should use KSPPREONLY\n"));
301: primme->preconditioner = NULL;
302: primme->applyPreconditioner = applyPreconditioner_PRIMME;
303: }
305: /* Prepare auxiliary vectors */
306: if (!ops->x) PetscCall(MatCreateVecsEmpty(ops->A,&ops->x,&ops->y));
307: PetscFunctionReturn(PETSC_SUCCESS);
308: }
310: static PetscErrorCode EPSSolve_PRIMME(EPS eps)
311: {
312: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
313: PetscScalar *a;
314: PetscInt i,ld,ierrprimme;
315: PetscReal *evals,*rnorms;
317: PetscFunctionBegin;
318: /* Reset some parameters left from previous runs */
319: #if defined(SLEPC_HAVE_PRIMME2p2)
320: ops->primme.aNorm = 0.0;
321: #else
322: /* Force PRIMME to stop by absolute error */
323: ops->primme.aNorm = 1.0;
324: #endif
325: ops->primme.initSize = (int)eps->nini;
326: ops->primme.iseed[0] = -1;
327: ops->primme.iseed[1] = -1;
328: ops->primme.iseed[2] = -1;
329: ops->primme.iseed[3] = -1;
330: PetscCall(BVGetLeadingDimension(eps->V,&ld));
331: ops->primme.ldevecs = (PRIMME_INT)ld;
333: /* Call PRIMME solver */
334: PetscCall(BVGetArray(eps->V,&a));
335: PetscCall(PetscMalloc2(eps->ncv,&evals,eps->ncv,&rnorms));
336: ierrprimme = PRIMME_DRIVER(evals,a,rnorms,&ops->primme);
337: for (i=0;i<eps->ncv;i++) eps->eigr[i] = evals[i];
338: for (i=0;i<eps->ncv;i++) eps->errest[i] = rnorms[i];
339: PetscCall(PetscFree2(evals,rnorms));
340: PetscCall(BVRestoreArray(eps->V,&a));
342: eps->nconv = ops->primme.initSize >= 0 ? (PetscInt)ops->primme.initSize : 0;
343: eps->reason = eps->nconv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
344: PetscCall(PetscIntCast(ops->primme.stats.numOuterIterations,&eps->its));
346: /* Process PRIMME error code */
347: switch (ierrprimme) {
348: case 0: /* no error */
349: break;
350: case -1:
351: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": unexpected error",ierrprimme);
352: case -2:
353: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": allocation error",ierrprimme);
354: case -3: /* stop due to maximum number of iterations or matvecs */
355: break;
356: default:
357: PetscCheck(ierrprimme<-39,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": configuration error; check PRIMME's manual",ierrprimme);
358: PetscCheck(ierrprimme>=-39,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": runtime error; check PRIMME's manual",ierrprimme);
359: }
360: PetscFunctionReturn(PETSC_SUCCESS);
361: }
363: static PetscErrorCode EPSReset_PRIMME(EPS eps)
364: {
365: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
367: PetscFunctionBegin;
368: primme_free(&ops->primme);
369: PetscCall(VecDestroy(&ops->x));
370: PetscCall(VecDestroy(&ops->y));
371: PetscFunctionReturn(PETSC_SUCCESS);
372: }
374: static PetscErrorCode EPSDestroy_PRIMME(EPS eps)
375: {
376: PetscFunctionBegin;
377: PetscCall(PetscFree(eps->data));
378: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL));
379: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL));
380: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL));
381: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL));
382: PetscFunctionReturn(PETSC_SUCCESS);
383: }
385: static PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
386: {
387: PetscBool isascii;
388: EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
389: PetscMPIInt rank;
391: PetscFunctionBegin;
392: PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
393: if (isascii) {
394: PetscCall(PetscViewerASCIIPrintf(viewer," block size=%" PetscInt_FMT "\n",ctx->bs));
395: PetscCall(PetscViewerASCIIPrintf(viewer," solver method: %s\n",EPSPRIMMEMethods[(EPSPRIMMEMethod)ctx->method]));
397: /* Display PRIMME params */
398: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank));
399: if (!rank) primme_display_params(ctx->primme);
400: }
401: PetscFunctionReturn(PETSC_SUCCESS);
402: }
404: static PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps,PetscOptionItems *PetscOptionsObject)
405: {
406: EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
407: PetscInt bs;
408: EPSPRIMMEMethod meth;
409: PetscBool flg;
411: PetscFunctionBegin;
412: PetscOptionsHeadBegin(PetscOptionsObject,"EPS PRIMME Options");
414: PetscCall(PetscOptionsInt("-eps_primme_blocksize","Maximum block size","EPSPRIMMESetBlockSize",ctx->bs,&bs,&flg));
415: if (flg) PetscCall(EPSPRIMMESetBlockSize(eps,bs));
417: PetscCall(PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg));
418: if (flg) PetscCall(EPSPRIMMESetMethod(eps,meth));
420: PetscOptionsHeadEnd();
421: PetscFunctionReturn(PETSC_SUCCESS);
422: }
424: static PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
425: {
426: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
428: PetscFunctionBegin;
429: if (bs == PETSC_DEFAULT || bs == PETSC_DECIDE) ops->bs = 0;
430: else {
431: PetscCheck(bs>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
432: ops->bs = bs;
433: }
434: PetscFunctionReturn(PETSC_SUCCESS);
435: }
437: /*@
438: EPSPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.
440: Logically Collective
442: Input Parameters:
443: + eps - the eigenproblem solver context
444: - bs - block size
446: Options Database Key:
447: . -eps_primme_blocksize - Sets the max allowed block size value
449: Notes:
450: If the block size is not set, the value established by primme_initialize
451: is used.
453: The user should set the block size based on the architecture specifics
454: of the target computer, as well as any a priori knowledge of multiplicities.
455: The code does NOT require bs > 1 to find multiple eigenvalues. For some
456: methods, keeping bs = 1 yields the best overall performance.
458: Level: advanced
460: .seealso: EPSPRIMMEGetBlockSize()
461: @*/
462: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
463: {
464: PetscFunctionBegin;
467: PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
468: PetscFunctionReturn(PETSC_SUCCESS);
469: }
471: static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
472: {
473: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
475: PetscFunctionBegin;
476: *bs = ops->bs;
477: PetscFunctionReturn(PETSC_SUCCESS);
478: }
480: /*@
481: EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
483: Not Collective
485: Input Parameter:
486: . eps - the eigenproblem solver context
488: Output Parameter:
489: . bs - returned block size
491: Level: advanced
493: .seealso: EPSPRIMMESetBlockSize()
494: @*/
495: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
496: {
497: PetscFunctionBegin;
499: PetscAssertPointer(bs,2);
500: PetscUseMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
501: PetscFunctionReturn(PETSC_SUCCESS);
502: }
504: static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
505: {
506: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
508: PetscFunctionBegin;
509: ops->method = (primme_preset_method)method;
510: PetscFunctionReturn(PETSC_SUCCESS);
511: }
513: /*@
514: EPSPRIMMESetMethod - Sets the method for the PRIMME library.
516: Logically Collective
518: Input Parameters:
519: + eps - the eigenproblem solver context
520: - method - method that will be used by PRIMME
522: Options Database Key:
523: . -eps_primme_method - Sets the method for the PRIMME library
525: Note:
526: If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.
528: Level: advanced
530: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
531: @*/
532: PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
533: {
534: PetscFunctionBegin;
537: PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
538: PetscFunctionReturn(PETSC_SUCCESS);
539: }
541: static PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
542: {
543: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
545: PetscFunctionBegin;
546: *method = (EPSPRIMMEMethod)ops->method;
547: PetscFunctionReturn(PETSC_SUCCESS);
548: }
550: /*@
551: EPSPRIMMEGetMethod - Gets the method for the PRIMME library.
553: Not Collective
555: Input Parameter:
556: . eps - the eigenproblem solver context
558: Output Parameter:
559: . method - method that will be used by PRIMME
561: Level: advanced
563: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
564: @*/
565: PetscErrorCode EPSPRIMMEGetMethod(EPS eps,EPSPRIMMEMethod *method)
566: {
567: PetscFunctionBegin;
569: PetscAssertPointer(method,2);
570: PetscUseMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
571: PetscFunctionReturn(PETSC_SUCCESS);
572: }
574: SLEPC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
575: {
576: EPS_PRIMME *primme;
578: PetscFunctionBegin;
579: PetscCall(PetscNew(&primme));
580: eps->data = (void*)primme;
582: primme_initialize(&primme->primme);
583: primme->primme.globalSumReal = par_GlobalSumReal;
584: #if defined(SLEPC_HAVE_PRIMME3)
585: primme->primme.broadcastReal = par_broadcastReal;
586: #endif
587: #if defined(SLEPC_HAVE_PRIMME2p2)
588: primme->primme.convTestFun = convTestFun;
589: primme->primme.monitorFun = monitorFun;
590: #endif
591: primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;
593: eps->categ = EPS_CATEGORY_PRECOND;
595: eps->ops->solve = EPSSolve_PRIMME;
596: eps->ops->setup = EPSSetUp_PRIMME;
597: eps->ops->setupsort = EPSSetUpSort_Basic;
598: eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
599: eps->ops->destroy = EPSDestroy_PRIMME;
600: eps->ops->reset = EPSReset_PRIMME;
601: eps->ops->view = EPSView_PRIMME;
602: eps->ops->backtransform = EPSBackTransform_Default;
603: eps->ops->setdefaultst = EPSSetDefaultST_GMRES;
605: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME));
606: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME));
607: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME));
608: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME));
609: PetscFunctionReturn(PETSC_SUCCESS);
610: }