Actual source code: primme.c
slepc-3.14.0 2020-09-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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 struc */
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: *MPIU_Allreduce(MPI_IN_PLACE,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));CHKERRABORT(PetscObjectComm((PetscObject)primme->commInfo),*ierr);
50: } else {
51: *MPIU_Allreduce(sendBuf,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));CHKERRABORT(PetscObjectComm((PetscObject)primme->commInfo),*ierr);
52: }
53: }
55: #if defined(SLEPC_HAVE_PRIMME3)
56: static void par_broadcastReal(void *buf,int *count,primme_params *primme,int *ierr)
57: {
58: *MPI_Bcast(buf,*count,MPIU_REAL,0/*root*/,PetscObjectComm((PetscObject)primme->commInfo));CHKERRABORT(PetscObjectComm((PetscObject)primme->commInfo),*ierr);
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: {
66: EPS eps = (EPS)primme->commInfo;
67: PetscScalar eigvr = eval?*eval:0.0;
68: PetscReal r = resNorm?*resNorm:HUGE_VAL,errest;
70: *err = 1;
71: (*eps->converged)(eps,eigvr,0.0,r,&errest,eps->convergedctx);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
72: *isconv = (errest<=eps->tol?1:0);
73: *err = 0;
74: }
76: 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,
77: #if defined(SLEPC_HAVE_PRIMME3)
78: const char *msg,double *time,
79: #endif
80: primme_event *event,struct primme_params *primme,int *err)
81: {
83: EPS eps = (EPS)primme->commInfo;
84: PetscInt i,k,nerrest;
86: *err = 1;
87: switch (*event) {
88: case primme_event_outer_iteration:
89: /* Update EPS */
90: eps->its = primme->stats.numOuterIterations;
91: eps->nconv = primme->initSize;
92: k=0;
93: if (lockedEvals && numLocked) for (i=0; i<*numLocked && k<eps->ncv; i++) eps->eigr[k++] = ((PetscReal*)lockedEvals)[i];
94: nerrest = k;
95: if (iblock && blockSize) {
96: for (i=0; i<*blockSize && k+iblock[i]<eps->ncv; i++) eps->errest[k+iblock[i]] = ((PetscReal*)basisNorms)[i];
97: nerrest = k+(*blockSize>0?1+iblock[*blockSize-1]:0);
98: }
99: if (basisEvals && basisSize) for (i=0; i<*basisSize && k<eps->ncv; i++) eps->eigr[k++] = ((PetscReal*)basisEvals)[i];
100: /* Show progress */
101: EPSMonitor(eps,eps->its,numConverged?*numConverged:0,eps->eigr,eps->eigi,eps->errest,nerrest);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
102: break;
103: #if defined(SLEPC_HAVE_PRIMME3)
104: case primme_event_message:
105: /* Print PRIMME information messages */
106: PetscInfo(eps,msg);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
107: break;
108: #endif
109: default:
110: break;
111: }
112: *err = 0;
113: }
114: #endif /* SLEPC_HAVE_PRIMME2p2 */
116: static void matrixMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
117: {
118: PetscInt i;
119: EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
120: Vec x = ops->x,y = ops->y;
121: Mat A = ops->A;
124: for (i=0;i<*blockSize;i++) {
125: *VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
126: *VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
127: *MatMult(A,x,y);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
128: *VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
129: *VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
130: }
131: PetscFunctionReturnVoid();
132: }
134: #if defined(SLEPC_HAVE_PRIMME3)
135: static void massMatrixMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
136: {
137: PetscInt i;
138: EPS_PRIMME *ops = (EPS_PRIMME*)primme->massMatrix;
139: Vec x = ops->x,y = ops->y;
140: Mat B = ops->B;
143: for (i=0;i<*blockSize;i++) {
144: *VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);CHKERRABORT(PetscObjectComm((PetscObject)B),*ierr);
145: *VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);CHKERRABORT(PetscObjectComm((PetscObject)B),*ierr);
146: *MatMult(B,x,y);CHKERRABORT(PetscObjectComm((PetscObject)B),*ierr);
147: *VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)B),*ierr);
148: *VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)B),*ierr);
149: }
150: PetscFunctionReturnVoid();
151: }
152: #endif
154: static void applyPreconditioner_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
155: {
156: PetscInt i;
157: EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
158: Vec x = ops->x,y = ops->y;
161: for (i=0;i<*blockSize;i++) {
162: *VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
163: *VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
164: *KSPSolve(ops->ksp,x,y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
165: *VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
166: *VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
167: }
168: PetscFunctionReturnVoid();
169: }
171: PetscErrorCode EPSSetUp_PRIMME(EPS eps)
172: {
174: PetscMPIInt numProcs,procID;
175: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
176: primme_params *primme = &ops->primme;
177: PetscBool flg;
180: EPSCheckHermitianDefinite(eps);
181: MPI_Comm_size(PetscObjectComm((PetscObject)eps),&numProcs);
182: MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&procID);
184: /* Check some constraints and set some default values */
185: if (eps->max_it==PETSC_DEFAULT) eps->max_it = PETSC_MAX_INT;
186: STGetMatrix(eps->st,0,&ops->A);
187: if (eps->isgeneralized) {
188: #if defined(SLEPC_HAVE_PRIMME3)
189: 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: EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
195: EPSCheckIgnored(eps,EPS_FEATURE_BALANCE);
196: if (!eps->which) eps->which = EPS_LARGEST_REAL;
197: #if !defined(SLEPC_HAVE_PRIMME2p2)
198: if (eps->converged != EPSConvergedAbsolute) { PetscInfo(eps,"Warning: using absolute convergence test\n"); }
199: #else
200: EPSCheckIgnored(eps,EPS_FEATURE_CONVERGENCE);
201: #endif
203: /* Transfer SLEPc options to PRIMME options */
204: primme_free(primme);
205: primme_initialize(primme);
206: primme->n = eps->n;
207: primme->nLocal = eps->nloc;
208: primme->numEvals = eps->nev;
209: primme->matrix = ops;
210: primme->matrixMatvec = matrixMatvec_PRIMME;
211: #if defined(SLEPC_HAVE_PRIMME3)
212: if (eps->isgeneralized) {
213: primme->massMatrix = ops;
214: primme->massMatrixMatvec = massMatrixMatvec_PRIMME;
215: }
216: #endif
217: primme->commInfo = eps;
218: primme->maxOuterIterations = eps->max_it;
219: #if !defined(SLEPC_HAVE_PRIMME2p2)
220: primme->eps = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
221: #endif
222: primme->numProcs = numProcs;
223: primme->procID = procID;
224: primme->printLevel = 1;
225: primme->correctionParams.precondition = 1;
226: primme->globalSumReal = par_GlobalSumReal;
227: #if defined(SLEPC_HAVE_PRIMME3)
228: primme->broadcastReal = par_broadcastReal;
229: #endif
230: #if defined(SLEPC_HAVE_PRIMME2p2)
231: primme->convTestFun = convTestFun;
232: primme->monitorFun = monitorFun;
233: #endif
234: if (ops->bs > 0) primme->maxBlockSize = ops->bs;
236: switch (eps->which) {
237: case EPS_LARGEST_REAL:
238: primme->target = primme_largest;
239: break;
240: case EPS_SMALLEST_REAL:
241: primme->target = primme_smallest;
242: break;
243: case EPS_LARGEST_MAGNITUDE:
244: primme->target = primme_largest_abs;
245: ops->target = 0.0;
246: primme->numTargetShifts = 1;
247: primme->targetShifts = &ops->target;
248: break;
249: case EPS_SMALLEST_MAGNITUDE:
250: primme->target = primme_closest_abs;
251: ops->target = 0.0;
252: primme->numTargetShifts = 1;
253: primme->targetShifts = &ops->target;
254: break;
255: case EPS_TARGET_MAGNITUDE:
256: case EPS_TARGET_REAL:
257: primme->target = primme_closest_abs;
258: primme->numTargetShifts = 1;
259: ops->target = PetscRealPart(eps->target);
260: primme->targetShifts = &ops->target;
261: break;
262: default:
263: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'which' value not supported by PRIMME");
264: break;
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: break;
280: }
282: /* If user sets mpd or ncv, maxBasisSize is modified */
283: if (eps->mpd!=PETSC_DEFAULT) {
284: primme->maxBasisSize = eps->mpd;
285: if (eps->ncv!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: 'ncv' is ignored by PRIMME\n"); }
286: } else if (eps->ncv!=PETSC_DEFAULT) primme->maxBasisSize = eps->ncv;
288: if (primme_set_method(ops->method,primme) < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME method not valid");
290: eps->mpd = primme->maxBasisSize;
291: eps->ncv = (primme->locking?eps->nev:0)+primme->maxBasisSize;
292: ops->bs = primme->maxBlockSize;
294: /* Set workspace */
295: EPSAllocateSolution(eps,0);
297: /* Setup the preconditioner */
298: if (primme->correctionParams.precondition) {
299: STGetKSP(eps->st,&ops->ksp);
300: PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&flg);
301: if (!flg) { 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) {
308: MatCreateVecsEmpty(ops->A,&ops->x,&ops->y);
309: PetscLogObjectParent((PetscObject)eps,(PetscObject)ops->x);
310: PetscLogObjectParent((PetscObject)eps,(PetscObject)ops->y);
311: }
312: return(0);
313: }
315: PetscErrorCode EPSSolve_PRIMME(EPS eps)
316: {
318: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
319: PetscScalar *a;
320: PetscInt i,ierrprimme;
321: PetscReal *evals,*rnorms;
324: /* Reset some parameters left from previous runs */
325: #if defined(SLEPC_HAVE_PRIMME2p2)
326: ops->primme.aNorm = 0.0;
327: #else
328: /* Force PRIMME to stop by absolute error */
329: ops->primme.aNorm = 1.0;
330: #endif
331: ops->primme.initSize = eps->nini;
332: ops->primme.iseed[0] = -1;
333: ops->primme.iseed[1] = -1;
334: ops->primme.iseed[2] = -1;
335: ops->primme.iseed[3] = -1;
337: /* Call PRIMME solver */
338: BVGetArray(eps->V,&a);
339: PetscMalloc2(eps->ncv,&evals,eps->ncv,&rnorms);
340: ierrprimme = PRIMME_DRIVER(evals,a,rnorms,&ops->primme);
341: for (i=0;i<eps->ncv;i++) eps->eigr[i] = evals[i];
342: for (i=0;i<eps->ncv;i++) eps->errest[i] = rnorms[i];
343: PetscFree2(evals,rnorms);
344: BVRestoreArray(eps->V,&a);
346: eps->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
347: eps->reason = eps->nconv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
348: eps->its = ops->primme.stats.numOuterIterations;
350: /* Process PRIMME error code */
351: if (ierrprimme == 0) {
352: /* no error */
353: } else if (ierrprimme == -1) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: unexpected error",ierrprimme);
354: else if (ierrprimme == -2) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: allocation error",ierrprimme);
355: else if (ierrprimme == -3) {
356: /* stop by maximum number of iteration or matvecs */
357: } 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);
358: else SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: runtime error; check PRIMME's manual",ierrprimme);
359: return(0);
360: }
362: PetscErrorCode EPSReset_PRIMME(EPS eps)
363: {
365: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
368: primme_free(&ops->primme);
369: VecDestroy(&ops->x);
370: VecDestroy(&ops->y);
371: return(0);
372: }
374: PetscErrorCode EPSDestroy_PRIMME(EPS eps)
375: {
379: PetscFree(eps->data);
380: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL);
381: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL);
382: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL);
383: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL);
384: return(0);
385: }
387: PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
388: {
390: PetscBool isascii;
391: EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
392: PetscMPIInt rank;
395: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
396: if (isascii) {
397: PetscViewerASCIIPrintf(viewer," block size=%D\n",ctx->bs);
398: PetscViewerASCIIPrintf(viewer," solver method: %s\n",EPSPRIMMEMethods[(EPSPRIMMEMethod)ctx->method]);
400: /* Display PRIMME params */
401: MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank);
402: if (!rank) primme_display_params(ctx->primme);
403: }
404: return(0);
405: }
407: PetscErrorCode EPSSetFromOptions_PRIMME(PetscOptionItems *PetscOptionsObject,EPS eps)
408: {
409: PetscErrorCode ierr;
410: EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
411: PetscInt bs;
412: EPSPRIMMEMethod meth;
413: PetscBool flg;
416: PetscOptionsHead(PetscOptionsObject,"EPS PRIMME Options");
418: PetscOptionsInt("-eps_primme_blocksize","Maximum block size","EPSPRIMMESetBlockSize",ctx->bs,&bs,&flg);
419: if (flg) { EPSPRIMMESetBlockSize(eps,bs); }
421: PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
422: if (flg) { EPSPRIMMESetMethod(eps,meth); }
424: PetscOptionsTail();
425: return(0);
426: }
428: static PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
429: {
430: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
433: if (bs == PETSC_DEFAULT) ops->bs = 0;
434: else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
435: else ops->bs = bs;
436: return(0);
437: }
439: /*@
440: EPSPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.
442: Logically Collective on eps
444: Input Parameters:
445: + eps - the eigenproblem solver context
446: - bs - block size
448: Options Database Key:
449: . -eps_primme_blocksize - Sets the max allowed block size value
451: Notes:
452: If the block size is not set, the value established by primme_initialize
453: is used.
455: The user should set the block size based on the architecture specifics
456: of the target computer, as well as any a priori knowledge of multiplicities.
457: The code does NOT require bs > 1 to find multiple eigenvalues. For some
458: methods, keeping bs = 1 yields the best overall performance.
460: Level: advanced
462: .seealso: EPSPRIMMEGetBlockSize()
463: @*/
464: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
465: {
471: PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
472: return(0);
473: }
475: static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
476: {
477: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
480: *bs = ops->bs;
481: return(0);
482: }
484: /*@
485: EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
487: Not Collective
489: Input Parameter:
490: . eps - the eigenproblem solver context
492: Output Parameter:
493: . bs - returned block size
495: Level: advanced
497: .seealso: EPSPRIMMESetBlockSize()
498: @*/
499: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
500: {
506: PetscUseMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
507: return(0);
508: }
510: static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
511: {
512: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
515: ops->method = (primme_preset_method)method;
516: return(0);
517: }
519: /*@
520: EPSPRIMMESetMethod - Sets the method for the PRIMME library.
522: Logically Collective on eps
524: Input Parameters:
525: + eps - the eigenproblem solver context
526: - method - method that will be used by PRIMME
528: Options Database Key:
529: . -eps_primme_method - Sets the method for the PRIMME library
531: Note:
532: If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.
534: Level: advanced
536: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
537: @*/
538: PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
539: {
545: PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
546: return(0);
547: }
549: static PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
550: {
551: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
554: *method = (EPSPRIMMEMethod)ops->method;
555: return(0);
556: }
558: /*@
559: EPSPRIMMEGetMethod - Gets the method for the PRIMME library.
561: Not Collective
563: Input Parameter:
564: . eps - the eigenproblem solver context
566: Output Parameter:
567: . method - method that will be used by PRIMME
569: Level: advanced
571: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
572: @*/
573: PetscErrorCode EPSPRIMMEGetMethod(EPS eps,EPSPRIMMEMethod *method)
574: {
580: PetscUseMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
581: return(0);
582: }
584: SLEPC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
585: {
587: EPS_PRIMME *primme;
590: PetscNewLog(eps,&primme);
591: eps->data = (void*)primme;
593: primme_initialize(&primme->primme);
594: primme->primme.globalSumReal = par_GlobalSumReal;
595: #if defined(SLEPC_HAVE_PRIMME3)
596: primme->primme.broadcastReal = par_broadcastReal;
597: #endif
598: #if defined(SLEPC_HAVE_PRIMME2p2)
599: primme->primme.convTestFun = convTestFun;
600: primme->primme.monitorFun = monitorFun;
601: #endif
602: primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;
604: eps->categ = EPS_CATEGORY_PRECOND;
606: eps->ops->solve = EPSSolve_PRIMME;
607: eps->ops->setup = EPSSetUp_PRIMME;
608: eps->ops->setupsort = EPSSetUpSort_Basic;
609: eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
610: eps->ops->destroy = EPSDestroy_PRIMME;
611: eps->ops->reset = EPSReset_PRIMME;
612: eps->ops->view = EPSView_PRIMME;
613: eps->ops->backtransform = EPSBackTransform_Default;
614: eps->ops->setdefaultst = EPSSetDefaultST_GMRES;
616: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME);
617: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME);
618: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME);
619: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME);
620: return(0);
621: }