Actual source code: primme.c

  1: /*
  2:        This file implements a wrapper to the PRIMME library

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.
  9:       
 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 <petscsys.h>
 25: #include <private/epsimpl.h>    /*I "slepceps.h" I*/
 26: #include <private/stimpl.h>

 28: PetscErrorCode EPSSolve_PRIMME(EPS);

 30: EXTERN_C_BEGIN
 31: #include <primme.h>
 32: EXTERN_C_END

 34: typedef struct {
 35:   primme_params primme;           /* param struc */
 36:   primme_preset_method method;    /* primme method */
 37:   Mat       A;                    /* problem matrix */
 38:   EPS       eps;                  /* EPS current context */
 39:   KSP       ksp;                  /* linear solver and preconditioner */
 40:   Vec       x,y;                  /* auxiliary vectors */
 41:   PetscReal target;               /* a copy of eps's target */
 42: } EPS_PRIMME;

 44: EPSPRIMMEMethod methodN[] = {
 45:   EPS_PRIMME_DYNAMIC,
 46:   EPS_PRIMME_DEFAULT_MIN_TIME,
 47:   EPS_PRIMME_DEFAULT_MIN_MATVECS,
 48:   EPS_PRIMME_ARNOLDI,
 49:   EPS_PRIMME_GD,
 50:   EPS_PRIMME_GD_PLUSK,
 51:   EPS_PRIMME_GD_OLSEN_PLUSK,
 52:   EPS_PRIMME_JD_OLSEN_PLUSK,
 53:   EPS_PRIMME_RQI,
 54:   EPS_PRIMME_JDQR,
 55:   EPS_PRIMME_JDQMR,
 56:   EPS_PRIMME_JDQMR_ETOL,
 57:   EPS_PRIMME_SUBSPACE_ITERATION,
 58:   EPS_PRIMME_LOBPCG_ORTHOBASIS,
 59:   EPS_PRIMME_LOBPCG_ORTHOBASISW
 60: };

 62: static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme);
 63: static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme);

 65: void par_GlobalSumDouble(void *sendBuf,void *recvBuf,int *count,primme_params *primme)
 66: {
 68:   MPI_Allreduce((double*)sendBuf,(double*)recvBuf,*count,MPI_DOUBLE,MPI_SUM,((PetscObject)(primme->commInfo))->comm);CHKERRABORT(((PetscObject)(primme->commInfo))->comm,ierr);
 69: }

 73: PetscErrorCode EPSSetUp_PRIMME(EPS eps)
 74: {
 76:   PetscMPIInt    numProcs,procID;
 77:   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;
 78:   primme_params  *primme = &ops->primme;
 79:   PetscBool      t;

 82:   MPI_Comm_size(((PetscObject)eps)->comm,&numProcs);
 83:   MPI_Comm_rank(((PetscObject)eps)->comm,&procID);
 84: 
 85:   /* Check some constraints and set some default values */
 86:   if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
 87:   STGetOperators(eps->OP,&ops->A,PETSC_NULL);
 88:   if (!ops->A) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"The problem matrix has to be specified first");
 89:   if (!eps->ishermitian)
 90:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME is only available for Hermitian problems");
 91:   if (eps->isgeneralized)
 92:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME is not available for generalized problems");
 93:   if (!eps->which) eps->which = EPS_LARGEST_REAL;

 95:   /* Change the default sigma to inf if necessary */
 96:   if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
 97:       eps->which == EPS_LARGEST_IMAGINARY) {
 98:     STSetDefaultShift(eps->OP,3e300);
 99:   }

101:   /* Avoid setting the automatic shift when a target is set */
102:   STSetDefaultShift(eps->OP,0.0);

104:   STSetUp(eps->OP);
105:   PetscTypeCompare((PetscObject)eps->OP,STPRECOND,&t);
106:   if (!t) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME only works with STPRECOND");

