Actual source code: feast.c

slepc-3.7.1 2016-05-27
Report Typos and Errors
  1: /*
  2:    This file implements a wrapper to the FEAST package

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, 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 <../src/eps/impls/external/feast/feastp.h>

 27: PetscErrorCode EPSSolve_FEAST(EPS);

 31: PetscErrorCode EPSSetUp_FEAST(EPS eps)
 32: {
 34:   PetscInt       ncv;
 35:   PetscBool      issinv,flg;
 36:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
 37:   PetscMPIInt    size;

 40:   MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size);
 41:   if (size!=1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The FEAST interface is supported for sequential runs only");
 42:   if (eps->ncv) {
 43:     if (eps->ncv<eps->nev+2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The value of ncv must be at least nev+2");
 44:   } else eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),eps->n); /* set default value of ncv */
 45:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 46:   if (!eps->max_it) eps->max_it = PetscMax(300,(PetscInt)(2*eps->n/eps->ncv));
 47:   if (!eps->which) eps->which = EPS_ALL;

 49:   ncv = eps->ncv;
 50:   PetscFree4(ctx->work1,ctx->work2,ctx->Aq,ctx->Bq);
 51:   PetscMalloc4(eps->nloc*ncv,&ctx->work1,eps->nloc*ncv,&ctx->work2,ncv*ncv,&ctx->Aq,ncv*ncv,&ctx->Bq);
 52:   PetscLogObjectMemory((PetscObject)eps,(2*eps->nloc*ncv+2*ncv*ncv)*sizeof(PetscScalar));

 54:   if (!((PetscObject)(eps->st))->type_name) { /* default to shift-and-invert */
 55:     STSetType(eps->st,STSINVERT);
 56:   }
 57:   PetscObjectTypeCompareAny((PetscObject)eps->st,&issinv,STSINVERT,STCAYLEY,"");
 58:   if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for FEAST");

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

 62:   if (eps->which!=EPS_ALL || (eps->inta==0.0 && eps->intb==0.0)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"FEAST must be used with a computational interval");
 63:   if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"FEAST only available for symmetric/Hermitian eigenproblems");
 64:   if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Balancing not supported in the FEAST interface");
 65:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
 66:   if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");

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

 70:   EPSAllocateSolution(eps,0);
 71:   PetscObjectTypeCompare((PetscObject)eps->V,BVVECS,&flg);
 72:   if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires a BV with contiguous storage");
 73:   EPSSetWorkVecs(eps,1);

 75:   /* dispatch solve method */
 76:   eps->ops->solve = EPSSolve_FEAST;
 77:   return(0);
 78: }

 82: PetscErrorCode EPSSolve_FEAST(EPS eps)
 83: {
 85:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
 86:   PetscBLASInt   n,fpm[64],ijob,info,nev,ncv,loop;
 87:   PetscReal      *evals,epsout;
 88:   PetscInt       i,k,nmat;
 89:   PetscScalar    *pV,Ze;
 90:   Vec            v0,x,y,w = eps->work[0];
 91:   Mat            A,B;

 94:   PetscBLASIntCast(eps->nev,&nev);
 95:   PetscBLASIntCast(eps->ncv,&ncv);
 96:   PetscBLASIntCast(eps->nloc,&n);

 98:   /* parameters */
 99:   FEASTinit_(fpm);
100:   fpm[0] = (eps->numbermonitors>0)? 1: 0;                      /* runtime comments */
101:   fpm[1] = ctx->npoints;                                       /* contour points */
102:   PetscBLASIntCast(eps->max_it,&fpm[3]);  /* refinement loops */
103: #if !defined(PETSC_HAVE_MPIUNI)
104:   PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&fpm[8]);
105: #endif

107:   PetscMalloc1(eps->ncv,&evals);
108:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&x);
109:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&y);
110:   BVGetColumn(eps->V,0,&v0);
111:   VecGetArray(v0,&pV);

113:   ijob = -1;           /* first call to reverse communication interface */
114:   STGetNumMatrices(eps->st,&nmat);
115:   STGetOperators(eps->st,0,&A);
116:   if (nmat>1) { STGetOperators(eps->st,1,&B); }
117:   else B = NULL;

