Actual source code: primme.c

slepc-3.20.2 2024-03-15
Report Typos and Errors
  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: {
 65:   PetscErrorCode ierr;
 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 = PETSC_SUCCESS;
 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:   PetscFunctionBegin;
125:   for (i=0;i<*blockSize;i++) {
126:     PetscCallAbort(PetscObjectComm((PetscObject)A),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i));
127:     PetscCallAbort(PetscObjectComm((PetscObject)A),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i));
128:     PetscCallAbort(PetscObjectComm((PetscObject)A),MatMult(A,x,y));
129:     PetscCallAbort(PetscObjectComm((PetscObject)A),VecResetArray(x));
130:     PetscCallAbort(PetscObjectComm((PetscObject)A),VecResetArray(y));
131:   }
132:   PetscFunctionReturnVoid();
133: }

135: #if defined(SLEPC_HAVE_PRIMME3)
136: static void massMatrixMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
137: {
138:   PetscInt   i;
139:   EPS_PRIMME *ops = (EPS_PRIMME*)primme->massMatrix;
140:   Vec        x = ops->x,y = ops->y;
141:   Mat        B = ops->B;

143:   PetscFunctionBegin;
144:   for (i=0;i<*blockSize;i++) {
145:     PetscCallAbort(PetscObjectComm((PetscObject)B),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i));
146:     PetscCallAbort(PetscObjectComm((PetscObject)B),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i));
147:     PetscCallAbort(PetscObjectComm((PetscObject)B),MatMult(B,x,y));
148:     PetscCallAbort(PetscObjectComm((PetscObject)B),VecResetArray(x));
149:     PetscCallAbort(PetscObjectComm((PetscObject)B),VecResetArray(y));
150:   }
151:   PetscFunctionReturnVoid();
152: }
153: #endif

155: static void applyPreconditioner_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
156: {
157:   PetscInt   i;
158:   EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
159:   Vec        x = ops->x,y = ops->y;

161:   PetscFunctionBegin;
162:   for (i=0;i<*blockSize;i++) {
163:     PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i));
164:     PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i));
165:     PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),KSPSolve(ops->ksp,x,y));
166:     PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecResetArray(x));
167:     PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecResetArray(y));
168:   }
169:   PetscFunctionReturnVoid();
170: }

172: static PetscErrorCode EPSSetUp_PRIMME(EPS eps)
173: {
174:   PetscMPIInt    numProcs,procID;
175:   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;
176:   primme_params  *primme = &ops->primme;
177:   PetscBool      flg;

179:   PetscFunctionBegin;
180:   EPSCheckHermitianDefinite(eps);
181:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&numProcs));
182:   PetscCallMPI(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:   PetscCall(STGetMatrix(eps->st,0,&ops->A));
187:   if (eps->isgeneralized) {
188: #if defined(SLEPC_HAVE_PRIMME3)
189:     PetscCall(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) PetscCall(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                           = SlepcDefaultTol(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:   }

266:   switch (eps->extraction) {
267:     case EPS_RITZ:
268:       primme->projectionParams.projection = primme_proj_RR;
269:       break;
270:     case EPS_HARMONIC:
271:       primme->projectionParams.projection = primme_proj_harmonic;
272:       break;
273:     case EPS_REFINED:
274:       primme->projectionParams.projection = primme_proj_refined;
275:       break;
276:     default:
277:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'extraction' value not supported by PRIMME");
278:   }

280:   /* If user sets mpd or ncv, maxBasisSize is modified */
281:   if (eps->mpd!=PETSC_DEFAULT) {
282:     primme->maxBasisSize = eps->mpd;
283:     if (eps->ncv!=PETSC_DEFAULT) PetscCall(PetscInfo(eps,"Warning: 'ncv' is ignored by PRIMME\n"));
284:   } else if (eps->ncv!=PETSC_DEFAULT) primme->maxBasisSize = eps->ncv;

286:   PetscCheck(primme_set_method(ops->method,primme)>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME method not valid");

288:   eps->mpd = primme->maxBasisSize;
289:   eps->ncv = (primme->locking?eps->nev:0)+primme->maxBasisSize;
290:   ops->bs  = primme->maxBlockSize;

292:   /* Set workspace */
293:   PetscCall(EPSAllocateSolution(eps,0));

295:   /* Setup the preconditioner */
296:   if (primme->correctionParams.precondition) {
297:     PetscCall(STGetKSP(eps->st,&ops->ksp));
298:     PetscCall(PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&flg));
299:     if (!flg) PetscCall(PetscInfo(eps,"Warning: ignoring KSP, should use KSPPREONLY\n"));
300:     primme->preconditioner = NULL;
301:     primme->applyPreconditioner = applyPreconditioner_PRIMME;
302:   }

304:   /* Prepare auxiliary vectors */
305:   if (!ops->x) PetscCall(MatCreateVecsEmpty(ops->A,&ops->x,&ops->y));
306:   PetscFunctionReturn(PETSC_SUCCESS);
307: }

309: static PetscErrorCode EPSSolve_PRIMME(EPS eps)
310: {
311:   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;
312:   PetscScalar    *a;
313:   PetscInt       i,ld,ierrprimme;
314:   PetscReal      *evals,*rnorms;

316:   PetscFunctionBegin;
317:   /* Reset some parameters left from previous runs */
318: #if defined(SLEPC_HAVE_PRIMME2p2)
319:   ops->primme.aNorm    = 0.0;
320: #else
321:   /* Force PRIMME to stop by absolute error */
322:   ops->primme.aNorm    = 1.0;
323: #endif
324:   ops->primme.initSize = eps->nini;
325:   ops->primme.iseed[0] = -1;
326:   ops->primme.iseed[1] = -1;
327:   ops->primme.iseed[2] = -1;
328:   ops->primme.iseed[3] = -1;
329:   PetscCall(BVGetLeadingDimension(eps->V,&ld));
330:   ops->primme.ldevecs  = ld;

332:   /* Call PRIMME solver */
333:   PetscCall(BVGetArray(eps->V,&a));
334:   PetscCall(PetscMalloc2(eps->ncv,&evals,eps->ncv,&rnorms));
335:   ierrprimme = PRIMME_DRIVER(evals,a,rnorms,&ops->primme);
336:   for (i=0;i<eps->ncv;i++) eps->eigr[i] = evals[i];
337:   for (i=0;i<eps->ncv;i++) eps->errest[i] = rnorms[i];
338:   PetscCall(PetscFree2(evals,rnorms));
339:   PetscCall(BVRestoreArray(eps->V,&a));

341:   eps->nconv  = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
342:   eps->reason = eps->nconv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
343:   eps->its    = ops->primme.stats.numOuterIterations;

345:   /* Process PRIMME error code */
346:   switch (ierrprimme) {
347:     case 0: /* no error */
348:       break;
349:     case -1:
350:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": unexpected error",ierrprimme);
351:     case -2:
352:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": allocation error",ierrprimme);
353:     case -3: /* stop due to maximum number of iterations or matvecs */
354:       break;
355:     default:
356:       PetscCheck(ierrprimme<-39,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": configuration error; check PRIMME's manual",ierrprimme);
357:       PetscCheck(ierrprimme>=-39,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": runtime error; check PRIMME's manual",ierrprimme);
358:   }
359:   PetscFunctionReturn(PETSC_SUCCESS);
360: }

