Actual source code: primme.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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>    /*I "slepceps.h" I*/

 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: typedef struct {
 33:   primme_params primme;           /* param struc */
 34:   primme_preset_method method;    /* primme method */
 35:   Mat       A;                    /* problem matrix */
 36:   KSP       ksp;                  /* linear solver and preconditioner */
 37:   Vec       x,y;                  /* auxiliary vectors */
 38:   double    target;               /* a copy of eps's target */
 39: } EPS_PRIMME;

 41: static void multMatvec_PRIMME(void*,PRIMME_INT*,void*,PRIMME_INT*,int*,struct primme_params*,int*);
 42: static void applyPreconditioner_PRIMME(void*,PRIMME_INT*,void*,PRIMME_INT*,int*,struct primme_params*,int*);

 44: static void par_GlobalSumReal(void *sendBuf,void *recvBuf,int *count,primme_params *primme,int *ierr)
 45: {
 46:   *MPI_Allreduce(sendBuf,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));CHKERRABORT(PetscObjectComm((PetscObject)primme->commInfo),*ierr);
 47: }

 49: PetscErrorCode EPSSetUp_PRIMME(EPS eps)
 50: {
 52:   PetscMPIInt    numProcs,procID;
 53:   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;
 54:   primme_params  *primme = &ops->primme;
 55:   PetscBool      istrivial,flg;

 58:   MPI_Comm_size(PetscObjectComm((PetscObject)eps),&numProcs);
 59:   MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&procID);

 61:   /* Check some constraints and set some default values */
 62:   if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
 63:   STGetMatrix(eps->st,0,&ops->A);
 64:   if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME is only available for Hermitian problems");
 65:   if (eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME is not available for generalized problems");
 66:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
 67:   if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");
 68:   if (!eps->which) eps->which = EPS_LARGEST_REAL;
 69:   if (eps->converged != EPSConvergedAbsolute) { PetscInfo(eps,"Warning: using absolute convergence test\n"); }
 70:   RGIsTrivial(eps->rg,&istrivial);
 71:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");

 73:   /* Transfer SLEPc options to PRIMME options */
 74:   primme->n          = eps->n;
 75:   primme->nLocal     = eps->nloc;
 76:   primme->numEvals   = eps->nev;
 77:   primme->matrix     = ops;
 78:   primme->commInfo   = eps;
 79:   primme->maxMatvecs = eps->max_it;
 80:   primme->eps        = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
 81:   primme->numProcs   = numProcs;
 82:   primme->procID     = procID;
 83:   primme->printLevel = 0;
 84:   primme->correctionParams.precondition = 1;

 86:   switch (eps->which) {
 87:     case EPS_LARGEST_REAL:
 88:       primme->target = primme_largest;
 89:       break;
 90:     case EPS_SMALLEST_REAL:
 91:       primme->target = primme_smallest;
 92:       break;
 93:     case EPS_TARGET_MAGNITUDE:
 94:     case EPS_TARGET_REAL:
 95:       primme->target = primme_closest_abs;
 96:       primme->numTargetShifts = 1;
 97:       ops->target = PetscRealPart(eps->target);
 98:       primme->targetShifts = &ops->target;
 99:       break;
100:     default:
101:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'which' value not supported by PRIMME");
102:       break;
103:   }

105:   if (primme_set_method(ops->method,primme) < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME method not valid");

107:   /* If user sets ncv, maxBasisSize is modified. If not, ncv is set as maxBasisSize */
108:   if (eps->ncv) primme->maxBasisSize = eps->ncv;
109:   else eps->ncv = primme->maxBasisSize;
110:   if (eps->ncv < eps->nev+primme->maxBlockSize) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME needs ncv >= nev+maxBlockSize");
111:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }

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

115:   /* Set workspace */
116:   EPSAllocateSolution(eps,0);
117:   PetscObjectTypeCompare((PetscObject)eps->V,BVVECS,&flg);
118:   if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires a BV with contiguous storage");

120:   /* Setup the preconditioner */
121:   if (primme->correctionParams.precondition) {
122:     STGetKSP(eps->st,&ops->ksp);
123:     PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&flg);
124:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME only works with KSPPREONLY");
125:     primme->preconditioner = NULL;
126:     primme->applyPreconditioner = applyPreconditioner_PRIMME;
127:   }

