Actual source code: primme.c
slepc-3.5.2 2014-10-10
1: /*
2: This file implements a wrapper to the PRIMME package
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
25: #include <slepc-private/stimpl.h>
27: PetscErrorCode EPSSolve_PRIMME(EPS);
29: EXTERN_C_BEGIN
30: #include <primme.h>
31: EXTERN_C_END
33: typedef struct {
34: primme_params primme; /* param struc */
35: primme_preset_method method; /* primme method */
36: Mat A; /* problem matrix */
37: EPS eps; /* EPS current context */
38: KSP ksp; /* linear solver and preconditioner */
39: Vec x,y; /* auxiliary vectors */
40: PetscReal target; /* a copy of eps's target */
41: } EPS_PRIMME;
43: EPSPRIMMEMethod methodN[] = {
44: EPS_PRIMME_DYNAMIC,
45: EPS_PRIMME_DEFAULT_MIN_TIME,
46: EPS_PRIMME_DEFAULT_MIN_MATVECS,
47: EPS_PRIMME_ARNOLDI,
48: EPS_PRIMME_GD,
49: EPS_PRIMME_GD_PLUSK,
50: EPS_PRIMME_GD_OLSEN_PLUSK,
51: EPS_PRIMME_JD_OLSEN_PLUSK,
52: EPS_PRIMME_RQI,
53: EPS_PRIMME_JDQR,
54: EPS_PRIMME_JDQMR,
55: EPS_PRIMME_JDQMR_ETOL,
56: EPS_PRIMME_SUBSPACE_ITERATION,
57: EPS_PRIMME_LOBPCG_ORTHOBASIS,
58: EPS_PRIMME_LOBPCG_ORTHOBASISW
59: };
61: static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme);
62: static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme);
64: static void par_GlobalSumDouble(void *sendBuf,void *recvBuf,int *count,primme_params *primme)
65: {
67: MPI_Allreduce((double*)sendBuf,(double*)recvBuf,*count,MPI_DOUBLE,MPI_SUM,PetscObjectComm((PetscObject)primme->commInfo));CHKERRABORT(PetscObjectComm((PetscObject)primme->commInfo),ierr);
68: }
72: PetscErrorCode EPSSetUp_PRIMME(EPS eps)
73: {
75: PetscMPIInt numProcs,procID;
76: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
77: primme_params *primme = &ops->primme;
78: PetscBool istrivial,flg;
81: MPI_Comm_size(PetscObjectComm((PetscObject)eps),&numProcs);
82: MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&procID);
84: /* Check some constraints and set some default values */
85: if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
86: STGetOperators(eps->st,0,&ops->A);
87: if (!ops->A) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"The problem matrix has to be specified first");
88: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME is only available for Hermitian problems");
89: if (eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME is not available for generalized problems");
90: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
91: if (!eps->which) eps->which = EPS_LARGEST_REAL;
92: if (eps->converged != EPSConvergedAbsolute) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME only supports absolute convergence test");
93: RGIsTrivial(eps->rg,&istrivial);
94: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
96: /* Set default sigma */
97: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL || eps->which == EPS_LARGEST_IMAGINARY) {
98: STSetDefaultShift(eps->st,3e300);
99: } else {
100: STSetDefaultShift(eps->st,0.0);
101: }
103: STSetUp(eps->st);
104: PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&flg);
105: if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME only works with STPRECOND");
107: /* Transfer SLEPc options to PRIMME options */
108: primme->n = eps->n;
109: primme->nLocal = eps->nloc;
110: primme->numEvals = eps->nev;
111: primme->matrix = ops;
112: primme->commInfo = eps;
113: primme->maxMatvecs = eps->max_it;
114: primme->eps = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
115: primme->numProcs = numProcs;
116: primme->procID = procID;
117: primme->printLevel = 0;
118: primme->correctionParams.precondition = 1;
120: switch (eps->which) {
121: case EPS_LARGEST_REAL:
122: primme->target = primme_largest;
123: break;
124: case EPS_SMALLEST_REAL:
125: primme->target = primme_smallest;
126: break;
127: case EPS_TARGET_MAGNITUDE:
128: case EPS_TARGET_REAL:
129: primme->target = primme_closest_abs;
130: primme->numTargetShifts = 1;
131: ops->target = PetscRealPart(eps->target);
132: primme->targetShifts = &ops->target;
133: break;
134: default:
135: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'which' value does not supported by PRIMME");
136: break;
137: }
139: if (primme_set_method(ops->method,primme) < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME method not valid");
141: /* If user sets ncv, maxBasisSize is modified. If not, ncv is set as maxBasisSize */
142: if (eps->ncv) primme->maxBasisSize = eps->ncv;
143: else eps->ncv = primme->maxBasisSize;
144: if (eps->ncv < eps->nev+primme->maxBlockSize) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME needs ncv >= nev+maxBlockSize");
145: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
147: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
149: /* Set workspace */
150: EPSAllocateSolution(eps,0);
151: PetscObjectTypeCompare((PetscObject)eps->V,BVVECS,&flg);
152: if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires a BV with contiguous storage");
154: /* Setup the preconditioner */
155: ops->eps = eps;
156: if (primme->correctionParams.precondition) {
157: STGetKSP(eps->st,&ops->ksp);
158: PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&flg);
159: if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME only works with KSPPREONLY");
160: primme->preconditioner = NULL;
161: primme->applyPreconditioner = applyPreconditioner_PRIMME;
162: }
164: /* Prepare auxiliary vectors */
165: VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,eps->n,NULL,&ops->x);
166: VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,eps->n,NULL,&ops->y);
167: PetscLogObjectParent((PetscObject)eps,(PetscObject)ops->x);
168: PetscLogObjectParent((PetscObject)eps,(PetscObject)ops->y);
170: /* dispatch solve method */
171: eps->ops->solve = EPSSolve_PRIMME;
172: return(0);
173: }
177: PetscErrorCode EPSSolve_PRIMME(EPS eps)
178: {
180: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
181: PetscScalar *a;
182: Vec v0;
183: #if defined(PETSC_USE_COMPLEX)
184: PetscInt i;
185: PetscReal *evals;
186: #endif
189: /* Reset some parameters left from previous runs */
190: ops->primme.aNorm = 1.0;
191: ops->primme.initSize = eps->nini;
192: ops->primme.iseed[0] = -1;
194: /* Call PRIMME solver */
195: BVGetColumn(eps->V,0,&v0);
196: VecGetArray(v0,&a);
197: #if !defined(PETSC_USE_COMPLEX)
198: dprimme(eps->eigr,a,eps->errest,&ops->primme);
199: #else
200: /* PRIMME returns real eigenvalues, but SLEPc works with complex ones */
201: PetscMalloc1(eps->ncv,&evals);
202: zprimme(evals,(Complex_Z*)a,eps->errest,&ops->primme);
203: for (i=0;i<eps->ncv;i++)
204: eps->eigr[i] = evals[i];
205: PetscFree(evals);
206: #endif
207: VecRestoreArray(v0,&a);
208: BVRestoreColumn(eps->V,0,&v0);
210: switch (ierr) {
211: case 0: /* Successful */
212: break;
213: case -1:
214: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME: Failed to open output file");
215: break;
216: case -2:
217: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME: Insufficient integer or real workspace allocated");
218: break;
219: case -3:
220: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME: main_iter encountered a problem");
221: break;
222: default:
223: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME: some parameters wrong configured");
224: break;
225: }
227: eps->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
228: eps->reason = eps->ncv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
229: eps->its = ops->primme.stats.numOuterIterations;
230: eps->st->applys = ops->primme.stats.numMatvecs;
231: return(0);
232: }
236: static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme)
237: {
239: PetscInt i,N = primme->n;
240: EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
241: Vec x = ops->x,y = ops->y;
242: Mat A = ops->A;
245: for (i=0;i<*blockSize;i++) {
246: /* build vectors using 'in' an 'out' workspace */
247: VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
248: VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
250: MatMult(A,x,y);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
252: VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
253: VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
254: }
255: PetscFunctionReturnVoid();
256: }
260: static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme)
261: {
263: PetscInt i,N = primme->n,lits;
264: EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
265: Vec x = ops->x,y = ops->y;
268: for (i=0;i<*blockSize;i++) {
269: /* build vectors using 'in' an 'out' workspace */
270: VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
271: VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
273: KSPSolve(ops->ksp,x,y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
274: KSPGetIterationNumber(ops->ksp,&lits);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
275: ops->eps->st->linearits += lits;
277: VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
278: VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
279: }
280: PetscFunctionReturnVoid();
281: }
285: PetscErrorCode EPSReset_PRIMME(EPS eps)
286: {
288: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
291: primme_Free(&ops->primme);
292: VecDestroy(&ops->x);
293: VecDestroy(&ops->y);
294: return(0);
295: }
299: PetscErrorCode EPSDestroy_PRIMME(EPS eps)
300: {
304: PetscFree(eps->data);
305: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL);
306: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL);
307: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL);
308: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL);
309: return(0);
310: }
314: PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
315: {
316: PetscErrorCode ierr;
317: PetscBool isascii;
318: primme_params *primme = &((EPS_PRIMME*)eps->data)->primme;
319: EPSPRIMMEMethod methodn;
320: PetscMPIInt rank;
323: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
324: if (isascii) {
325: PetscViewerASCIIPrintf(viewer," PRIMME: block size=%d\n",primme->maxBlockSize);
326: EPSPRIMMEGetMethod(eps,&methodn);
327: PetscViewerASCIIPrintf(viewer," PRIMME: solver method: %s\n",EPSPRIMMEMethods[methodn]);
329: /* Display PRIMME params */
330: MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank);
331: if (!rank) primme_display_params(*primme);
332: }
333: return(0);
334: }
338: PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps)
339: {
340: PetscErrorCode ierr;
341: EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
342: PetscInt bs;
343: EPSPRIMMEMethod meth;
344: PetscBool flg;
345: KSP ksp;
348: PetscOptionsHead("EPS PRIMME Options");
349: PetscOptionsInt("-eps_primme_block_size","Maximum block size","EPSPRIMMESetBlockSize",ctx->primme.maxBlockSize,&bs,&flg);
350: if (flg) {
351: EPSPRIMMESetBlockSize(eps,bs);
352: }
353: PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
354: if (flg) {
355: EPSPRIMMESetMethod(eps,meth);
356: }
358: /* Set STPrecond as the default ST */
359: if (!((PetscObject)eps->st)->type_name) {
360: STSetType(eps->st,STPRECOND);
361: }
362: STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);
364: /* Set the default options of the KSP */
365: STGetKSP(eps->st,&ksp);
366: if (!((PetscObject)ksp)->type_name) {
367: KSPSetType(ksp,KSPPREONLY);
368: }
369: PetscOptionsTail();
370: return(0);
371: }
375: static PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
376: {
377: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
380: if (bs == PETSC_DEFAULT) ops->primme.maxBlockSize = 1;
381: else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
382: else ops->primme.maxBlockSize = bs;
383: return(0);
384: }
388: /*@
389: EPSPRIMMESetBlockSize - The maximum block size the code will try to use.
390: The user should set
391: this based on the architecture specifics of the target computer,
392: as well as any a priori knowledge of multiplicities. The code does
393: NOT require BlockSize > 1 to find multiple eigenvalues. For some
394: methods, keeping BlockSize = 1 yields the best overall performance.
396: Logically Collective on EPS
398: Input Parameters:
399: + eps - the eigenproblem solver context
400: - bs - block size
402: Options Database Key:
403: . -eps_primme_block_size - Sets the max allowed block size value
405: Notes:
406: If the block size is not set, the value established by primme_initialize
407: is used.
409: Level: advanced
411: .seealso: EPSPRIMMEGetBlockSize()
412: @*/
413: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
414: {
420: PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
421: return(0);
422: }
426: static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
427: {
428: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
431: if (bs) *bs = ops->primme.maxBlockSize;
432: return(0);
433: }
437: /*@
438: EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
440: Not Collective
442: Input Parameter:
443: . eps - the eigenproblem solver context
445: Output Parameter:
446: . bs - returned block size
448: Level: advanced
450: .seealso: EPSPRIMMESetBlockSize()
451: @*/
452: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
453: {
458: PetscTryMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
459: return(0);
460: }
464: static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
465: {
466: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
469: if (method == PETSC_DEFAULT) ops->method = DEFAULT_MIN_TIME;
470: else ops->method = (primme_preset_method)method;
471: return(0);
472: }
476: /*@
477: EPSPRIMMESetMethod - Sets the method for the PRIMME library.
479: Logically Collective on EPS
481: Input Parameters:
482: + eps - the eigenproblem solver context
483: - method - method that will be used by PRIMME. It must be one of:
484: EPS_PRIMME_DYNAMIC, EPS_PRIMME_DEFAULT_MIN_TIME(EPS_PRIMME_JDQMR_ETOL),
485: EPS_PRIMME_DEFAULT_MIN_MATVECS(EPS_PRIMME_GD_OLSEN_PLUSK), EPS_PRIMME_ARNOLDI,
486: EPS_PRIMME_GD, EPS_PRIMME_GD_PLUSK, EPS_PRIMME_GD_OLSEN_PLUSK,
487: EPS_PRIMME_JD_OLSEN_PLUSK, EPS_PRIMME_RQI, EPS_PRIMME_JDQR, EPS_PRIMME_JDQMR,
488: EPS_PRIMME_JDQMR_ETOL, EPS_PRIMME_SUBSPACE_ITERATION,
489: EPS_PRIMME_LOBPCG_ORTHOBASIS, EPS_PRIMME_LOBPCG_ORTHOBASISW
491: Options Database Key:
492: . -eps_primme_method - Sets the method for the PRIMME library (one of
493: 'dynamic', 'default_min_time', 'default_min_matvecs', 'arnoldi',
494: 'gd', 'gd_plusk', 'gd_olsen_plusk', 'jd_olsen_plusk', 'rqi', 'jdqr', 'jdqmr',
495: 'jdqmr_etol', 'subspace_iteration', 'lobpcg_orthobasis', 'lobpcg_orthobasisw').
497: Note:
498: If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.
500: Level: advanced
502: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
503: @*/
504: PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
505: {
511: PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
512: return(0);
513: }
517: static PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
518: {
519: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
522: if (method) *method = (EPSPRIMMEMethod)ops->method;
523: return(0);
524: }
528: /*@C
529: EPSPRIMMEGetMethod - Gets the method for the PRIMME library.
531: Not Collective
533: Input Parameter:
534: . eps - the eigenproblem solver context
536: Output Parameter:
537: . method - method that will be used by PRIMME, one of
538: EPS_PRIMME_DYNAMIC, EPS_PRIMME_DEFAULT_MIN_TIME(EPS_PRIMME_JDQMR_ETOL),
539: EPS_PRIMME_DEFAULT_MIN_MATVECS(EPS_PRIMME_GD_OLSEN_PLUSK), EPS_PRIMME_ARNOLDI,
540: EPS_PRIMME_GD, EPS_PRIMME_GD_PLUSK, EPS_PRIMME_GD_OLSEN_PLUSK,
541: EPS_PRIMME_JD_OLSEN_PLUSK, EPS_PRIMME_RQI, EPS_PRIMME_JDQR, EPS_PRIMME_JDQMR,
542: EPS_PRIMME_JDQMR_ETOL, EPS_PRIMME_SUBSPACE_ITERATION,
543: EPS_PRIMME_LOBPCG_ORTHOBASIS, EPS_PRIMME_LOBPCG_ORTHOBASISW
545: Level: advanced
547: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
548: @*/
549: PetscErrorCode EPSPRIMMEGetMethod(EPS eps,EPSPRIMMEMethod *method)
550: {
555: PetscTryMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
556: return(0);
557: }
561: PETSC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
562: {
564: EPS_PRIMME *primme;
567: PetscNewLog(eps,&primme);
568: eps->data = (void*)primme;
570: eps->ops->setup = EPSSetUp_PRIMME;
571: eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
572: eps->ops->destroy = EPSDestroy_PRIMME;
573: eps->ops->reset = EPSReset_PRIMME;
574: eps->ops->view = EPSView_PRIMME;
575: eps->ops->backtransform = EPSBackTransform_Default;
577: primme_initialize(&primme->primme);
578: primme->primme.matrixMatvec = multMatvec_PRIMME;
579: primme->primme.globalSumDouble = par_GlobalSumDouble;
580: primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;
582: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME);
583: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME);
584: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME);
585: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME);
586: return(0);
587: }