362: static PetscErrorCode EPSReset_PRIMME(EPS eps)
363: {
364:   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;

366:   PetscFunctionBegin;
367:   primme_free(&ops->primme);
368:   PetscCall(VecDestroy(&ops->x));
369:   PetscCall(VecDestroy(&ops->y));
370:   PetscFunctionReturn(PETSC_SUCCESS);
371: }

373: static PetscErrorCode EPSDestroy_PRIMME(EPS eps)
374: {
375:   PetscFunctionBegin;
376:   PetscCall(PetscFree(eps->data));
377:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL));
378:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL));
379:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL));
380:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL));
381:   PetscFunctionReturn(PETSC_SUCCESS);
382: }

384: static PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
385: {
386:   PetscBool      isascii;
387:   EPS_PRIMME     *ctx = (EPS_PRIMME*)eps->data;
388:   PetscMPIInt    rank;

390:   PetscFunctionBegin;
391:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
392:   if (isascii) {
393:     PetscCall(PetscViewerASCIIPrintf(viewer,"  block size=%" PetscInt_FMT "\n",ctx->bs));
394:     PetscCall(PetscViewerASCIIPrintf(viewer,"  solver method: %s\n",EPSPRIMMEMethods[(EPSPRIMMEMethod)ctx->method]));

396:     /* Display PRIMME params */
397:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank));
398:     if (!rank) primme_display_params(ctx->primme);
399:   }
400:   PetscFunctionReturn(PETSC_SUCCESS);
401: }

403: static PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps,PetscOptionItems *PetscOptionsObject)
404: {
405:   EPS_PRIMME      *ctx = (EPS_PRIMME*)eps->data;
406:   PetscInt        bs;
407:   EPSPRIMMEMethod meth;
408:   PetscBool       flg;

410:   PetscFunctionBegin;
411:   PetscOptionsHeadBegin(PetscOptionsObject,"EPS PRIMME Options");

413:     PetscCall(PetscOptionsInt("-eps_primme_blocksize","Maximum block size","EPSPRIMMESetBlockSize",ctx->bs,&bs,&flg));
414:     if (flg) PetscCall(EPSPRIMMESetBlockSize(eps,bs));

416:     PetscCall(PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg));
417:     if (flg) PetscCall(EPSPRIMMESetMethod(eps,meth));

419:   PetscOptionsHeadEnd();
420:   PetscFunctionReturn(PETSC_SUCCESS);
421: }

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

