Actual source code: primme.c
slepc-3.16.1 2021-11-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: #if defined(PETSC_HAVE_MPI_IN_PLACE)
50: *MPI_Allreduce(MPI_IN_PLACE,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));
51: #else
52: #error Cannot use SLEPc-PRIMME interface if your MPI does not support MPI_IN_PLACE
53: #endif
54: } else {
55: *MPI_Allreduce(sendBuf,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));
56: }
57: }
59: #if defined(SLEPC_HAVE_PRIMME3)
60: static void par_broadcastReal(void *buf,int *count,primme_params *primme,int *ierr)
61: {
62: *MPI_Bcast(buf,*count,MPIU_REAL,0/*root*/,PetscObjectComm((PetscObject)primme->commInfo));
63: }
64: #endif
66: #if defined(SLEPC_HAVE_PRIMME2p2)
67: static void convTestFun(double *eval,void *evec,double *resNorm,int *isconv,primme_params *primme,int *err)
68: {
70: EPS eps = (EPS)primme->commInfo;
71: PetscScalar eigvr = eval?*eval:0.0;
72: PetscReal r = resNorm?*resNorm:HUGE_VAL,errest;
74: (*eps->converged)(eps,eigvr,0.0,r,&errest,eps->convergedctx);
75: if (ierr) *err = 1;
76: else {
77: *isconv = (errest<=eps->tol?1:0);
78: *err = 0;
79: }
80: }
82: 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,
83: #if defined(SLEPC_HAVE_PRIMME3)
84: const char *msg,double *time,
85: #endif
86: primme_event *event,struct primme_params *primme,int *err)
87: {
88: PetscErrorCode 0;
89: EPS eps = (EPS)primme->commInfo;
90: PetscInt i,k,nerrest;
92: switch (*event) {
93: case primme_event_outer_iteration:
94: /* Update EPS */
95: eps->its = primme->stats.numOuterIterations;
96: eps->nconv = primme->initSize;
97: k=0;
98: if (lockedEvals && numLocked) for (i=0; i<*numLocked && k<eps->ncv; i++) eps->eigr[k++] = ((PetscReal*)lockedEvals)[i];
99: nerrest = k;
100: if (iblock && blockSize) {
101: for (i=0; i<*blockSize && k+iblock[i]<eps->ncv; i++) eps->errest[k+iblock[i]] = ((PetscReal*)basisNorms)[i];
102: nerrest = k+(*blockSize>0?1+iblock[*blockSize-1]:0);
103: }
104: if (basisEvals && basisSize) for (i=0; i<*basisSize && k<eps->ncv; i++) eps->eigr[k++] = ((PetscReal*)basisEvals)[i];
105: /* Show progress */
106: EPSMonitor(eps,eps->its,numConverged?*numConverged:0,eps->eigr,eps->eigi,eps->errest,nerrest);
107: break;
108: #if defined(SLEPC_HAVE_PRIMME3)
109: case primme_event_message:
110: /* Print PRIMME information messages */
111: PetscInfo(eps,msg);
112: break;
113: #endif
114: default:
115: break;
116: }
117: *err = (ierr!=0)? 1: 0;
118: }
119: #endif /* SLEPC_HAVE_PRIMME2p2 */
121: static void matrixMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
122: {
123: PetscInt i;
124: EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
125: Vec x = ops->x,y = ops->y;
126: Mat A = ops->A;
129: for (i=0;i<*blockSize;i++) {
130: *VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
131: *VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
132: *MatMult(A,x,y);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
133: *VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
134: *VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
135: }
136: PetscFunctionReturnVoid();
137: }
139: #if defined(SLEPC_HAVE_PRIMME3)
140: static void massMatrixMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
141: {
142: PetscInt i;
143: EPS_PRIMME *ops = (EPS_PRIMME*)primme->massMatrix;
144: Vec x = ops->x,y = ops->y;
145: Mat B = ops->B;
148: for (i=0;i<*blockSize;i++) {
149: *VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);CHKERRABORT(PetscObjectComm((PetscObject)B),*ierr);
150: *VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);CHKERRABORT(PetscObjectComm((PetscObject)B),*ierr);
151: *MatMult(B,x,y);CHKERRABORT(PetscObjectComm((PetscObject)B),*ierr);
152: *VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)B),*ierr);
153: *VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)B),*ierr);
154: }
155: PetscFunctionReturnVoid();
156: }
157: #endif
159: static void applyPreconditioner_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
160: {
161: PetscInt i;
162: EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
163: Vec x = ops->x,y = ops->y;
166: for (i=0;i<*blockSize;i++) {
167: *VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
168: *VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
169: *KSPSolve(ops->ksp,x,y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
170: *VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
171: *VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
172: }
173: PetscFunctionReturnVoid();
174: }
176: PetscErrorCode EPSSetUp_PRIMME(EPS eps)
177: {
179: PetscMPIInt numProcs,procID;
180: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
181: primme_params *primme = &ops->primme;
182: PetscBool flg;
185: EPSCheckHermitianDefinite(eps);
186: MPI_Comm_size(PetscObjectComm((PetscObject)eps),&numProcs);
187: MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&procID);
189: /* Check some constraints and set some default values */
190: if (eps->max_it==PETSC_DEFAULT) eps->max_it = PETSC_MAX_INT;
191: STGetMatrix(eps->st,0,&ops->A);
192: if (eps->isgeneralized) {
193: #if defined(SLEPC_HAVE_PRIMME3)
194: STGetMatrix(eps->st,1,&ops->B);
195: #else
196: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This version of PRIMME is not available for generalized problems");
197: #endif
198: }
199: EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
200: EPSCheckIgnored(eps,EPS_FEATURE_BALANCE);
201: if (!eps->which) eps->which = EPS_LARGEST_REAL;
202: #if !defined(SLEPC_HAVE_PRIMME2p2)
203: if (eps->converged != EPSConvergedAbsolute) { PetscInfo(eps,"Warning: using absolute convergence test\n"); }
204: #else
205: EPSCheckIgnored(eps,EPS_FEATURE_CONVERGENCE);
206: #endif
208: /* Transfer SLEPc options to PRIMME options */
209: primme_free(primme);
210: primme_initialize(primme);
211: primme->n = eps->n;
212: primme->nLocal = eps->nloc;
213: primme->numEvals = eps->nev;
214: primme->matrix = ops;
215: primme->matrixMatvec = matrixMatvec_PRIMME;
216: #if defined(SLEPC_HAVE_PRIMME3)
217: if (eps->isgeneralized) {
218: primme->massMatrix = ops;
219: primme->massMatrixMatvec = massMatrixMatvec_PRIMME;
220: }
221: #endif
222: primme->commInfo = eps;
223: primme->maxOuterIterations = eps->max_it;
224: #if !defined(SLEPC_HAVE_PRIMME2p2)
225: primme->eps = SlepcDefaultTol(eps->tol);
226: #endif
227: primme->numProcs = numProcs;
228: primme->procID = procID;
229: primme->printLevel = 1;
230: primme->correctionParams.precondition = 1;
231: primme->globalSumReal = par_GlobalSumReal;
232: #if defined(SLEPC_HAVE_PRIMME3)
233: primme->broadcastReal = par_broadcastReal;
234: #endif
235: #if defined(SLEPC_HAVE_PRIMME2p2)
236: primme->convTestFun = convTestFun;
237: primme->monitorFun = monitorFun;
238: #endif
239: if (ops->bs > 0) primme->maxBlockSize = ops->bs;
241: switch (eps->which) {
242: case EPS_LARGEST_REAL:
243: primme->target = primme_largest;
244: break;
245: case EPS_SMALLEST_REAL:
246: primme->target = primme_smallest;
247: break;
248: case EPS_LARGEST_MAGNITUDE:
249: primme->target = primme_largest_abs;
250: ops->target = 0.0;
251: primme->numTargetShifts = 1;
252: primme->targetShifts = &ops->target;
253: break;
254: case EPS_SMALLEST_MAGNITUDE:
255: primme->target = primme_closest_abs;
256: ops->target = 0.0;
257: primme->numTargetShifts = 1;
258: primme->targetShifts = &ops->target;
259: break;
260: case EPS_TARGET_MAGNITUDE:
261: case EPS_TARGET_REAL:
262: primme->target = primme_closest_abs;
263: primme->numTargetShifts = 1;
264: ops->target = PetscRealPart(eps->target);
265: primme->targetShifts = &ops->target;
266: break;
267: default:
268: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'which' value not supported by PRIMME");
269: }
271: switch (eps->extraction) {
272: case EPS_RITZ:
273: primme->projectionParams.projection = primme_proj_RR;
274: break;
275: case EPS_HARMONIC:
276: primme->projectionParams.projection = primme_proj_harmonic;
277: break;
278: case EPS_REFINED:
279: primme->projectionParams.projection = primme_proj_refined;
280: break;
281: default:
282: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'extraction' value not supported by PRIMME");
283: }
285: /* If user sets mpd or ncv, maxBasisSize is modified */
286: if (eps->mpd!=PETSC_DEFAULT) {
287: primme->maxBasisSize = eps->mpd;
288: if (eps->ncv!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: 'ncv' is ignored by PRIMME\n"); }
289: } else if (eps->ncv!=PETSC_DEFAULT) primme->maxBasisSize = eps->ncv;
291: if (primme_set_method(ops->method,primme) < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME method not valid");
293: eps->mpd = primme->maxBasisSize;
294: eps->ncv = (primme->locking?eps->nev:0)+primme->maxBasisSize;
295: ops->bs = primme->maxBlockSize;
297: /* Set workspace */
298: EPSAllocateSolution(eps,0);
300: /* Setup the preconditioner */
301: if (primme->correctionParams.precondition) {
302: STGetKSP(eps->st,&ops->ksp);
303: PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&flg);
304: if (!flg) { PetscInfo(eps,"Warning: ignoring KSP, should use KSPPREONLY\n"); }
305: primme->preconditioner = NULL;
306: primme->applyPreconditioner = applyPreconditioner_PRIMME;
307: }
309: /* Prepare auxiliary vectors */
310: if (!ops->x) {
311: MatCreateVecsEmpty(ops->A,&ops->x,&ops->y);
312: PetscLogObjectParent((PetscObject)eps,(PetscObject)ops->x);
313: PetscLogObjectParent((PetscObject)eps,(PetscObject)ops->y);
314: }
315: return(0);
316: }
318: PetscErrorCode EPSSolve_PRIMME(EPS eps)
319: {
321: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
322: PetscScalar *a;
323: PetscInt i,ierrprimme;
324: PetscReal *evals,*rnorms;
327: /* Reset some parameters left from previous runs */
328: #if defined(SLEPC_HAVE_PRIMME2p2)
329: ops->primme.aNorm = 0.0;
330: #else
331: /* Force PRIMME to stop by absolute error */
332: ops->primme.aNorm = 1.0;
333: #endif
334: ops->primme.initSize = eps->nini;
335: ops->primme.iseed[0] = -1;
336: ops->primme.iseed[1] = -1;
337: ops->primme.iseed[2] = -1;
338: ops->primme.iseed[3] = -1;
340: /* Call PRIMME solver */
341: BVGetArray(eps->V,&a);
342: PetscMalloc2(eps->ncv,&evals,eps->ncv,&rnorms);
343: ierrprimme = PRIMME_DRIVER(evals,a,rnorms,&ops->primme);
344: for (i=0;i<eps->ncv;i++) eps->eigr[i] = evals[i];
345: for (i=0;i<eps->ncv;i++) eps->errest[i] = rnorms[i];
346: PetscFree2(evals,rnorms);
347: BVRestoreArray(eps->V,&a);
349: eps->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
350: eps->reason = eps->nconv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
351: eps->its = ops->primme.stats.numOuterIterations;
353: /* Process PRIMME error code */
354: if (ierrprimme == 0) {
355: /* no error */
356: } else if (ierrprimme == -1) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: unexpected error",ierrprimme);
357: else if (ierrprimme == -2) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: allocation error",ierrprimme);
358: else if (ierrprimme == -3) {
359: /* stop by maximum number of iteration or matvecs */
360: } else if (ierrprimme >= -39) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: configuration error; check PRIMME's manual",ierrprimme);
361: else SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: runtime error; check PRIMME's manual",ierrprimme);
362: return(0);
363: }
365: PetscErrorCode EPSReset_PRIMME(EPS eps)
366: {
368: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
371: primme_free(&ops->primme);
372: VecDestroy(&ops->x);
373: VecDestroy(&ops->y);
374: return(0);
375: }
377: PetscErrorCode EPSDestroy_PRIMME(EPS eps)
378: {
382: PetscFree(eps->data);
383: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL);
384: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL);
385: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL);
386: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL);
387: return(0);
388: }
390: PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
391: {
393: PetscBool isascii;
394: EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
395: PetscMPIInt rank;
398: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
399: if (isascii) {
400: PetscViewerASCIIPrintf(viewer," block size=%D\n",ctx->bs);
401: PetscViewerASCIIPrintf(viewer," solver method: %s\n",EPSPRIMMEMethods[(EPSPRIMMEMethod)ctx->method]);
403: /* Display PRIMME params */
404: MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank);
405: if (!rank) primme_display_params(ctx->primme);
406: }
407: return(0);
408: }
410: PetscErrorCode EPSSetFromOptions_PRIMME(PetscOptionItems *PetscOptionsObject,EPS eps)
411: {
412: PetscErrorCode ierr;
413: EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
414: PetscInt bs;
415: EPSPRIMMEMethod meth;
416: PetscBool flg;
419: PetscOptionsHead(PetscOptionsObject,"EPS PRIMME Options");
421: PetscOptionsInt("-eps_primme_blocksize","Maximum block size","EPSPRIMMESetBlockSize",ctx->bs,&bs,&flg);
422: if (flg) { EPSPRIMMESetBlockSize(eps,bs); }
424: PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
425: if (flg) { EPSPRIMMESetMethod(eps,meth); }
427: PetscOptionsTail();
428: return(0);
429: }
431: static PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
432: {
433: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
436: if (bs == PETSC_DEFAULT) ops->bs = 0;
437: else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
438: else ops->bs = bs;
439: return(0);
440: }
442: /*@
443: EPSPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.
445: Logically Collective on eps
447: Input Parameters:
448: + eps - the eigenproblem solver context
449: - bs - block size
451: Options Database Key:
452: . -eps_primme_blocksize - Sets the max allowed block size value
454: Notes:
455: If the block size is not set, the value established by primme_initialize
456: is used.
458: The user should set the block size based on the architecture specifics
459: of the target computer, as well as any a priori knowledge of multiplicities.
460: The code does NOT require bs > 1 to find multiple eigenvalues. For some
461: methods, keeping bs = 1 yields the best overall performance.
463: Level: advanced
465: .seealso: EPSPRIMMEGetBlockSize()
466: @*/
467: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
468: {
474: PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
475: return(0);
476: }
478: static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
479: {
480: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
483: *bs = ops->bs;
484: return(0);
485: }
487: /*@
488: EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
490: Not Collective
492: Input Parameter:
493: . eps - the eigenproblem solver context
495: Output Parameter:
496: . bs - returned block size
498: Level: advanced
500: .seealso: EPSPRIMMESetBlockSize()
501: @*/
502: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
503: {
509: PetscUseMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
510: return(0);
511: }
513: static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
514: {
515: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
518: ops->method = (primme_preset_method)method;
519: return(0);
520: }
522: /*@
523: EPSPRIMMESetMethod - Sets the method for the PRIMME library.
525: Logically Collective on eps
527: Input Parameters:
528: + eps - the eigenproblem solver context
529: - method - method that will be used by PRIMME
531: Options Database Key:
532: . -eps_primme_method - Sets the method for the PRIMME library
534: Note:
535: If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.
537: Level: advanced
539: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
540: @*/
541: PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
542: {
548: PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
549: return(0);
550: }
552: static PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
553: {
554: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
557: *method = (EPSPRIMMEMethod)ops->method;
558: return(0);
559: }
561: /*@
562: EPSPRIMMEGetMethod - Gets the method for the PRIMME library.
564: Not Collective
566: Input Parameter:
567: . eps - the eigenproblem solver context
569: Output Parameter:
570: . method - method that will be used by PRIMME
572: Level: advanced
574: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
575: @*/
576: PetscErrorCode EPSPRIMMEGetMethod(EPS eps,EPSPRIMMEMethod *method)
577: {
583: PetscUseMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
584: return(0);
585: }
587: SLEPC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
588: {
590: EPS_PRIMME *primme;
593: PetscNewLog(eps,&primme);
594: eps->data = (void*)primme;
596: primme_initialize(&primme->primme);
597: primme->primme.globalSumReal = par_GlobalSumReal;
598: #if defined(SLEPC_HAVE_PRIMME3)
599: primme->primme.broadcastReal = par_broadcastReal;
600: #endif
601: #if defined(SLEPC_HAVE_PRIMME2p2)
602: primme->primme.convTestFun = convTestFun;
603: primme->primme.monitorFun = monitorFun;
604: #endif
605: primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;
607: eps->categ = EPS_CATEGORY_PRECOND;
609: eps->ops->solve = EPSSolve_PRIMME;
610: eps->ops->setup = EPSSetUp_PRIMME;
611: eps->ops->setupsort = EPSSetUpSort_Basic;
612: eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
613: eps->ops->destroy = EPSDestroy_PRIMME;
614: eps->ops->reset = EPSReset_PRIMME;
615: eps->ops->view = EPSView_PRIMME;
616: eps->ops->backtransform = EPSBackTransform_Default;
617: eps->ops->setdefaultst = EPSSetDefaultST_GMRES;
619: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME);
620: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME);
621: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME);
622: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME);
623: return(0);
624: }