108:   /* Transfer SLEPc options to PRIMME options */
109:   primme->n = eps->n;
110:   primme->nLocal = eps->nloc;
111:   primme->numEvals = eps->nev;
112:   primme->matrix = ops;
113:   primme->commInfo = eps;
114:   primme->maxMatvecs = eps->max_it;
115:   primme->eps = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
116:   primme->numProcs = numProcs;
117:   primme->procID = procID;
118:   primme->printLevel = 0;
119:   primme->correctionParams.precondition = 1;

121:   if (!eps->which) eps->which = EPS_LARGEST_REAL;
122:   switch(eps->which) {
123:     case EPS_LARGEST_REAL:
124:       primme->target = primme_largest;
125:       break;
126:     case EPS_SMALLEST_REAL:
127:       primme->target = primme_smallest;
128:       break;
129:     case EPS_TARGET_MAGNITUDE:
130:     case EPS_TARGET_REAL:
131:       primme->target = primme_closest_abs;
132:       primme->numTargetShifts = 1;
133:       ops->target = PetscRealPart(eps->target);
134:       primme->targetShifts = &ops->target;
135:       break;
136:     default:
137:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"'which' value does not supported by PRIMME");
138:       break;
139:   }
140: 
141:   if (primme_set_method(ops->method,primme) < 0)
142:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME method not valid");
143: 
144:   /* If user sets ncv, maxBasisSize is modified. If not, ncv is set as maxBasisSize */
145:   if (eps->ncv) primme->maxBasisSize = eps->ncv;
146:   else eps->ncv = primme->maxBasisSize;
147:   if (eps->ncv < eps->nev+primme->maxBlockSize)
148:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME needs ncv >= nev+maxBlockSize");
149:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }

151:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }

153:   /* Set workspace */
154:   EPSAllocateSolution(eps);

156:   /* Setup the preconditioner */
157:   ops->eps = eps;
158:   if (primme->correctionParams.precondition) {
159:     STGetKSP(eps->OP,&ops->ksp);
160:     PetscTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&t);
161:     if (!t) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME only works with KSPPREONLY");
162:     primme->preconditioner = PETSC_NULL;
163:     primme->applyPreconditioner = applyPreconditioner_PRIMME;
164:   }

166:   /* Prepare auxiliary vectors */
167:   VecCreateMPIWithArray(PETSC_COMM_WORLD,eps->nloc,eps->n,PETSC_NULL,&ops->x);
168:   VecCreateMPIWithArray(PETSC_COMM_WORLD,eps->nloc,eps->n,PETSC_NULL,&ops->y);
169: 
170:   /* dispatch solve method */
171:   if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
172:   eps->ops->solve = EPSSolve_PRIMME;
173:   return(0);
174: }

178: PetscErrorCode EPSSolve_PRIMME(EPS eps)
179: {
181:   EPS_PRIMME     *ops = (EPS_PRIMME *)eps->data;
182:   PetscScalar    *a;
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    = 0.0;
191:   ops->primme.initSize = eps->nini;
192:   ops->primme.iseed[0] = -1;

194:   /* Call PRIMME solver */
195:   VecGetArray(eps->V[0],&a);
196: #if !defined(PETSC_USE_COMPLEX)
197:   dprimme(eps->eigr,a,eps->errest,&ops->primme);
198: #else
199:   /* PRIMME returns real eigenvalues, but SLEPc works with complex ones */
200:   PetscMalloc(eps->ncv*sizeof(PetscReal),&evals);
201:   zprimme(evals,(Complex_Z*)a,eps->errest,&ops->primme);
202:   for (i=0;i<eps->ncv;i++)
203:     eps->eigr[i] = evals[i];
204:   PetscFree(evals);
205: #endif
206:   VecRestoreArray(eps->V[0],&a);
207: 
208:   switch(ierr) {
209:     case 0: /* Successful */
210:       break;

212:     case -1:
213:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: Failed to open output file");
214:       break;

216:     case -2:
217:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: Insufficient integer or real workspace allocated");
218:       break;

220:     case -3:
221:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: main_iter encountered a problem");
222:       break;

224:     default:
225:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: some parameters wrong configured");
226:       break;
227:   }

229:   eps->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
230:   eps->reason = eps->ncv >= eps->nev ? EPS_CONVERGED_TOL : EPS_DIVERGED_ITS;
231:   eps->its = ops->primme.stats.numOuterIterations;
232:   eps->OP->applys = ops->primme.stats.numMatvecs;
233:   return(0);
234: }