427:   PetscFunctionBegin;
428:   if (bs == PETSC_DEFAULT) ops->bs = 0;
429:   else {
430:     PetscCheck(bs>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
431:     ops->bs = bs;
432:   }
433:   PetscFunctionReturn(PETSC_SUCCESS);
434: }

436: /*@
437:    EPSPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.

439:    Logically Collective

441:    Input Parameters:
442: +  eps - the eigenproblem solver context
443: -  bs - block size

445:    Options Database Key:
446: .  -eps_primme_blocksize - Sets the max allowed block size value

448:    Notes:
449:    If the block size is not set, the value established by primme_initialize
450:    is used.

452:    The user should set the block size based on the architecture specifics
453:    of the target computer, as well as any a priori knowledge of multiplicities.
454:    The code does NOT require bs > 1 to find multiple eigenvalues. For some
455:    methods, keeping bs = 1 yields the best overall performance.

457:    Level: advanced

459: .seealso: EPSPRIMMEGetBlockSize()
460: @*/
461: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
462: {
463:   PetscFunctionBegin;
466:   PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
467:   PetscFunctionReturn(PETSC_SUCCESS);
468: }

470: static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
471: {
472:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

474:   PetscFunctionBegin;
475:   *bs = ops->bs;
476:   PetscFunctionReturn(PETSC_SUCCESS);
477: }

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

482:    Not Collective

484:    Input Parameter:
485: .  eps - the eigenproblem solver context

487:    Output Parameter:
488: .  bs - returned block size

490:    Level: advanced

492: .seealso: EPSPRIMMESetBlockSize()
493: @*/
494: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
495: {
496:   PetscFunctionBegin;
498:   PetscAssertPointer(bs,2);
499:   PetscUseMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
500:   PetscFunctionReturn(PETSC_SUCCESS);
501: }

503: static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
504: {
505:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

507:   PetscFunctionBegin;
508:   ops->method = (primme_preset_method)method;
509:   PetscFunctionReturn(PETSC_SUCCESS);
510: }

512: /*@
513:    EPSPRIMMESetMethod - Sets the method for the PRIMME library.

515:    Logically Collective

517:    Input Parameters:
518: +  eps - the eigenproblem solver context
519: -  method - method that will be used by PRIMME

521:    Options Database Key:
522: .  -eps_primme_method - Sets the method for the PRIMME library

524:    Note:
525:    If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.

527:    Level: advanced

529: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
530: @*/
531: PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
532: {
533:   PetscFunctionBegin;
536:   PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
537:   PetscFunctionReturn(PETSC_SUCCESS);
538: }

540: static PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
541: {
542:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

544:   PetscFunctionBegin;
545:   *method = (EPSPRIMMEMethod)ops->method;
546:   PetscFunctionReturn(PETSC_SUCCESS);
547: }

549: /*@
550:    EPSPRIMMEGetMethod - Gets the method for the PRIMME library.

552:    Not Collective

554:    Input Parameter:
555: .  eps - the eigenproblem solver context

557:    Output Parameter:
558: .  method - method that will be used by PRIMME

560:    Level: advanced

562: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
563: @*/
564: PetscErrorCode EPSPRIMMEGetMethod(EPS eps,EPSPRIMMEMethod *method)
565: {
566:   PetscFunctionBegin;
568:   PetscAssertPointer(method,2);
569:   PetscUseMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
570:   PetscFunctionReturn(PETSC_SUCCESS);
571: }

573: SLEPC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
574: {
575:   EPS_PRIMME     *primme;

577:   PetscFunctionBegin;
578:   PetscCall(PetscNew(&primme));
579:   eps->data = (void*)primme;

581:   primme_initialize(&primme->primme);
582:   primme->primme.globalSumReal = par_GlobalSumReal;
583: #if defined(SLEPC_HAVE_PRIMME3)
584:   primme->primme.broadcastReal = par_broadcastReal;
585: #endif
586: #if defined(SLEPC_HAVE_PRIMME2p2)
587:   primme->primme.convTestFun = convTestFun;
588:   primme->primme.monitorFun = monitorFun;
589: #endif
590:   primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;

592:   eps->categ = EPS_CATEGORY_PRECOND;

594:   eps->ops->solve          = EPSSolve_PRIMME;
595:   eps->ops->setup          = EPSSetUp_PRIMME;
596:   eps->ops->setupsort      = EPSSetUpSort_Basic;
597:   eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
598:   eps->ops->destroy        = EPSDestroy_PRIMME;
599:   eps->ops->reset          = EPSReset_PRIMME;
600:   eps->ops->view           = EPSView_PRIMME;
601:   eps->ops->backtransform  = EPSBackTransform_Default;
602:   eps->ops->setdefaultst   = EPSSetDefaultST_GMRES;

604:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME));
605:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME));
606:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME));
607:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME));
608:   PetscFunctionReturn(PETSC_SUCCESS);
609: }