129:   /* Prepare auxiliary vectors */
130:   if (!ops->x) {
131:     MatCreateVecsEmpty(ops->A,&ops->x,&ops->y);
132:     PetscLogObjectParent((PetscObject)eps,(PetscObject)ops->x);
133:     PetscLogObjectParent((PetscObject)eps,(PetscObject)ops->y);
134:   }
135:   return(0);
136: }

138: PetscErrorCode EPSSolve_PRIMME(EPS eps)
139: {
141:   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;
142:   PetscScalar    *a;
143:   Vec            v0;
144: #if defined(PETSC_USE_COMPLEX)
145:   PetscInt       i;
146:   PetscReal      *evals;
147: #endif

150:   /* Reset some parameters left from previous runs */
151:   ops->primme.aNorm    = 1.0;
152:   ops->primme.initSize = eps->nini;
153:   ops->primme.iseed[0] = -1;

155:   /* Call PRIMME solver */
156:   BVGetColumn(eps->V,0,&v0);
157:   VecGetArray(v0,&a);
158: #if !defined(PETSC_USE_COMPLEX)
159:   PRIMME_DRIVER(eps->eigr,a,eps->errest,&ops->primme);
160:   if (ierr) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d",ierr);
161: #else
162:   /* PRIMME returns real eigenvalues, but SLEPc works with complex ones */
163:   PetscMalloc1(eps->ncv,&evals);
164:   PRIMME_DRIVER(evals,a,eps->errest,&ops->primme);
165:   if (ierr) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d",ierr);
166:   for (i=0;i<eps->ncv;i++) eps->eigr[i] = evals[i];
167:   PetscFree(evals);
168: #endif
169:   VecRestoreArray(v0,&a);
170:   BVRestoreColumn(eps->V,0,&v0);

172:   eps->nconv  = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
173:   eps->reason = eps->ncv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
174:   eps->its    = ops->primme.stats.numOuterIterations;
175:   return(0);
176: }

178: static void multMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
179: {
180:   PetscInt   i;
181:   EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
182:   Vec        x = ops->x,y = ops->y;
183:   Mat        A = ops->A;

186:   for (i=0;i<*blockSize;i++) {
187:     *VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
188:     *VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
189:     *MatMult(A,x,y);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
190:     *VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
191:     *VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)A),*ierr);
192:   }
193:   PetscFunctionReturnVoid();
194: }

196: static void applyPreconditioner_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
197: {
198:   PetscInt   i;
199:   EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
200:   Vec        x = ops->x,y = ops->y;

203:   for (i=0;i<*blockSize;i++) {
204:     *VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
205:     *VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
206:     *KSPSolve(ops->ksp,x,y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
207:     *VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
208:     *VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),*ierr);
209:   }
210:   PetscFunctionReturnVoid();
211: }

213: PetscErrorCode EPSReset_PRIMME(EPS eps)
214: {
216:   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;

219:   primme_free(&ops->primme);
220:   VecDestroy(&ops->x);
221:   VecDestroy(&ops->y);
222:   return(0);
223: }

225: PetscErrorCode EPSDestroy_PRIMME(EPS eps)
226: {

230:   PetscFree(eps->data);
231:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL);
232:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL);
233:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL);
234:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL);
235:   return(0);
236: }

238: PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
239: {
240:   PetscErrorCode  ierr;
241:   PetscBool       isascii;
242:   primme_params   *primme = &((EPS_PRIMME*)eps->data)->primme;
243:   EPSPRIMMEMethod methodn;
244:   PetscMPIInt     rank;

247:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
248:   if (isascii) {
249:     PetscViewerASCIIPrintf(viewer,"  block size=%D\n",primme->maxBlockSize);
250:     EPSPRIMMEGetMethod(eps,&methodn);
251:     PetscViewerASCIIPrintf(viewer,"  solver method: %s\n",EPSPRIMMEMethods[methodn]);

253:     /* Display PRIMME params */
254:     MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank);
255:     if (!rank) primme_display_params(*primme);
256:   }
257:   return(0);
258: }