238: static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme)
239: {
241:   PetscInt       i,N = primme->n;
242:   EPS_PRIMME     *ops = (EPS_PRIMME *)primme->matrix;
243:   Vec            x = ops->x,y = ops->y;
244:   Mat            A = ops->A;

247:   for (i=0;i<*blockSize;i++) {
248:     /* build vectors using 'in' an 'out' workspace */
249:     VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);
250:     VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);

252:     MatMult(A,x,y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
253: 
254:     VecResetArray(x);CHKERRABORT(PETSC_COMM_WORLD,ierr);
255:     VecResetArray(y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
256:   }
257:   PetscFunctionReturnVoid();
258: }

262: static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme)
263: {
265:   PetscInt       i,N = primme->n,lits;
266:   EPS_PRIMME     *ops = (EPS_PRIMME *)primme->matrix;
267:   Vec            x = ops->x,y = ops->y;
268: 
270:   for (i=0;i<*blockSize;i++) {
271:     /* build vectors using 'in' an 'out' workspace */
272:     VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);
273:     VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);

275:     KSPSolve(ops->ksp,x,y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
276:     KSPGetIterationNumber(ops->ksp,&lits);CHKERRABORT(PETSC_COMM_WORLD,ierr);
277:     ops->eps->OP->lineariterations+= lits;
278: 
279:     VecResetArray(x);CHKERRABORT(PETSC_COMM_WORLD,ierr);
280:     VecResetArray(y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
281:   }
282:   PetscFunctionReturnVoid();
283: }

287: PetscErrorCode EPSReset_PRIMME(EPS eps)
288: {
290:   EPS_PRIMME     *ops = (EPS_PRIMME *)eps->data;

293:   primme_Free(&ops->primme);
294:   VecDestroy(&ops->x);
295:   VecDestroy(&ops->y);
296:   EPSFreeSolution(eps);
297:   return(0);
298: }

302: PetscErrorCode EPSDestroy_PRIMME(EPS eps)
303: {

307:   PetscFree(eps->data);
308:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetBlockSize_C","",PETSC_NULL);
309:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetMethod_C","",PETSC_NULL);
310:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetBlockSize_C","",PETSC_NULL);
311:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetMethod_C","",PETSC_NULL);
312:   return(0);
313: }

317: PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
318: {
319:   PetscErrorCode  ierr;
320:   PetscBool       isascii;
321:   primme_params   *primme = &((EPS_PRIMME *)eps->data)->primme;
322:   EPSPRIMMEMethod methodn;
323:   PetscMPIInt     rank;

326:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
327:   if (!isascii) {
328:     SETERRQ1(((PetscObject)eps)->comm,1,"Viewer type %s not supported for EPSPRIMME",((PetscObject)viewer)->type_name);
329:   }
330: 
331:   PetscViewerASCIIPrintf(viewer,"  PRIMME: block size=%d\n",primme->maxBlockSize);
332:   EPSPRIMMEGetMethod(eps,&methodn);
333:   PetscViewerASCIIPrintf(viewer,"  PRIMME: solver method: %s\n",EPSPRIMMEMethods[methodn]);

335:   /* Display PRIMME params */
336:   MPI_Comm_rank(((PetscObject)eps)->comm,&rank);
337:   if (!rank) primme_display_params(*primme);
338:   return(0);
339: }

