Actual source code: feast.c

slepc-3.14.0 2020-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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 FEAST solver in MKL
 12: */

 14: #include <petscsys.h>
 15: #if defined(PETSC_HAVE_MKL_INTEL_ILP64)
 16: #define MKL_ILP64
 17: #endif
 18: #include <mkl.h>
 19: #include <slepc/private/epsimpl.h>

 21: #if defined(PETSC_USE_COMPLEX)
 22: #  if defined(PETSC_USE_REAL_SINGLE)
 23: #    define FEAST_RCI cfeast_hrci
 24: #    define SCALAR_CAST (MKL_Complex8*)
 25: #  else
 26: #    define FEAST_RCI zfeast_hrci
 27: #    define SCALAR_CAST (MKL_Complex16*)
 28: #  endif
 29: #else
 30: #  if defined(PETSC_USE_REAL_SINGLE)
 31: #    define FEAST_RCI sfeast_srci
 32: #  else
 33: #    define FEAST_RCI dfeast_srci
 34: #  endif
 35: #  define SCALAR_CAST
 36: #endif

 38: typedef struct {
 39:   PetscInt      npoints;          /* number of contour points */
 40:   PetscScalar   *work1,*Aq,*Bq;   /* workspace */
 41: #if defined(PETSC_USE_REAL_SINGLE)
 42:   MKL_Complex8  *work2;
 43: #else
 44:   MKL_Complex16 *work2;
 45: #endif
 46: } EPS_FEAST;

 48: PetscErrorCode EPSSetUp_FEAST(EPS eps)
 49: {
 51:   PetscInt       ncv;
 52:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
 53:   PetscMPIInt    size;

 56:   MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size);
 57:   if (size!=1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The FEAST interface is supported for sequential runs only");
 58:   EPSCheckSinvertCayley(eps);
 59:   if (eps->ncv!=PETSC_DEFAULT) {
 60:     if (eps->ncv<eps->nev+2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The value of ncv must be at least nev+2");
 61:   } else eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),eps->n); /* set default value of ncv */
 62:   if (eps->mpd!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 63:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = 20;
 64:   if (!eps->which) eps->which = EPS_ALL;
 65:   if (eps->which!=EPS_ALL || eps->inta==eps->intb) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver must be used with a computational interval");
 66:   EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_STOPPING | EPS_FEATURE_TWOSIDED);
 67:   EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION);

 69:   if (!ctx->npoints) ctx->npoints = 8;

 71:   ncv = eps->ncv;
 72:   PetscFree4(ctx->work1,ctx->work2,ctx->Aq,ctx->Bq);
 73:   PetscMalloc4(eps->nloc*ncv,&ctx->work1,eps->nloc*ncv,&ctx->work2,ncv*ncv,&ctx->Aq,ncv*ncv,&ctx->Bq);
 74:   PetscLogObjectMemory((PetscObject)eps,(2*eps->nloc*ncv+2*ncv*ncv)*sizeof(PetscScalar));

 76:   EPSAllocateSolution(eps,0);
 77:   EPSSetWorkVecs(eps,2);
 78:   return(0);
 79: }

 81: PetscErrorCode EPSSolve_FEAST(EPS eps)
 82: {
 84:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
 85:   MKL_INT        fpm[128],ijob,n,ncv,nconv,loop,info;
 86:   PetscReal      *evals,epsout=0.0;
 87:   PetscInt       i,k,nmat;
 88:   PetscScalar    *pV,*pz;
 89:   Vec            x,y,w=eps->work[0],z=eps->work[1];
 90:   Mat            A,B;
 91: #if defined(PETSC_USE_REAL_SINGLE)
 92:   MKL_Complex8   Ze;
 93: #else
 94:   MKL_Complex16  Ze;
 95: #endif

 98:   ncv = eps->ncv;
 99:   n   = eps->nloc;

101:   /* parameters */
102:   feastinit(fpm);
103:   fpm[0] = (eps->numbermonitors>0)? 1: 0;   /* runtime comments */
104:   fpm[1] = ctx->npoints;                    /* contour points */
105: #if !defined(PETSC_USE_REAL_SINGLE)
106:   fpm[2] = -PetscLog10Real(eps->tol);       /* tolerance for trace */
107: #endif
108:   fpm[3] = eps->max_it;                     /* refinement loops */
109:   fpm[5] = 1;                               /* second stopping criterion */
110: #if defined(PETSC_USE_REAL_SINGLE)
111:   fpm[6] = -PetscLog10Real(eps->tol);       /* tolerance for trace */
112: #endif

114:   PetscMalloc1(eps->ncv,&evals);
115:   BVGetArray(eps->V,&pV);

117:   ijob = -1;           /* first call to reverse communication interface */
118:   STGetNumMatrices(eps->st,&nmat);
119:   STGetMatrix(eps->st,0,&A);
120:   if (nmat>1) { STGetMatrix(eps->st,1,&B); }
121:   else B = NULL;
122:   MatCreateVecsEmpty(A,&x,&y);

124:   do {

126:     FEAST_RCI(&ijob,&n,&Ze,SCALAR_CAST ctx->work1,ctx->work2,SCALAR_CAST ctx->Aq,SCALAR_CAST ctx->Bq,fpm,&epsout,&loop,&eps->inta,&eps->intb,&ncv,evals,SCALAR_CAST pV,&nconv,eps->errest,&info);

128:     if (ncv!=eps->ncv) SETERRQ1(PetscObjectComm((PetscObject)eps),1,"FEAST changed value of ncv to %d",ncv);
129:     if (ijob == 10) {
130:       /* set new quadrature point */
131:       STSetShift(eps->st,Ze.real);
132:     } else if (ijob == 20) {
133:       /* use same quadrature point and factorization for transpose solve */
134:     } else if (ijob == 11 || ijob == 21) {
135:       /* linear solve (A-sigma*B)\work2, overwrite work2 */
136:       for (k=0;k<ncv;k++) {
137:         VecGetArray(z,&pz);
138: #if defined(PETSC_USE_COMPLEX)
139:         for (i=0;i<eps->nloc;i++) pz[i] = PetscCMPLX(ctx->work2[eps->nloc*k+i].real,ctx->work2[eps->nloc*k+i].imag);
140: #else
141:         for (i=0;i<eps->nloc;i++) pz[i] = ctx->work2[eps->nloc*k+i].real;
142: #endif
143:         VecRestoreArray(z,&pz);
144:         if (ijob == 11) {
145:           STMatSolve(eps->st,z,w);
146:         } else {
147:           VecConjugate(z);
148:           STMatSolveTranspose(eps->st,z,w);
149:           VecConjugate(w);
150:         }
151:         VecGetArray(w,&pz);
152: #if defined(PETSC_USE_COMPLEX)
153:         for (i=0;i<eps->nloc;i++) {
154:           ctx->work2[eps->nloc*k+i].real = PetscRealPart(pz[i]);
155:           ctx->work2[eps->nloc*k+i].imag = PetscImaginaryPart(pz[i]);
156:         }
157: #else
158:         for (i=0;i<eps->nloc;i++) ctx->work2[eps->nloc*k+i].real = pz[i];
159: #endif
160:         VecRestoreArray(w,&pz);
161:       }
162:     } else if (ijob == 30 || ijob == 40) {
163:       /* multiplication A*V or B*V, result in work1 */
164:       for (k=fpm[23]-1;k<fpm[23]+fpm[24]-1;k++) {
165:         VecPlaceArray(x,&pV[k*eps->nloc]);
166:         VecPlaceArray(y,&ctx->work1[k*eps->nloc]);
167:         if (ijob == 30) {
168:           MatMult(A,x,y);
169:         } else if (nmat>1) {
170:           MatMult(B,x,y);
171:         } else {
172:           VecCopy(x,y);
173:         }
174:         VecResetArray(x);
175:         VecResetArray(y);
176:       }
177:     } else if (ijob && ijob!=-2) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Internal error in FEAST reverse comunication interface (ijob=%d)",ijob);

179:   } while (ijob);

181:   eps->reason = EPS_CONVERGED_TOL;
182:   eps->its    = loop;
183:   eps->nconv  = nconv;
184:   if (info) {
185:     if (info==1) { /* No eigenvalue has been found in the proposed search interval */
186:       eps->nconv = 0;
187:     } else if (info==2) { /* FEAST did not converge "yet" */
188:       eps->reason = EPS_DIVERGED_ITS;
189:     } else SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by FEAST (%d)",info);
190:   }

192:   for (i=0;i<eps->nconv;i++) eps->eigr[i] = evals[i];

194:   BVRestoreArray(eps->V,&pV);
195:   VecDestroy(&x);
196:   VecDestroy(&y);
197:   PetscFree(evals);
198:   return(0);
199: }