260: PetscErrorCode EPSSetFromOptions_PRIMME(PetscOptionItems *PetscOptionsObject,EPS eps)
261: {
262:   PetscErrorCode  ierr;
263:   EPS_PRIMME      *ctx = (EPS_PRIMME*)eps->data;
264:   PetscInt        bs;
265:   EPSPRIMMEMethod meth;
266:   PetscBool       flg;

269:   PetscOptionsHead(PetscOptionsObject,"EPS PRIMME Options");

271:     PetscOptionsInt("-eps_primme_blocksize","Maximum block size","EPSPRIMMESetBlockSize",ctx->primme.maxBlockSize,&bs,&flg);
272:     if (flg) { EPSPRIMMESetBlockSize(eps,bs); }

274:     PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
275:     if (flg) { EPSPRIMMESetMethod(eps,meth); }

277:   PetscOptionsTail();
278:   return(0);
279: }

281: static PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
282: {
283:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

286:   if (bs == PETSC_DEFAULT) ops->primme.maxBlockSize = 1;
287:   else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
288:   else ops->primme.maxBlockSize = bs;
289:   return(0);
290: }

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

295:    Logically Collective on EPS

297:    Input Parameters:
298: +  eps - the eigenproblem solver context
299: -  bs - block size

301:    Options Database Key:
302: .  -eps_primme_blocksize - Sets the max allowed block size value

304:    Notes:
305:    If the block size is not set, the value established by primme_initialize
306:    is used.

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

313:    Level: advanced

315: .seealso: EPSPRIMMEGetBlockSize()
316: @*/
317: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
318: {

324:   PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
325:   return(0);
326: }

328: static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
329: {
330:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

333:   *bs = ops->primme.maxBlockSize;
334:   return(0);
335: }

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

340:    Not Collective

342:    Input Parameter:
343: .  eps - the eigenproblem solver context

345:    Output Parameter:
346: .  bs - returned block size

348:    Level: advanced

350: .seealso: EPSPRIMMESetBlockSize()
351: @*/
352: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
353: {

359:   PetscUseMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
360:   return(0);
361: }

363: static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
364: {
365:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

368:   ops->method = (primme_preset_method)method;
369:   return(0);
370: }

372: /*@
373:    EPSPRIMMESetMethod - Sets the method for the PRIMME library.

375:    Logically Collective on EPS

377:    Input Parameters:
378: +  eps - the eigenproblem solver context
379: -  method - method that will be used by PRIMME

381:    Options Database Key:
382: .  -eps_primme_method - Sets the method for the PRIMME library

384:    Note:
385:    If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.

387:    Level: advanced

389: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
390: @*/
391: PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
392: {

398:   PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
399:   return(0);
400: }

402: static PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
403: {
404:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

407:   *method = (EPSPRIMMEMethod)ops->method;
408:   return(0);
409: }

411: /*@
412:    EPSPRIMMEGetMethod - Gets the method for the PRIMME library.

414:    Not Collective

416:    Input Parameter:
417: .  eps - the eigenproblem solver context

419:    Output Parameter:
420: .  method - method that will be used by PRIMME

422:    Level: advanced

424: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
425: @*/
426: PetscErrorCode EPSPRIMMEGetMethod(EPS eps,EPSPRIMMEMethod *method)
427: {

433:   PetscUseMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
434:   return(0);
435: }

437: SLEPC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
438: {
440:   EPS_PRIMME     *primme;

443:   PetscNewLog(eps,&primme);
444:   eps->data = (void*)primme;

446:   primme_initialize(&primme->primme);
447:   primme->primme.matrixMatvec = multMatvec_PRIMME;
448:   primme->primme.globalSumReal = par_GlobalSumReal;
449:   primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;

451:   eps->categ = EPS_CATEGORY_PRECOND;

453:   eps->ops->solve          = EPSSolve_PRIMME;
454:   eps->ops->setup          = EPSSetUp_PRIMME;
455:   eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
456:   eps->ops->destroy        = EPSDestroy_PRIMME;
457:   eps->ops->reset          = EPSReset_PRIMME;
458:   eps->ops->view           = EPSView_PRIMME;
459:   eps->ops->backtransform  = EPSBackTransform_Default;
460:   eps->ops->setdefaultst   = EPSSetDefaultST_Precond;

462:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME);
463:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME);
464:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME);
465:   PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME);
466:   return(0);
467: }