343: PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps)
344: {
345:   PetscErrorCode  ierr;
346:   EPS_PRIMME      *ctx = (EPS_PRIMME *)eps->data;
347:   PetscInt        bs;
348:   EPSPRIMMEMethod meth;
349:   PetscBool       flg;

352:   PetscOptionsHead("EPS PRIMME Options");
353:   PetscOptionsInt("-eps_primme_block_size","Maximum block size","EPSPRIMMESetBlockSize",ctx->primme.maxBlockSize,&bs,&flg);
354:   if (flg) { EPSPRIMMESetBlockSize(eps,bs); }
355:   PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
356:   if (flg) { EPSPRIMMESetMethod(eps,meth); }
357:   PetscOptionsTail();
358:   return(0);
359: }

361: EXTERN_C_BEGIN
364: PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
365: {
366:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

369:   if (bs == PETSC_DEFAULT) ops->primme.maxBlockSize = 1;
370:   else if (bs <= 0) {
371:     SETERRQ(((PetscObject)eps)->comm,1,"PRIMME: wrong block size");
372:   } else ops->primme.maxBlockSize = bs;
373:   return(0);
374: }
375: EXTERN_C_END

379: /*@
380:    EPSPRIMMESetBlockSize - The maximum block size the code will try to use. 
381:    The user should set
382:    this based on the architecture specifics of the target computer, 
383:    as well as any a priori knowledge of multiplicities. The code does 
384:    NOT require BlockSize > 1 to find multiple eigenvalues.  For some 
385:    methods, keeping BlockSize = 1 yields the best overall performance.

387:    Collective on EPS

389:    Input Parameters:
390: +  eps - the eigenproblem solver context
391: -  bs - block size

393:    Options Database Key:
394: .  -eps_primme_block_size - Sets the max allowed block size value

396:    Notes:
397:    If the block size is not set, the value established by primme_initialize
398:    is used.

400:    Level: advanced
401: .seealso: EPSPRIMMEGetBlockSize()
402: @*/
403: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
404: {

410:   PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
411:   return(0);
412: }

414: EXTERN_C_BEGIN
417: PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
418: {
419:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

422:   if (bs) *bs = ops->primme.maxBlockSize;
423:   return(0);
424: }
425: EXTERN_C_END

429: /*@
430:    EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use. 

432:    Collective on EPS

434:    Input Parameters:
435: .  eps - the eigenproblem solver context
436:     
437:    Output Parameters:  
438: .  bs - returned block size 

440:    Level: advanced
441: .seealso: EPSPRIMMESetBlockSize()
442: @*/
443: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
444: {

449:   PetscTryMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
450:   return(0);
451: }

453: EXTERN_C_BEGIN
456: PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
457: {
458:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

461:   if (method == PETSC_DEFAULT) ops->method = DEFAULT_MIN_TIME;
462:   else ops->method = (primme_preset_method)method;
463:   return(0);
464: }
465: EXTERN_C_END

469: /*@
470:    EPSPRIMMESetMethod - Sets the method for the PRIMME library.

472:    Collective on EPS

474:    Input Parameters:
475: +  eps - the eigenproblem solver context
476: -  method - method that will be used by PRIMME. It must be one of:
477:     EPS_PRIMME_DYNAMIC, EPS_PRIMME_DEFAULT_MIN_TIME(EPS_PRIMME_JDQMR_ETOL),
478:     EPS_PRIMME_DEFAULT_MIN_MATVECS(EPS_PRIMME_GD_OLSEN_PLUSK), EPS_PRIMME_ARNOLDI,
479:     EPS_PRIMME_GD, EPS_PRIMME_GD_PLUSK, EPS_PRIMME_GD_OLSEN_PLUSK, 
480:     EPS_PRIMME_JD_OLSEN_PLUSK, EPS_PRIMME_RQI, EPS_PRIMME_JDQR, EPS_PRIMME_JDQMR, 
481:     EPS_PRIMME_JDQMR_ETOL, EPS_PRIMME_SUBSPACE_ITERATION,
482:     EPS_PRIMME_LOBPCG_ORTHOBASIS, EPS_PRIMME_LOBPCG_ORTHOBASISW

484:    Options Database Key:
485: .  -eps_primme_set_method - Sets the method for the PRIMME library.

487:    Note:
488:    If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.

490:    Level: advanced

492: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
493: @*/
494: PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
495: {

501:   PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
502:   return(0);
503: }

