Actual source code: primme.c

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

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2013, 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,primme->commInfo);CHKERRABORT(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      t;

 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");

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

 99:   /* Avoid setting the automatic shift when a target is set */
100:   STSetDefaultShift(eps->st,0.0);

102:   STSetUp(eps->st);
103:   PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&t);
104:   if (!t) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME only works with STPRECOND");

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

119:   if (!eps->which) eps->which = EPS_LARGEST_REAL;
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);

152:   /* Setup the preconditioner */
153:   ops->eps = eps;
154:   if (primme->correctionParams.precondition) {
155:     STGetKSP(eps->st,&ops->ksp);
156:     PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&t);
157:     if (!t) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME only works with KSPPREONLY");
158:     primme->preconditioner = NULL;
159:     primme->applyPreconditioner = applyPreconditioner_PRIMME;
160:   }

162:   /* Prepare auxiliary vectors */
163:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,eps->n,NULL,&ops->x);
164:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,eps->n,NULL,&ops->y);
165:   PetscLogObjectParent(eps,ops->x);
166:   PetscLogObjectParent(eps,ops->y);

168:   /* dispatch solve method */
169:   if (eps->leftvecs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Left vectors not supported in this solver");
170:   eps->ops->solve = EPSSolve_PRIMME;
171:   return(0);
172: }

176: PetscErrorCode EPSSolve_PRIMME(EPS eps)
177: {
179:   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;
180:   PetscScalar    *a;
181: #if defined(PETSC_USE_COMPLEX)
182:   PetscInt       i;
183:   PetscReal      *evals;
184: #endif

187:   /* Reset some parameters left from previous runs */
188:   ops->primme.aNorm    = 1.0;
189:   ops->primme.initSize = eps->nini;
190:   ops->primme.iseed[0] = -1;

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

206:   switch (ierr) {
207:     case 0: /* Successful */
208:       break;
209:     case -1:
210:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME: Failed to open output file");
211:       break;
212:     case -2:
213:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME: Insufficient integer or real workspace allocated");
214:       break;
215:     case -3:
216:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME: main_iter encountered a problem");
217:       break;
218:     default:
219:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME: some parameters wrong configured");
220:       break;
221:   }

223:   eps->nconv      = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
224:   eps->reason     = eps->ncv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
225:   eps->its        = ops->primme.stats.numOuterIterations;
226:   eps->st->applys = ops->primme.stats.numMatvecs;
227:   return(0);
228: }

232: static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme)
233: {
235:   PetscInt       i,N = primme->n;
236:   EPS_PRIMME     *ops = (EPS_PRIMME*)primme->matrix;
237:   Vec            x = ops->x,y = ops->y;
238:   Mat            A = ops->A;

241:   for (i=0;i<*blockSize;i++) {
242:     /* build vectors using 'in' an 'out' workspace */
243:     VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
244:     VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);

246:     MatMult(A,x,y);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);

248:     VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
249:     VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
250:   }
251:   PetscFunctionReturnVoid();
252: }

256: static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme)
257: {
259:   PetscInt       i,N = primme->n,lits;
260:   EPS_PRIMME     *ops = (EPS_PRIMME*)primme->matrix;
261:   Vec            x = ops->x,y = ops->y;

264:   for (i=0;i<*blockSize;i++) {
265:     /* build vectors using 'in' an 'out' workspace */
266:     VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
267:     VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);

269:     KSPSolve(ops->ksp,x,y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
270:     KSPGetIterationNumber(ops->ksp,&lits);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
271:     ops->eps->st->lineariterations+= lits;

273:     VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
274:     VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
275:   }
276:   PetscFunctionReturnVoid();
277: }

281: PetscErrorCode EPSReset_PRIMME(EPS eps)
282: {
284:   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;

287:   primme_Free(&ops->primme);
288:   VecDestroy(&ops->x);
289:   VecDestroy(&ops->y);
290:   EPSFreeSolution(eps);
291:   return(0);
292: }

296: PetscErrorCode EPSDestroy_PRIMME(EPS eps)
297: {

301:   PetscFree(eps->data);
302:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL);
303:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL);
304:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL);
305:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL);
306:   return(0);
307: }

