Actual source code: primme.c
slepc-3.18.3 2023-03-24
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: {
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 = 0;
85: EPS eps = (EPS)primme->commInfo;
86: PetscInt i,k,nerrest;
88: switch (*event) {
89: case primme_event_outer_iteration:
90: /* Update EPS */
91: eps->its = primme->stats.numOuterIterations;
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: for (i=0;i<*blockSize;i++) {
125: PetscObjectComm((PetscObject)A),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);
126: PetscObjectComm((PetscObject)A),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);
127: PetscObjectComm((PetscObject)A),MatMult(A,x,y);
128: PetscObjectComm((PetscObject)A),VecResetArray(x);
129: PetscObjectComm((PetscObject)A),VecResetArray(y);
130: }
131: return;
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;
142: for (i=0;i<*blockSize;i++) {
143: PetscObjectComm((PetscObject)B),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);
144: PetscObjectComm((PetscObject)B),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);
145: PetscObjectComm((PetscObject)B),MatMult(B,x,y);
146: PetscObjectComm((PetscObject)B),VecResetArray(x);
147: PetscObjectComm((PetscObject)B),VecResetArray(y);
148: }
149: return;
150: }
151: #endif
153: static void applyPreconditioner_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
154: {
155: PetscInt i;
156: EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
157: Vec x = ops->x,y = ops->y;
159: for (i=0;i<*blockSize;i++) {
160: PetscObjectComm((PetscObject)ops->ksp),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);
161: PetscObjectComm((PetscObject)ops->ksp),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);
162: PetscObjectComm((PetscObject)ops->ksp),KSPSolve(ops->ksp,x,y);
163: PetscObjectComm((PetscObject)ops->ksp),VecResetArray(x);
164: PetscObjectComm((PetscObject)ops->ksp),VecResetArray(y);
165: }
166: return;
167: }
169: PetscErrorCode EPSSetUp_PRIMME(EPS eps)
170: {
171: PetscMPIInt numProcs,procID;
172: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
173: primme_params *primme = &ops->primme;
174: PetscBool flg;
176: EPSCheckHermitianDefinite(eps);
177: MPI_Comm_size(PetscObjectComm((PetscObject)eps),&numProcs);
178: MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&procID);
180: /* Check some constraints and set some default values */
181: if (eps->max_it==PETSC_DEFAULT) eps->max_it = PETSC_MAX_INT;
182: STGetMatrix(eps->st,0,&ops->A);
183: if (eps->isgeneralized) {
184: #if defined(SLEPC_HAVE_PRIMME3)
185: STGetMatrix(eps->st,1,&ops->B);
186: #else
187: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This version of PRIMME is not available for generalized problems");
188: #endif
189: }
190: EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
191: EPSCheckIgnored(eps,EPS_FEATURE_BALANCE);
192: if (!eps->which) eps->which = EPS_LARGEST_REAL;
193: #if !defined(SLEPC_HAVE_PRIMME2p2)
194: if (eps->converged != EPSConvergedAbsolute) PetscInfo(eps,"Warning: using absolute convergence test\n");
195: #else
196: EPSCheckIgnored(eps,EPS_FEATURE_CONVERGENCE);
197: #endif
199: /* Transfer SLEPc options to PRIMME options */
200: primme_free(primme);
201: primme_initialize(primme);
202: primme->n = eps->n;
203: primme->nLocal = eps->nloc;
204: primme->numEvals = eps->nev;
205: primme->matrix = ops;
206: primme->matrixMatvec = matrixMatvec_PRIMME;
207: #if defined(SLEPC_HAVE_PRIMME3)
208: if (eps->isgeneralized) {
209: primme->massMatrix = ops;
210: primme->massMatrixMatvec = massMatrixMatvec_PRIMME;
211: }
212: #endif
213: primme->commInfo = eps;
214: primme->maxOuterIterations = eps->max_it;
215: #if !defined(SLEPC_HAVE_PRIMME2p2)
216: primme->eps = SlepcDefaultTol(eps->tol);
217: #endif
218: primme->numProcs = numProcs;
219: primme->procID = procID;
220: primme->printLevel = 1;
221: primme->correctionParams.precondition = 1;
222: primme->globalSumReal = par_GlobalSumReal;
223: #if defined(SLEPC_HAVE_PRIMME3)
224: primme->broadcastReal = par_broadcastReal;
225: #endif
226: #if defined(SLEPC_HAVE_PRIMME2p2)
227: primme->convTestFun = convTestFun;
228: primme->monitorFun = monitorFun;
229: #endif
230: if (ops->bs > 0) primme->maxBlockSize = ops->bs;
232: switch (eps->which) {
233: case EPS_LARGEST_REAL:
234: primme->target = primme_largest;
235: break;
236: case EPS_SMALLEST_REAL:
237: primme->target = primme_smallest;
238: break;
239: case EPS_LARGEST_MAGNITUDE:
240: primme->target = primme_largest_abs;
241: ops->target = 0.0;
242: primme->numTargetShifts = 1;
243: primme->targetShifts = &ops->target;
244: break;
245: case EPS_SMALLEST_MAGNITUDE:
246: primme->target = primme_closest_abs;
247: ops->target = 0.0;
248: primme->numTargetShifts = 1;
249: primme->targetShifts = &ops->target;
250: break;
251: case EPS_TARGET_MAGNITUDE:
252: case EPS_TARGET_REAL:
253: primme->target = primme_closest_abs;
254: primme->numTargetShifts = 1;
255: ops->target = PetscRealPart(eps->target);
256: primme->targetShifts = &ops->target;
257: break;
258: default:
259: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'which' value not supported by PRIMME");
260: }
262: switch (eps->extraction) {
263: case EPS_RITZ:
264: primme->projectionParams.projection = primme_proj_RR;
265: break;
266: case EPS_HARMONIC:
267: primme->projectionParams.projection = primme_proj_harmonic;
268: break;
269: case EPS_REFINED:
270: primme->projectionParams.projection = primme_proj_refined;
271: break;
272: default:
273: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'extraction' value not supported by PRIMME");
274: }
276: /* If user sets mpd or ncv, maxBasisSize is modified */
277: if (eps->mpd!=PETSC_DEFAULT) {
278: primme->maxBasisSize = eps->mpd;
279: if (eps->ncv!=PETSC_DEFAULT) PetscInfo(eps,"Warning: 'ncv' is ignored by PRIMME\n");
280: } else if (eps->ncv!=PETSC_DEFAULT) primme->maxBasisSize = eps->ncv;
284: eps->mpd = primme->maxBasisSize;
285: eps->ncv = (primme->locking?eps->nev:0)+primme->maxBasisSize;
286: ops->bs = primme->maxBlockSize;
288: /* Set workspace */
289: EPSAllocateSolution(eps,0);
291: /* Setup the preconditioner */
292: if (primme->correctionParams.precondition) {
293: STGetKSP(eps->st,&ops->ksp);
294: PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&flg);
295: if (!flg) PetscInfo(eps,"Warning: ignoring KSP, should use KSPPREONLY\n");
296: primme->preconditioner = NULL;
297: primme->applyPreconditioner = applyPreconditioner_PRIMME;
298: }
300: /* Prepare auxiliary vectors */
301: if (!ops->x) MatCreateVecsEmpty(ops->A,&ops->x,&ops->y);
302: return 0;
303: }
305: PetscErrorCode EPSSolve_PRIMME(EPS eps)
306: {
307: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
308: PetscScalar *a;
309: PetscInt i,ierrprimme;
310: PetscReal *evals,*rnorms;
312: /* Reset some parameters left from previous runs */
313: #if defined(SLEPC_HAVE_PRIMME2p2)
314: ops->primme.aNorm = 0.0;
315: #else
316: /* Force PRIMME to stop by absolute error */
317: ops->primme.aNorm = 1.0;
318: #endif
319: ops->primme.initSize = eps->nini;
320: ops->primme.iseed[0] = -1;
321: ops->primme.iseed[1] = -1;
322: ops->primme.iseed[2] = -1;
323: ops->primme.iseed[3] = -1;
325: /* Call PRIMME solver */
326: BVGetArray(eps->V,&a);
327: PetscMalloc2(eps->ncv,&evals,eps->ncv,&rnorms);
328: ierrprimme = PRIMME_DRIVER(evals,a,rnorms,&ops->primme);
329: for (i=0;i<eps->ncv;i++) eps->eigr[i] = evals[i];
330: for (i=0;i<eps->ncv;i++) eps->errest[i] = rnorms[i];
331: PetscFree2(evals,rnorms);
332: BVRestoreArray(eps->V,&a);
334: eps->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
335: eps->reason = eps->nconv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
336: eps->its = ops->primme.stats.numOuterIterations;
338: /* Process PRIMME error code */
339: switch (ierrprimme) {
340: case 0: /* no error */
341: break;
342: case -1:
343: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": unexpected error",ierrprimme);
344: case -2:
345: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": allocation error",ierrprimme);
346: case -3: /* stop due to maximum number of iterations or matvecs */
347: break;
348: default:
351: }
352: return 0;
353: }
355: PetscErrorCode EPSReset_PRIMME(EPS eps)
356: {
357: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
359: primme_free(&ops->primme);
360: VecDestroy(&ops->x);
361: VecDestroy(&ops->y);
362: return 0;
363: }
365: PetscErrorCode EPSDestroy_PRIMME(EPS eps)
366: {
367: PetscFree(eps->data);
368: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL);
369: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL);
370: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL);
371: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL);
372: return 0;
373: }
375: PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
376: {
377: PetscBool isascii;
378: EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
379: PetscMPIInt rank;
381: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
382: if (isascii) {
383: PetscViewerASCIIPrintf(viewer," block size=%" PetscInt_FMT "\n",ctx->bs);
384: PetscViewerASCIIPrintf(viewer," solver method: %s\n",EPSPRIMMEMethods[(EPSPRIMMEMethod)ctx->method]);
386: /* Display PRIMME params */
387: MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank);
388: if (!rank) primme_display_params(ctx->primme);
389: }
390: return 0;
391: }
393: PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps,PetscOptionItems *PetscOptionsObject)
394: {
395: EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
396: PetscInt bs;
397: EPSPRIMMEMethod meth;
398: PetscBool flg;
400: PetscOptionsHeadBegin(PetscOptionsObject,"EPS PRIMME Options");
402: PetscOptionsInt("-eps_primme_blocksize","Maximum block size","EPSPRIMMESetBlockSize",ctx->bs,&bs,&flg);
403: if (flg) EPSPRIMMESetBlockSize(eps,bs);
405: PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
406: if (flg) EPSPRIMMESetMethod(eps,meth);
408: PetscOptionsHeadEnd();
409: return 0;
410: }
412: static PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
413: {
414: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
416: if (bs == PETSC_DEFAULT) ops->bs = 0;
417: else {
419: ops->bs = bs;
420: }
421: return 0;
422: }
424: /*@
425: EPSPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.
427: Logically Collective on eps
429: Input Parameters:
430: + eps - the eigenproblem solver context
431: - bs - block size
433: Options Database Key:
434: . -eps_primme_blocksize - Sets the max allowed block size value
436: Notes:
437: If the block size is not set, the value established by primme_initialize
438: is used.
440: The user should set the block size based on the architecture specifics
441: of the target computer, as well as any a priori knowledge of multiplicities.
442: The code does NOT require bs > 1 to find multiple eigenvalues. For some
443: methods, keeping bs = 1 yields the best overall performance.
445: Level: advanced
447: .seealso: EPSPRIMMEGetBlockSize()
448: @*/
449: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
450: {
453: PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
454: return 0;
455: }
457: static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
458: {
459: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
461: *bs = ops->bs;
462: return 0;
463: }
465: /*@
466: EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
468: Not Collective
470: Input Parameter:
471: . eps - the eigenproblem solver context
473: Output Parameter:
474: . bs - returned block size
476: Level: advanced
478: .seealso: EPSPRIMMESetBlockSize()
479: @*/
480: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
481: {
484: PetscUseMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
485: return 0;
486: }
488: static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
489: {
490: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
492: ops->method = (primme_preset_method)method;
493: return 0;
494: }
496: /*@
497: EPSPRIMMESetMethod - Sets the method for the PRIMME library.
499: Logically Collective on eps
501: Input Parameters:
502: + eps - the eigenproblem solver context
503: - method - method that will be used by PRIMME
505: Options Database Key:
506: . -eps_primme_method - Sets the method for the PRIMME library
508: Note:
509: If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.
511: Level: advanced
513: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
514: @*/
515: PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
516: {
519: PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
520: return 0;
521: }
523: static PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
524: {
525: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
527: *method = (EPSPRIMMEMethod)ops->method;
528: return 0;
529: }
531: /*@
532: EPSPRIMMEGetMethod - Gets the method for the PRIMME library.
534: Not Collective
536: Input Parameter:
537: . eps - the eigenproblem solver context
539: Output Parameter:
540: . method - method that will be used by PRIMME
542: Level: advanced
544: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
545: @*/
546: PetscErrorCode EPSPRIMMEGetMethod(EPS eps,EPSPRIMMEMethod *method)
547: {
550: PetscUseMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
551: return 0;
552: }
554: SLEPC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
555: {
556: EPS_PRIMME *primme;
558: PetscNew(&primme);
559: eps->data = (void*)primme;
561: primme_initialize(&primme->primme);
562: primme->primme.globalSumReal = par_GlobalSumReal;
563: #if defined(SLEPC_HAVE_PRIMME3)
564: primme->primme.broadcastReal = par_broadcastReal;
565: #endif
566: #if defined(SLEPC_HAVE_PRIMME2p2)
567: primme->primme.convTestFun = convTestFun;
568: primme->primme.monitorFun = monitorFun;
569: #endif
570: primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;
572: eps->categ = EPS_CATEGORY_PRECOND;
574: eps->ops->solve = EPSSolve_PRIMME;
575: eps->ops->setup = EPSSetUp_PRIMME;
576: eps->ops->setupsort = EPSSetUpSort_Basic;
577: eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
578: eps->ops->destroy = EPSDestroy_PRIMME;
579: eps->ops->reset = EPSReset_PRIMME;
580: eps->ops->view = EPSView_PRIMME;
581: eps->ops->backtransform = EPSBackTransform_Default;
582: eps->ops->setdefaultst = EPSSetDefaultST_GMRES;
584: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME);
585: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME);
586: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME);
587: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME);
588: return 0;
589: }