Actual source code: feast.c

slepc-3.19.2 2023-09-05
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 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: {
 50:   PetscInt       ncv;
 51:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
 52:   PetscMPIInt    size;

 54:   PetscFunctionBegin;
 55:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size));
 56:   PetscCheck(size==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The FEAST interface is supported for sequential runs only");
 57:   EPSCheckHermitianDefinite(eps);
 58:   EPSCheckSinvertCayley(eps);
 59:   if (eps->ncv!=PETSC_DEFAULT) {
 60:     PetscCheck(eps->ncv>=eps->nev+2,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) PetscCall(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:   PetscCheck(eps->which==EPS_ALL && eps->inta!=eps->intb,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:   PetscCall(PetscFree4(ctx->work1,ctx->work2,ctx->Aq,ctx->Bq));
 73:   PetscCall(PetscMalloc4(eps->nloc*ncv,&ctx->work1,eps->nloc*ncv,&ctx->work2,ncv*ncv,&ctx->Aq,ncv*ncv,&ctx->Bq));

 75:   PetscCall(EPSAllocateSolution(eps,0));
 76:   PetscCall(EPSSetWorkVecs(eps,2));
 77:   PetscFunctionReturn(PETSC_SUCCESS);
 78: }

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

 95:   PetscFunctionBegin;
 96:   ncv = eps->ncv;
 97:   n   = eps->nloc;

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

112:   PetscCall(PetscMalloc1(eps->ncv,&evals));
113:   PetscCall(BVGetArray(eps->V,&pV));

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

122:   do {

124:     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);

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

172:   } while (ijob);

174:   eps->reason = EPS_CONVERGED_TOL;
175:   eps->its    = loop;
176:   eps->nconv  = nconv;
177:   if (info) {
178:     switch (info) {
179:       case 1:  /* No eigenvalue has been found in the proposed search interval */
180:         eps->nconv = 0;
181:         break;
182:       case 2:   /* FEAST did not converge "yet" */
183:         eps->reason = EPS_DIVERGED_ITS;
184:         break;
185:       default:
186:         SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by FEAST (%d)",(int)info);
187:     }
188:   }

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

192:   PetscCall(BVRestoreArray(eps->V,&pV));
193:   PetscCall(VecDestroy(&x));
194:   PetscCall(VecDestroy(&y));
195:   PetscCall(PetscFree(evals));
196:   PetscFunctionReturn(PETSC_SUCCESS);
197: }

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

203:   PetscFunctionBegin;
204:   PetscCall(PetscFree4(ctx->work1,ctx->work2,ctx->Aq,ctx->Bq));
205:   PetscFunctionReturn(PETSC_SUCCESS);
206: }

208: PetscErrorCode EPSDestroy_FEAST(EPS eps)
209: {
210:   PetscFunctionBegin;
211:   PetscCall(PetscFree(eps->data));
212:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTSetNumPoints_C",NULL));
213:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTGetNumPoints_C",NULL));
214:   PetscFunctionReturn(PETSC_SUCCESS);
215: }

217: PetscErrorCode EPSSetFromOptions_FEAST(EPS eps,PetscOptionItems *PetscOptionsObject)
218: {
219:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
220:   PetscInt       n;
221:   PetscBool      flg;

223:   PetscFunctionBegin;
224:   PetscOptionsHeadBegin(PetscOptionsObject,"EPS FEAST Options");

226:     n = ctx->npoints;
227:     PetscCall(PetscOptionsInt("-eps_feast_num_points","Number of contour integration points","EPSFEASTSetNumPoints",n,&n,&flg));
228:     if (flg) PetscCall(EPSFEASTSetNumPoints(eps,n));

230:   PetscOptionsHeadEnd();
231:   PetscFunctionReturn(PETSC_SUCCESS);
232: }

234: PetscErrorCode EPSView_FEAST(EPS eps,PetscViewer viewer)
235: {
236:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
237:   PetscBool      isascii;

239:   PetscFunctionBegin;
240:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
241:   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer,"  number of contour integration points=%" PetscInt_FMT "\n",ctx->npoints));
242:   PetscFunctionReturn(PETSC_SUCCESS);
243: }