201: PetscErrorCode EPSReset_FEAST(EPS eps)
202: {
204:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;

207:   PetscFree4(ctx->work1,ctx->work2,ctx->Aq,ctx->Bq);
208:   return(0);
209: }

211: PetscErrorCode EPSDestroy_FEAST(EPS eps)
212: {

216:   PetscFree(eps->data);
217:   PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTSetNumPoints_C",NULL);
218:   PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTGetNumPoints_C",NULL);
219:   return(0);
220: }

222: PetscErrorCode EPSSetFromOptions_FEAST(PetscOptionItems *PetscOptionsObject,EPS eps)
223: {
225:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
226:   PetscInt       n;
227:   PetscBool      flg;

230:   PetscOptionsHead(PetscOptionsObject,"EPS FEAST Options");

232:     n = ctx->npoints;
233:     PetscOptionsInt("-eps_feast_num_points","Number of contour integration points","EPSFEASTSetNumPoints",n,&n,&flg);
234:     if (flg) { EPSFEASTSetNumPoints(eps,n); }

236:   PetscOptionsTail();
237:   return(0);
238: }

240: PetscErrorCode EPSView_FEAST(EPS eps,PetscViewer viewer)
241: {
243:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
244:   PetscBool      isascii;

247:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
248:   if (isascii) {
249:     PetscViewerASCIIPrintf(viewer,"  number of contour integration points=%D\n",ctx->npoints);
250:   }
251:   return(0);
252: }