505: EXTERN_C_BEGIN
508: PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
509: {
510:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

513:   if (method) *method = (EPSPRIMMEMethod)ops->method;
514:   return(0);
515: }
516: EXTERN_C_END

520: /*@C
521:     EPSPRIMMEGetMethod - Gets the method for the PRIMME library.

523:     Mon Collective on EPS

525:    Input Parameters:
526: .  eps - the eigenproblem solver context
527:     
528:    Output Parameters: 
529: .  method - method that will be used by PRIMME. It must be one of:
530:     EPS_PRIMME_DYNAMIC, EPS_PRIMME_DEFAULT_MIN_TIME(EPS_PRIMME_JDQMR_ETOL),
531:     EPS_PRIMME_DEFAULT_MIN_MATVECS(EPS_PRIMME_GD_OLSEN_PLUSK), EPS_PRIMME_ARNOLDI,
532:     EPS_PRIMME_GD, EPS_PRIMME_GD_PLUSK, EPS_PRIMME_GD_OLSEN_PLUSK, 
533:     EPS_PRIMME_JD_OLSEN_PLUSK, EPS_PRIMME_RQI, EPS_PRIMME_JDQR, EPS_PRIMME_JDQMR, 
534:     EPS_PRIMME_JDQMR_ETOL, EPS_PRIMME_SUBSPACE_ITERATION,
535:     EPS_PRIMME_LOBPCG_ORTHOBASIS, EPS_PRIMME_LOBPCG_ORTHOBASISW

537:     Level: advanced

539: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
540: @*/
541: PetscErrorCode EPSPRIMMEGetMethod(EPS eps, EPSPRIMMEMethod *method)
542: {

547:   PetscTryMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
548:   return(0);
549: }

551: EXTERN_C_BEGIN
554: PetscErrorCode EPSCreate_PRIMME(EPS eps)
555: {
557:   EPS_PRIMME     *primme;

560:   STSetType(eps->OP,STPRECOND);
561:   STPrecondSetKSPHasMat(eps->OP,PETSC_TRUE);

563:   PetscNewLog(eps,EPS_PRIMME,&primme);
564:   eps->data                      = (void*)primme;
565:   eps->ops->setup                = EPSSetUp_PRIMME;
566:   eps->ops->setfromoptions       = EPSSetFromOptions_PRIMME;
567:   eps->ops->destroy              = EPSDestroy_PRIMME;
568:   eps->ops->reset                = EPSReset_PRIMME;
569:   eps->ops->view                 = EPSView_PRIMME;
570:   eps->ops->backtransform        = EPSBackTransform_Default;
571:   eps->ops->computevectors       = EPSComputeVectors_Default;

573:   primme_initialize(&primme->primme);
574:   primme->primme.matrixMatvec = multMatvec_PRIMME;
575:   primme->primme.globalSumDouble = par_GlobalSumDouble;
576:   primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;

578:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetBlockSize_C","EPSPRIMMESetBlockSize_PRIMME",EPSPRIMMESetBlockSize_PRIMME);
579:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetMethod_C","EPSPRIMMESetMethod_PRIMME",EPSPRIMMESetMethod_PRIMME);
580:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetBlockSize_C","EPSPRIMMEGetBlockSize_PRIMME",EPSPRIMMEGetBlockSize_PRIMME);
581:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetMethod_C","EPSPRIMMEGetMethod_PRIMME",EPSPRIMMEGetMethod_PRIMME);
582:   return(0);
583: }
584: EXTERN_C_END