245: PetscErrorCode EPSSetDefaultST_FEAST(EPS eps)
246: {
247:   PetscFunctionBegin;
248:   if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STSINVERT));
249:   PetscFunctionReturn(PETSC_SUCCESS);
250: }

252: static PetscErrorCode EPSFEASTSetNumPoints_FEAST(EPS eps,PetscInt npoints)
253: {
254:   EPS_FEAST *ctx = (EPS_FEAST*)eps->data;

256:   PetscFunctionBegin;
257:   if (npoints == PETSC_DEFAULT) ctx->npoints = 8;
258:   else ctx->npoints = npoints;
259:   PetscFunctionReturn(PETSC_SUCCESS);
260: }

262: /*@
263:    EPSFEASTSetNumPoints - Sets the number of contour integration points for
264:    the FEAST package.

266:    Logically Collective

268:    Input Parameters:
269: +  eps     - the eigenproblem solver context
270: -  npoints - number of contour integration points

272:    Options Database Key:
273: .  -eps_feast_num_points - Sets the number of points

275:    Level: advanced

277: .seealso: EPSFEASTGetNumPoints()
278: @*/
279: PetscErrorCode EPSFEASTSetNumPoints(EPS eps,PetscInt npoints)
280: {
281:   PetscFunctionBegin;
284:   PetscTryMethod(eps,"EPSFEASTSetNumPoints_C",(EPS,PetscInt),(eps,npoints));
285:   PetscFunctionReturn(PETSC_SUCCESS);
286: }

288: static PetscErrorCode EPSFEASTGetNumPoints_FEAST(EPS eps,PetscInt *npoints)
289: {
290:   EPS_FEAST *ctx = (EPS_FEAST*)eps->data;

292:   PetscFunctionBegin;
293:   *npoints = ctx->npoints;
294:   PetscFunctionReturn(PETSC_SUCCESS);
295: }

297: /*@
298:    EPSFEASTGetNumPoints - Gets the number of contour integration points for
299:    the FEAST package.

301:    Not Collective

303:    Input Parameter:
304: .  eps     - the eigenproblem solver context

306:    Output Parameter:
307: .  npoints - number of contour integration points

309:    Level: advanced

311: .seealso: EPSFEASTSetNumPoints()
312: @*/
313: PetscErrorCode EPSFEASTGetNumPoints(EPS eps,PetscInt *npoints)
314: {
315:   PetscFunctionBegin;
318:   PetscUseMethod(eps,"EPSFEASTGetNumPoints_C",(EPS,PetscInt*),(eps,npoints));
319:   PetscFunctionReturn(PETSC_SUCCESS);
320: }

322: SLEPC_EXTERN PetscErrorCode EPSCreate_FEAST(EPS eps)
323: {
324:   EPS_FEAST      *ctx;

326:   PetscFunctionBegin;
327:   PetscCall(PetscNew(&ctx));
328:   eps->data = (void*)ctx;

330:   eps->categ = EPS_CATEGORY_CONTOUR;

332:   eps->ops->solve          = EPSSolve_FEAST;
333:   eps->ops->setup          = EPSSetUp_FEAST;
334:   eps->ops->setupsort      = EPSSetUpSort_Basic;
335:   eps->ops->setfromoptions = EPSSetFromOptions_FEAST;
336:   eps->ops->destroy        = EPSDestroy_FEAST;
337:   eps->ops->reset          = EPSReset_FEAST;
338:   eps->ops->view           = EPSView_FEAST;
339:   eps->ops->setdefaultst   = EPSSetDefaultST_FEAST;

341:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTSetNumPoints_C",EPSFEASTSetNumPoints_FEAST));
342:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTGetNumPoints_C",EPSFEASTGetNumPoints_FEAST));
343:   PetscFunctionReturn(PETSC_SUCCESS);
344: }