254: PetscErrorCode EPSSetDefaultST_FEAST(EPS eps)
255: {

259:   if (!((PetscObject)eps->st)->type_name) {
260:     STSetType(eps->st,STSINVERT);
261:   }
262:   return(0);
263: }

265: static PetscErrorCode EPSFEASTSetNumPoints_FEAST(EPS eps,PetscInt npoints)
266: {
267:   EPS_FEAST *ctx = (EPS_FEAST*)eps->data;

270:   if (npoints == PETSC_DEFAULT) ctx->npoints = 8;
271:   else ctx->npoints = npoints;
272:   return(0);
273: }

275: /*@
276:    EPSFEASTSetNumPoints - Sets the number of contour integration points for
277:    the FEAST package.

279:    Collective on EPS

281:    Input Parameters:
282: +  eps     - the eigenproblem solver context
283: -  npoints - number of contour integration points

285:    Options Database Key:
286: .  -eps_feast_num_points - Sets the number of points

288:    Level: advanced

290: .seealso: EPSFEASTGetNumPoints()
291: @*/
292: PetscErrorCode EPSFEASTSetNumPoints(EPS eps,PetscInt npoints)
293: {

299:   PetscTryMethod(eps,"EPSFEASTSetNumPoints_C",(EPS,PetscInt),(eps,npoints));
300:   return(0);
301: }

303: static PetscErrorCode EPSFEASTGetNumPoints_FEAST(EPS eps,PetscInt *npoints)
304: {
305:   EPS_FEAST *ctx = (EPS_FEAST*)eps->data;

308:   *npoints = ctx->npoints;
309:   return(0);
310: }

312: /*@
313:    EPSFEASTGetNumPoints - Gets the number of contour integration points for
314:    the FEAST package.

316:    Collective on EPS

318:    Input Parameter:
319: .  eps     - the eigenproblem solver context

321:    Output Parameter:
322: -  npoints - number of contour integration points

324:    Level: advanced

326: .seealso: EPSFEASTSetNumPoints()
327: @*/
328: PetscErrorCode EPSFEASTGetNumPoints(EPS eps,PetscInt *npoints)
329: {

335:   PetscUseMethod(eps,"EPSFEASTGetNumPoints_C",(EPS,PetscInt*),(eps,npoints));
336:   return(0);
337: }

339: SLEPC_EXTERN PetscErrorCode EPSCreate_FEAST(EPS eps)
340: {
341:   EPS_FEAST      *ctx;

345:   PetscNewLog(eps,&ctx);
346:   eps->data = (void*)ctx;

348:   eps->categ = EPS_CATEGORY_CONTOUR;

350:   eps->ops->solve          = EPSSolve_FEAST;
351:   eps->ops->setup          = EPSSetUp_FEAST;
352:   eps->ops->setupsort      = EPSSetUpSort_Basic;
353:   eps->ops->setfromoptions = EPSSetFromOptions_FEAST;
354:   eps->ops->destroy        = EPSDestroy_FEAST;
355:   eps->ops->reset          = EPSReset_FEAST;
356:   eps->ops->view           = EPSView_FEAST;
357:   eps->ops->setdefaultst   = EPSSetDefaultST_FEAST;

359:   PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTSetNumPoints_C",EPSFEASTSetNumPoints_FEAST);
360:   PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTGetNumPoints_C",EPSFEASTGetNumPoints_FEAST);
361:   return(0);
362: }