119:   do {

121:     PetscStackCall("FEASTrci",FEASTrci_(&ijob,&n,&Ze,ctx->work1,ctx->work2,ctx->Aq,ctx->Bq,fpm,&epsout,&loop,&eps->inta,&eps->intb,&eps->ncv,evals,pV,&eps->nconv,eps->errest,&info));

123:     if (ncv!=eps->ncv) SETERRQ1(PetscObjectComm((PetscObject)eps),1,"FEAST changed value of ncv to %d",ncv);
124:     if (ijob == 10 || ijob == 20) {
125:       /* set new quadrature point */
126:       STSetShift(eps->st,-Ze);
127:     } else if (ijob == 11 || ijob == 21) {
128:       /* linear solve (A-sigma*B)\work2, overwrite work2 */
129:       for (k=0;k<ncv;k++) {
130:         VecPlaceArray(x,ctx->work2+eps->nloc*k);
131:         if (ijob == 11) {
132:           STMatSolve(eps->st,x,w);
133:         } else {
134:           STMatSolveTranspose(eps->st,x,w);
135:         }
136:         VecCopy(w,x);
137:         VecScale(x,-1.0);
138:         VecResetArray(x);
139:       }
140:     } else if (ijob == 30 || ijob == 40) {
141:       /* multiplication A*V or B*V, result in work1 */
142:       for (k=0;k<fpm[24];k++) {
143:         VecPlaceArray(x,&pV[(fpm[23]+k-1)*eps->nloc]);
144:         VecPlaceArray(y,&ctx->work1[(fpm[23]+k-1)*eps->nloc]);
145:         if (ijob == 30) {
146:           MatMult(A,x,y);
147:         } else if (nmat>1) {
148:           MatMult(B,x,y);
149:         } else {
150:           VecCopy(x,y);
151:         }
152:         VecResetArray(x);
153:         VecResetArray(y);
154:       }
155:     } else if (ijob != 0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Internal error in FEAST reverse comunication interface (ijob=%d)",ijob);

157:   } while (ijob != 0);

159:   eps->reason = EPS_CONVERGED_TOL;
160:   eps->its = loop;
161:   if (info!=0) {
162:     if (info==1) { /* No eigenvalue has been found in the proposed search interval */
163:       eps->nconv = 0;
164:     } else if (info==2) { /* FEAST did not converge "yet" */
165:       eps->reason = EPS_DIVERGED_ITS;
166:     } else SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by FEAST (%d)",info);
167:   }

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

171:   VecRestoreArray(v0,&pV);
172:   BVRestoreColumn(eps->V,0,&v0);
173:   VecDestroy(&x);
174:   VecDestroy(&y);
175:   PetscFree(evals);
176:   return(0);
177: }

181: PetscErrorCode EPSReset_FEAST(EPS eps)
182: {
184:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;

187:   PetscFree4(ctx->work1,ctx->work2,ctx->Aq,ctx->Bq);
188:   return(0);
189: }

193: PetscErrorCode EPSDestroy_FEAST(EPS eps)
194: {

198:   PetscFree(eps->data);
199:   PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTSetNumPoints_C",NULL);
200:   PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTGetNumPoints_C",NULL);
201:   return(0);
202: }

206: PetscErrorCode EPSSetFromOptions_FEAST(PetscOptionItems *PetscOptionsObject,EPS eps)
207: {
209:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
210:   PetscInt       n;
211:   PetscBool      flg;

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

216:   n = ctx->npoints;
217:   PetscOptionsInt("-eps_feast_num_points","Number of contour integration points","EPSFEASTSetNumPoints",n,&n,&flg);
218:   if (flg) {
219:     EPSFEASTSetNumPoints(eps,n);
220:   }

222:   PetscOptionsTail();
223:   return(0);
224: }

228: PetscErrorCode EPSView_FEAST(EPS eps,PetscViewer viewer)
229: {
231:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;
232:   PetscBool      isascii;

235:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
236:   if (isascii) {
237:     PetscViewerASCIIPrintf(viewer,"  FEAST: number of contour integration points=%D\n",ctx->npoints);
238:   }
239:   return(0);
240: }

244: static PetscErrorCode EPSFEASTSetNumPoints_FEAST(EPS eps,PetscInt npoints)
245: {
247:   EPS_FEAST      *ctx = (EPS_FEAST*)eps->data;

250:   if (npoints == PETSC_DEFAULT) ctx->npoints = 8;
251:   else {
252:     PetscBLASIntCast(npoints,&ctx->npoints);
253:   }
254:   return(0);
255: }

259: /*@
260:    EPSFEASTSetNumPoints - Sets the number of contour integration points for
261:    the FEAST package.

263:    Collective on EPS

265:    Input Parameters:
266: +  eps     - the eigenproblem solver context
267: -  npoints - number of contour integration points

269:    Options Database Key:
270: .  -eps_feast_num_points - Sets the number of points

272:    Level: advanced

274: .seealso: EPSFEASTGetNumPoints()
275: @*/
276: PetscErrorCode EPSFEASTSetNumPoints(EPS eps,PetscInt npoints)
277: {

283:   PetscTryMethod(eps,"EPSFEASTSetNumPoints_C",(EPS,PetscInt),(eps,npoints));
284:   return(0);
285: }

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

294:   *npoints = ctx->npoints;
295:   return(0);
296: }

300: /*@
301:    EPSFEASTGetNumPoints - Gets the number of contour integration points for
302:    the FEAST package.

304:    Collective on EPS

306:    Input Parameter:
307: .  eps     - the eigenproblem solver context

309:    Output Parameter:
310: -  npoints - number of contour integration points

312:    Level: advanced

314: .seealso: EPSFEASTSetNumPoints()
315: @*/
316: PetscErrorCode EPSFEASTGetNumPoints(EPS eps,PetscInt *npoints)
317: {

323:   PetscUseMethod(eps,"EPSFEASTGetNumPoints_C",(EPS,PetscInt*),(eps,npoints));
324:   return(0);
325: }

329: PETSC_EXTERN PetscErrorCode EPSCreate_FEAST(EPS eps)
330: {
331:   EPS_FEAST      *ctx;

335:   PetscNewLog(eps,&ctx);
336:   eps->data = (void*)ctx;

338:   eps->ops->setup                = EPSSetUp_FEAST;
339:   eps->ops->setfromoptions       = EPSSetFromOptions_FEAST;
340:   eps->ops->destroy              = EPSDestroy_FEAST;
341:   eps->ops->reset                = EPSReset_FEAST;
342:   eps->ops->view                 = EPSView_FEAST;
343:   PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTSetNumPoints_C",EPSFEASTSetNumPoints_FEAST);
344:   PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTGetNumPoints_C",EPSFEASTGetNumPoints_FEAST);
345:   return(0);
346: }