311: PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
312: {
313:   PetscErrorCode  ierr;
314:   PetscBool       isascii;
315:   primme_params   *primme = &((EPS_PRIMME*)eps->data)->primme;
316:   EPSPRIMMEMethod methodn;
317:   PetscMPIInt     rank;

320:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
321:   if (isascii) {
322:     PetscViewerASCIIPrintf(viewer,"  PRIMME: block size=%d\n",primme->maxBlockSize);
323:     EPSPRIMMEGetMethod(eps,&methodn);
324:     PetscViewerASCIIPrintf(viewer,"  PRIMME: solver method: %s\n",EPSPRIMMEMethods[methodn]);

326:     /* Display PRIMME params */
327:     MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank);
328:     if (!rank) primme_display_params(*primme);
329:   }
330:   return(0);
331: }

335: PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps)
336: {
337:   PetscErrorCode  ierr;
338:   EPS_PRIMME      *ctx = (EPS_PRIMME*)eps->data;
339:   PetscInt        bs;
340:   EPSPRIMMEMethod meth;
341:   PetscBool       flg;
342:   KSP             ksp;

345:   PetscOptionsHead("EPS PRIMME Options");
346:   PetscOptionsInt("-eps_primme_block_size","Maximum block size","EPSPRIMMESetBlockSize",ctx->primme.maxBlockSize,&bs,&flg);
347:   if (flg) {
348:     EPSPRIMMESetBlockSize(eps,bs);
349:   }
350:   PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
351:   if (flg) {
352:     EPSPRIMMESetMethod(eps,meth);
353:   }

355:   /* Set STPrecond as the default ST */
356:   if (!((PetscObject)eps->st)->type_name) {
357:     STSetType(eps->st,STPRECOND);
358:   }
359:   STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);

361:   /* Set the default options of the KSP */
362:   STGetKSP(eps->st,&ksp);
363:   if (!((PetscObject)ksp)->type_name) {
364:     KSPSetType(ksp,KSPPREONLY);
365:   }
366:   PetscOptionsTail();
367:   return(0);
368: }

372: static PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
373: {
374:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

377:   if (bs == PETSC_DEFAULT) ops->primme.maxBlockSize = 1;
378:   else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
379:   else ops->primme.maxBlockSize = bs;
380:   return(0);
381: }

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

393:    Collective on EPS

395:    Input Parameters:
396: +  eps - the eigenproblem solver context
397: -  bs - block size

399:    Options Database Key:
400: .  -eps_primme_block_size - Sets the max allowed block size value

402:    Notes:
403:    If the block size is not set, the value established by primme_initialize
404:    is used.

406:    Level: advanced
407: .seealso: EPSPRIMMEGetBlockSize()
408: @*/
409: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
410: {

416:   PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
417:   return(0);
418: }

422: static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
423: {
424:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

427:   if (bs) *bs = ops->primme.maxBlockSize;
428:   return(0);
429: }

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

436:    Collective on EPS

438:    Input Parameters:
439: .  eps - the eigenproblem solver context

441:    Output Parameters:
442: .  bs - returned block size

444:    Level: advanced
445: .seealso: EPSPRIMMESetBlockSize()
446: @*/
447: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
448: {

453:   PetscTryMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
454:   return(0);
455: }

459: static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
460: {
461:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

464:   if (method == PETSC_DEFAULT) ops->method = DEFAULT_MIN_TIME;
465:   else ops->method = (primme_preset_method)method;
466:   return(0);
467: }

471: /*@
472:    EPSPRIMMESetMethod - Sets the method for the PRIMME library.

474:    Collective on EPS

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

486:    Options Database Key:
487: .  -eps_primme_set_method - Sets the method for the PRIMME library.

489:    Note:
490:    If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.

492:    Level: advanced

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

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

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

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

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

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: }

553: PETSC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
554: {
556:   EPS_PRIMME     *primme;

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

569:   primme_initialize(&primme->primme);
570:   primme->primme.matrixMatvec = multMatvec_PRIMME;
571:   primme->primme.globalSumDouble = par_GlobalSumDouble;
572:   primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;

574:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME);
575:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME);
576:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME);
577:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME);
578:   return(0);
579: }