Actual source code: evsl.c

slepc-3.21.0 2024-03-30
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 eigensolvers in EVSL.
 12: */

 14: #include <slepc/private/epsimpl.h>
 15: #include <evsl.h>

 17: #define PetscCallEVSL(func, ...) do {                                                   \
 18:     PetscStackPushExternal(PetscStringize(func));                                                      \
 19:     int evsl_ierr_ = func(__VA_ARGS__);                                              \
 20:     PetscStackPop;                                                                             \
 21:     PetscCheck(!evsl_ierr_,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error calling %s: error code %d",PetscStringize(func(__VA_ARGS__)),evsl_ierr_); \
 22:   } while (0)

 24: typedef struct {
 25:   PetscBool         initialized;
 26:   Mat               A;           /* problem matrix */
 27:   Vec               x,y;         /* auxiliary vectors */
 28:   PetscReal         *sli;        /* slice bounds */
 29:   PetscInt          nev;         /* approximate number of wanted eigenvalues in each slice */
 30:   PetscLayout       map;         /* used to distribute slices among MPI processes */
 31:   PetscBool         estimrange;  /* the filter range was not set by the user */
 32:   /* user parameters */
 33:   PetscInt          nslices;     /* number of slices */
 34:   PetscReal         lmin,lmax;   /* numerical range (min and max eigenvalue) */
 35:   EPSEVSLDOSMethod  dos;         /* DOS method, either KPM or Lanczos */
 36:   PetscInt          nvec;        /* number of sample vectors used for DOS */
 37:   PetscInt          deg;         /* polynomial degree used for DOS (KPM only) */
 38:   PetscInt          steps;       /* number of Lanczos steps used for DOS (Lanczos only) */
 39:   PetscInt          npoints;     /* number of sample points used for DOS (Lanczos only) */
 40:   PetscInt          max_deg;     /* maximum degree allowed for the polynomial */
 41:   PetscReal         thresh;      /* threshold for accepting polynomial */
 42:   EPSEVSLDamping    damping;     /* type of damping (for polynomial and for DOS-KPM) */
 43: } EPS_EVSL;

 45: static void AMatvec_EVSL(double *xa,double *ya,void *data)
 46: {
 47:   EPS_EVSL       *ctx = (EPS_EVSL*)data;
 48:   Vec            x = ctx->x,y = ctx->y;
 49:   Mat            A = ctx->A;

 51:   PetscFunctionBegin;
 52:   PetscCallAbort(PetscObjectComm((PetscObject)A),VecPlaceArray(x,(PetscScalar*)xa));
 53:   PetscCallAbort(PetscObjectComm((PetscObject)A),VecPlaceArray(y,(PetscScalar*)ya));
 54:   PetscCallAbort(PetscObjectComm((PetscObject)A),MatMult(A,x,y));
 55:   PetscCallAbort(PetscObjectComm((PetscObject)A),VecResetArray(x));
 56:   PetscCallAbort(PetscObjectComm((PetscObject)A),VecResetArray(y));
 57:   PetscFunctionReturnVoid();
 58: }

 60: static PetscErrorCode EPSSetUp_EVSL(EPS eps)
 61: {
 62:   EPS_EVSL       *ctx = (EPS_EVSL*)eps->data;
 63:   PetscMPIInt    size,rank;
 64:   PetscBool      isshift;
 65:   PetscScalar    *vinit;
 66:   PetscReal      *mu,ecount,xintv[4],*xdos,*ydos;
 67:   Vec            v0;
 68:   Mat            A;
 69:   PetscRandom    rnd;

 71:   PetscFunctionBegin;
 72:   EPSCheckStandard(eps);
 73:   EPSCheckHermitian(eps);
 74:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift));
 75:   PetscCheck(isshift,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations");

 77:   if (ctx->initialized) EVSLFinish();
 78:   EVSLStart();
 79:   ctx->initialized=PETSC_TRUE;

 81:   /* get number of slices per process */
 82:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size));
 83:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank));
 84:   if (!ctx->nslices) ctx->nslices = size;
 85:   PetscCall(PetscLayoutDestroy(&ctx->map));
 86:   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)eps),PETSC_DECIDE,ctx->nslices,1,&ctx->map));

 88:   /* get matrix and prepare auxiliary vectors */
 89:   PetscCall(MatDestroy(&ctx->A));
 90:   PetscCall(STGetMatrix(eps->st,0,&A));
 91:   if (size==1) {
 92:     PetscCall(PetscObjectReference((PetscObject)A));
 93:     ctx->A = A;
 94:   } else PetscCall(MatCreateRedundantMatrix(A,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&ctx->A));
 95:   SetAMatvec(eps->n,&AMatvec_EVSL,(void*)ctx);
 96:   if (!ctx->x) PetscCall(MatCreateVecsEmpty(ctx->A,&ctx->x,&ctx->y));
 97:   EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
 98:   EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE);

100:   if (!eps->which) eps->which=EPS_ALL;
101:   PetscCheck(eps->which==EPS_ALL && eps->inta!=eps->intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires setting an interval with EPSSetInterval()");

103:   /* estimate numerical range */
104:   if (ctx->estimrange || ctx->lmin == PETSC_MIN_REAL || ctx->lmax == PETSC_MAX_REAL) {
105:     PetscCall(MatCreateVecs(ctx->A,&v0,NULL));
106:     if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
107:     PetscCall(BVGetRandomContext(eps->V,&rnd));
108:     PetscCall(VecSetRandom(v0,rnd));
109:     PetscCall(VecGetArray(v0,&vinit));
110:     PetscCallEVSL(LanTrbounds,50,200,eps->tol,vinit,1,&ctx->lmin,&ctx->lmax,NULL);
111:     PetscCall(VecRestoreArray(v0,&vinit));
112:     PetscCall(VecDestroy(&v0));
113:     ctx->estimrange = PETSC_TRUE;   /* estimate if called again with another matrix */
114:   }
115:   PetscCheck(ctx->lmin<=eps->inta && ctx->lmax>=eps->intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The requested interval [%g,%g] must be contained in the numerical range [%g,%g]",(double)eps->inta,(double)eps->intb,(double)ctx->lmin,(double)ctx->lmax);
116:   xintv[0] = eps->inta;
117:   xintv[1] = eps->intb;
118:   xintv[2] = ctx->lmin;
119:   xintv[3] = ctx->lmax;

121:   /* estimate number of eigenvalues in the interval */
122:   switch (ctx->dos) {
123:     case EPS_EVSL_DOS_KPM:
124:       PetscCall(PetscMalloc1(ctx->deg+1,&mu));
125:       if (!rank) PetscCallEVSL(kpmdos,ctx->deg,(int)ctx->damping,ctx->nvec,xintv,mu,&ecount);
126:       PetscCallMPI(MPI_Bcast(mu,ctx->deg+1,MPIU_REAL,0,PetscObjectComm((PetscObject)eps)));
127:       break;
128:     case EPS_EVSL_DOS_LANCZOS:
129:       PetscCall(PetscMalloc2(ctx->npoints,&xdos,ctx->npoints,&ydos));
130:       if (!rank) PetscCallEVSL(LanDos,ctx->nvec,PetscMin(ctx->steps,eps->n/2),ctx->npoints,xdos,ydos,&ecount,xintv);
131:       PetscCallMPI(MPI_Bcast(xdos,ctx->npoints,MPIU_REAL,0,PetscObjectComm((PetscObject)eps)));
132:       PetscCallMPI(MPI_Bcast(ydos,ctx->npoints,MPIU_REAL,0,PetscObjectComm((PetscObject)eps)));
133:       break;
134:     default:
135:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid DOS method");
136:   }
137:   PetscCallMPI(MPI_Bcast(&ecount,1,MPIU_REAL,0,PetscObjectComm((PetscObject)eps)));

139:   PetscCall(PetscInfo(eps,"Estimated eigenvalue count in the interval: %g\n",ecount));
140:   eps->ncv = (PetscInt)PetscCeilReal(1.5*ecount);

142:   /* slice the spectrum */
143:   PetscCall(PetscFree(ctx->sli));
144:   PetscCall(PetscMalloc1(ctx->nslices+1,&ctx->sli));
145:   if (ctx->dos == EPS_EVSL_DOS_KPM) {
146:     PetscCallEVSL(spslicer,ctx->sli,mu,ctx->deg,xintv,ctx->nslices,10*(PetscInt)ecount);
147:     PetscCall(PetscFree(mu));
148:   } else if (ctx->dos == EPS_EVSL_DOS_LANCZOS) {
149:     spslicer2(xdos,ydos,ctx->nslices,ctx->npoints,ctx->sli);
150:     PetscCall(PetscFree2(xdos,ydos));
151:   }

153:   /* approximate number of eigenvalues wanted in each slice */
154:   ctx->nev = (PetscInt)(1.0 + ecount/(PetscReal)ctx->nslices) + 2;

156:   if (eps->mpd!=PETSC_DEFAULT) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
157:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
158:   PetscCall(EPSAllocateSolution(eps,0));
159:   PetscFunctionReturn(PETSC_SUCCESS);
160: }

162: static PetscErrorCode EPSSolve_EVSL(EPS eps)
163: {
164:   EPS_EVSL       *ctx = (EPS_EVSL*)eps->data;
165:   PetscInt       i,j,k=0,sl,mlan,nevout,*ind,nevmax,rstart,rend,*nevloc,*disp,N;
166:   PetscReal      *res,xintv[4],*errest;
167:   PetscScalar    *lam,*X,*Y,*vinit,*eigr;
168:   PetscMPIInt    size,rank;
169:   PetscRandom    rnd;
170:   Vec            v,w,v0,x;
171:   VecScatter     vs;
172:   IS             is;
173:   polparams      pol;

175:   PetscFunctionBegin;
176:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size));
177:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank));
178:   PetscCall(PetscLayoutGetRange(ctx->map,&rstart,&rend));
179:   nevmax = (rend-rstart)*ctx->nev;
180:   PetscCall(MatCreateVecs(ctx->A,&v0,NULL));
181:   PetscCall(BVGetRandomContext(eps->V,&rnd));
182:   PetscCall(VecSetRandom(v0,rnd));
183:   PetscCall(VecGetArray(v0,&vinit));
184:   PetscCall(PetscMalloc5(size,&nevloc,size+1,&disp,nevmax,&eigr,nevmax,&errest,nevmax*eps->n,&X));
185:   mlan = PetscMin(PetscMax(5*ctx->nev,300),eps->n);
186:   for (sl=rstart; sl<rend; sl++) {
187:     xintv[0] = ctx->sli[sl];
188:     xintv[1] = ctx->sli[sl+1];
189:     xintv[2] = ctx->lmin;
190:     xintv[3] = ctx->lmax;
191:     PetscCall(PetscInfo(ctx->A,"Subinterval %" PetscInt_FMT ": [%.4e, %.4e]\n",sl+1,xintv[0],xintv[1]));
192:     set_pol_def(&pol);
193:     pol.max_deg    = ctx->max_deg;
194:     pol.damping    = (int)ctx->damping;
195:     pol.thresh_int = ctx->thresh;
196:     find_pol(xintv,&pol);
197:     PetscCall(PetscInfo(ctx->A,"Polynomial [type = %" PetscInt_FMT "], deg %" PetscInt_FMT ", bar %e gam %e\n",pol.type,pol.deg,pol.bar,pol.gam));
198:     PetscCallEVSL(ChebLanNr,xintv,mlan,eps->tol,vinit,&pol,&nevout,&lam,&Y,&res,NULL);
199:     PetscCheck(k+nevout<=nevmax,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Too low estimation of eigenvalue count, try modifying the sampling parameters");
200:     free_pol(&pol);
201:     PetscCall(PetscInfo(ctx->A,"Computed %" PetscInt_FMT " eigenvalues\n",nevout));
202:     PetscCall(PetscMalloc1(nevout,&ind));
203:     sort_double(nevout,lam,ind);
204:     for (i=0;i<nevout;i++) {
205:       eigr[i+k]   = lam[i];
206:       errest[i+k] = res[ind[i]];
207:       PetscCall(PetscArraycpy(X+(i+k)*eps->n,Y+ind[i]*eps->n,eps->n));
208:     }
209:     k += nevout;
210:     if (lam) evsl_Free(lam);
211:     if (Y)   evsl_Free_device(Y);
212:     if (res) evsl_Free(res);
213:     PetscCall(PetscFree(ind));
214:   }
215:   PetscCall(VecRestoreArray(v0,&vinit));
216:   PetscCall(VecDestroy(&v0));

218:   /* gather eigenvalues computed by each MPI process */
219:   PetscCallMPI(MPI_Allgather(&k,1,MPIU_INT,nevloc,1,MPIU_INT,PetscObjectComm((PetscObject)eps)));
220:   eps->nev = nevloc[0];
221:   disp[0]  = 0;
222:   for (i=1;i<size;i++) {
223:     eps->nev += nevloc[i];
224:     disp[i]   = disp[i-1]+nevloc[i-1];
225:   }
226:   disp[size] = disp[size-1]+nevloc[size-1];
227:   PetscCheck(eps->nev<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Too low estimation of eigenvalue count, try modifying the sampling parameters");
228:   PetscCallMPI(MPI_Allgatherv(eigr,k,MPIU_SCALAR,eps->eigr,nevloc,disp,MPIU_SCALAR,PetscObjectComm((PetscObject)eps)));
229:   PetscCallMPI(MPI_Allgatherv(errest,k,MPIU_REAL,eps->errest,nevloc,disp,MPIU_REAL,PetscObjectComm((PetscObject)eps)));
230:   eps->nconv  = eps->nev;
231:   eps->its    = 1;
232:   eps->reason = EPS_CONVERGED_TOL;

234:   /* scatter computed eigenvectors and store them in eps->V */
235:   PetscCall(BVCreateVec(eps->V,&w));
236:   for (i=0;i<size;i++) {
237:     N = (rank==i)? eps->n: 0;
238:     PetscCall(VecCreateSeq(PETSC_COMM_SELF,N,&x));
239:     PetscCall(VecSetFromOptions(x));
240:     PetscCall(ISCreateStride(PETSC_COMM_SELF,N,0,1,&is));
241:     PetscCall(VecScatterCreate(x,is,w,is,&vs));
242:     PetscCall(ISDestroy(&is));
243:     for (j=disp[i];j<disp[i+1];j++) {
244:       PetscCall(BVGetColumn(eps->V,j,&v));
245:       if (rank==i) PetscCall(VecPlaceArray(x,X+(j-disp[i])*eps->n));
246:       PetscCall(VecScatterBegin(vs,x,v,INSERT_VALUES,SCATTER_FORWARD));
247:       PetscCall(VecScatterEnd(vs,x,v,INSERT_VALUES,SCATTER_FORWARD));
248:       if (rank==i) PetscCall(VecResetArray(x));
249:       PetscCall(BVRestoreColumn(eps->V,j,&v));
250:     }
251:     PetscCall(VecScatterDestroy(&vs));
252:     PetscCall(VecDestroy(&x));
253:   }
254:   PetscCall(VecDestroy(&w));
255:   PetscCall(PetscFree5(nevloc,disp,eigr,errest,X));
256:   PetscFunctionReturn(PETSC_SUCCESS);
257: }

259: static PetscErrorCode EPSEVSLSetSlices_EVSL(EPS eps,PetscInt nslices)
260: {
261:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

263:   PetscFunctionBegin;
264:   if (nslices == PETSC_DECIDE || nslices == PETSC_DEFAULT) nslices = 0;
265:   else PetscCheck(nslices>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Number of slices must be 1 at least");
266:   if (ctx->nslices != nslices) {
267:     ctx->nslices = nslices;
268:     eps->state   = EPS_STATE_INITIAL;
269:   }
270:   PetscFunctionReturn(PETSC_SUCCESS);
271: }

273: /*@
274:    EPSEVSLSetSlices - Set the number of slices in which the interval must be
275:    subdivided.

277:    Logically Collective

279:    Input Parameters:
280: +  eps     - the eigensolver context
281: -  nslices - the number of slices

283:    Options Database Key:
284: .  -eps_evsl_slices <n> - set the number of slices to n

286:    Notes:
287:    By default, one slice per MPI process is used. Depending on the number of
288:    eigenvalues, using more slices may be beneficial, but very narrow subintervals
289:    imply higher polynomial degree.

291:    Level: intermediate

293: .seealso: EPSEVSLGetSlices()
294: @*/
295: PetscErrorCode EPSEVSLSetSlices(EPS eps,PetscInt nslices)
296: {
297:   PetscFunctionBegin;
300:   PetscTryMethod(eps,"EPSEVSLSetSlices_C",(EPS,PetscInt),(eps,nslices));
301:   PetscFunctionReturn(PETSC_SUCCESS);
302: }

304: static PetscErrorCode EPSEVSLGetSlices_EVSL(EPS eps,PetscInt *nslices)
305: {
306:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

308:   PetscFunctionBegin;
309:   *nslices = ctx->nslices;
310:   PetscFunctionReturn(PETSC_SUCCESS);
311: }

313: /*@
314:    EPSEVSLGetSlices - Gets the number of slices in which the interval must be
315:    subdivided.

317:    Not Collective

319:    Input Parameter:
320: .  eps - the eigensolver context

322:    Output Parameter:
323: .  nslices - the number of slices

325:    Level: intermediate

327: .seealso: EPSEVSLSetSlices()
328: @*/
329: PetscErrorCode EPSEVSLGetSlices(EPS eps,PetscInt *nslices)
330: {
331:   PetscFunctionBegin;
333:   PetscAssertPointer(nslices,2);
334:   PetscUseMethod(eps,"EPSEVSLGetSlices_C",(EPS,PetscInt*),(eps,nslices));
335:   PetscFunctionReturn(PETSC_SUCCESS);
336: }

338: static PetscErrorCode EPSEVSLSetRange_EVSL(EPS eps,PetscReal lmin,PetscReal lmax)
339: {
340:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

342:   PetscFunctionBegin;
343:   PetscCheck(lmin<lmax,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be lmin<lmax");
344:   if (ctx->lmin != lmin || ctx->lmax != lmax) {
345:     ctx->lmin  = lmin;
346:     ctx->lmax  = lmax;
347:     eps->state = EPS_STATE_INITIAL;
348:   }
349:   PetscFunctionReturn(PETSC_SUCCESS);
350: }

352: /*@
353:    EPSEVSLSetRange - Defines the numerical range (or field of values) of the problem,
354:    that is, the interval containing all eigenvalues.

356:    Logically Collective

358:    Input Parameters:
359: +  eps  - the eigensolver context
360: .  lmin - left end of the interval
361: -  lmax - right end of the interval

363:    Options Database Key:
364: .  -eps_evsl_range <a,b> - set [a,b] as the numerical range

366:    Notes:
367:    The filter will be most effective if the numerical range is tight, that is, lmin
368:    and lmax are good approximations to the leftmost and rightmost eigenvalues,
369:    respectively. If not set by the user, an approximation is computed internally.

371:    The wanted computational interval specified via EPSSetInterval() must be
372:    contained in the numerical range.

374:    Level: intermediate

376: .seealso: EPSEVSLGetRange(), EPSSetInterval()
377: @*/
378: PetscErrorCode EPSEVSLSetRange(EPS eps,PetscReal lmin,PetscReal lmax)
379: {
380:   PetscFunctionBegin;
384:   PetscTryMethod(eps,"EPSEVSLSetRange_C",(EPS,PetscReal,PetscReal),(eps,lmin,lmax));
385:   PetscFunctionReturn(PETSC_SUCCESS);
386: }

388: static PetscErrorCode EPSEVSLGetRange_EVSL(EPS eps,PetscReal *lmin,PetscReal *lmax)
389: {
390:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

392:   PetscFunctionBegin;
393:   if (lmin) *lmin = ctx->lmin;
394:   if (lmax) *lmax = ctx->lmax;
395:   PetscFunctionReturn(PETSC_SUCCESS);
396: }

398: /*@
399:    EPSEVSLGetRange - Gets the interval containing all eigenvalues.

401:    Not Collective

403:    Input Parameter:
404: .  eps - the eigensolver context

406:    Output Parameters:
407: +  lmin - left end of the interval
408: -  lmax - right end of the interval

410:    Level: intermediate

412: .seealso: EPSEVSLSetRange()
413: @*/
414: PetscErrorCode EPSEVSLGetRange(EPS eps,PetscReal *lmin,PetscReal *lmax)
415: {
416:   PetscFunctionBegin;
418:   PetscUseMethod(eps,"EPSEVSLGetRange_C",(EPS,PetscReal*,PetscReal*),(eps,lmin,lmax));
419:   PetscFunctionReturn(PETSC_SUCCESS);
420: }

422: static PetscErrorCode EPSEVSLSetDOSParameters_EVSL(EPS eps,EPSEVSLDOSMethod dos,PetscInt nvec,PetscInt deg,PetscInt steps,PetscInt npoints)
423: {
424:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

426:   PetscFunctionBegin;
427:   ctx->dos = dos;
428:   if (nvec == PETSC_DECIDE || nvec == PETSC_DEFAULT) ctx->nvec = 80;
429:   else {
430:     PetscCheck(nvec>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The nvec argument must be > 0");
431:     ctx->nvec = nvec;
432:   }
433:   switch (dos) {
434:     case EPS_EVSL_DOS_KPM:
435:       if (deg == PETSC_DECIDE || deg == PETSC_DEFAULT) ctx->deg = 300;
436:       else {
437:         PetscCheck(deg>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The deg argument must be > 0");
438:         ctx->deg = deg;
439:       }
440:       break;
441:     case EPS_EVSL_DOS_LANCZOS:
442:       if (steps == PETSC_DECIDE || steps == PETSC_DEFAULT) ctx->steps = 40;
443:       else {
444:         PetscCheck(steps>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The steps argument must be > 0");
445:         ctx->steps = steps;
446:       }
447:       if (npoints == PETSC_DECIDE || npoints == PETSC_DEFAULT) ctx->npoints = 200;
448:       else {
449:         PetscCheck(npoints>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The npoints argument must be > 0");
450:         ctx->npoints = npoints;
451:       }
452:       break;
453:   }
454:   eps->state = EPS_STATE_INITIAL;
455:   PetscFunctionReturn(PETSC_SUCCESS);
456: }

458: /*@
459:    EPSEVSLSetDOSParameters - Defines the parameters used for computing the
460:    density of states (DOS) in the EVSL solver.

462:    Logically Collective

464:    Input Parameters:
465: +  eps     - the eigensolver context
466: .  dos     - DOS method, either KPM or Lanczos
467: .  nvec    - number of sample vectors
468: .  deg     - polynomial degree (KPM only)
469: .  steps   - number of Lanczos steps (Lanczos only)
470: -  npoints - number of sample points (Lanczos only)

472:    Options Database Keys:
473: +  -eps_evsl_dos_method <dos> - set the DOS method, either kpm or lanczos
474: .  -eps_evsl_dos_nvec <n> - set the number of sample vectors
475: .  -eps_evsl_dos_degree <n> - set the polynomial degree
476: .  -eps_evsl_dos_steps <n> - set the number of Lanczos steps
477: -  -eps_evsl_dos_npoints <n> - set the number of sample points

479:    Notes:
480:    The density of states (or spectral density) can be approximated with two
481:    methods, kernel polynomial method (KPM) or Lanczos. Some parameters for
482:    these methods can be set by the user with this function, with some of
483:    them being relevant for one of the methods only.

485:    Level: intermediate

487: .seealso: EPSEVSLGetDOSParameters()
488: @*/
489: PetscErrorCode EPSEVSLSetDOSParameters(EPS eps,EPSEVSLDOSMethod dos,PetscInt nvec,PetscInt deg,PetscInt steps,PetscInt npoints)
490: {
491:   PetscFunctionBegin;
498:   PetscTryMethod(eps,"EPSEVSLSetDOSParameters_C",(EPS,EPSEVSLDOSMethod,PetscInt,PetscInt,PetscInt,PetscInt),(eps,dos,nvec,deg,steps,npoints));
499:   PetscFunctionReturn(PETSC_SUCCESS);
500: }

502: static PetscErrorCode EPSEVSLGetDOSParameters_EVSL(EPS eps,EPSEVSLDOSMethod *dos,PetscInt *nvec,PetscInt *deg,PetscInt *steps,PetscInt *npoints)
503: {
504:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

506:   PetscFunctionBegin;
507:   if (dos)     *dos     = ctx->dos;
508:   if (nvec)    *nvec    = ctx->nvec;
509:   if (deg)     *deg     = ctx->deg;
510:   if (steps)   *steps   = ctx->steps;
511:   if (npoints) *npoints = ctx->npoints;
512:   PetscFunctionReturn(PETSC_SUCCESS);
513: }

515: /*@
516:    EPSEVSLGetDOSParameters - Gets the parameters used for computing the
517:    density of states (DOS) in the EVSL solver.

519:    Not Collective

521:    Input Parameter:
522: .  eps - the eigensolver context

524:    Output Parameters:
525: +  dos     - DOS method, either KPM or Lanczos
526: .  nvec    - number of sample vectors
527: .  deg     - polynomial degree (KPM only)
528: .  steps   - number of Lanczos steps (Lanczos only)
529: -  npoints - number of sample points (Lanczos only)

531:    Level: intermediate

533: .seealso: EPSEVSLSetDOSParameters()
534: @*/
535: PetscErrorCode EPSEVSLGetDOSParameters(EPS eps,EPSEVSLDOSMethod *dos,PetscInt *nvec,PetscInt *deg,PetscInt *steps,PetscInt *npoints)
536: {
537:   PetscFunctionBegin;
539:   PetscUseMethod(eps,"EPSEVSLGetDOSParameters_C",(EPS,EPSEVSLDOSMethod*,PetscInt*,PetscInt*,PetscInt*,PetscInt*),(eps,dos,nvec,deg,steps,npoints));
540:   PetscFunctionReturn(PETSC_SUCCESS);
541: }

543: static PetscErrorCode EPSEVSLSetPolParameters_EVSL(EPS eps,PetscInt max_deg,PetscReal thresh)
544: {
545:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

547:   PetscFunctionBegin;
548:   if (max_deg == PETSC_DECIDE || max_deg == PETSC_DEFAULT) ctx->max_deg = 10000;
549:   else {
550:     PetscCheck(max_deg>2,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The max_deg argument must be > 2");
551:     ctx->max_deg = max_deg;
552:   }
553:   if (thresh == (PetscReal)PETSC_DECIDE || thresh == (PetscReal)PETSC_DEFAULT) ctx->thresh = 0.8;
554:   else {
555:     PetscCheck(thresh>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The thresh argument must be > 0.0");
556:     ctx->thresh = thresh;
557:   }
558:   eps->state = EPS_STATE_INITIAL;
559:   PetscFunctionReturn(PETSC_SUCCESS);
560: }

562: /*@
563:    EPSEVSLSetPolParameters - Defines the parameters used for building the
564:    building the polynomial in the EVSL solver.

566:    Logically Collective

568:    Input Parameters:
569: +  eps     - the eigensolver context
570: .  max_deg - maximum degree allowed for the polynomial
571: -  thresh  - threshold for accepting polynomial

573:    Options Database Keys:
574: +  -eps_evsl_pol_max_deg <d> - set maximum polynomial degree
575: -  -eps_evsl_pol_thresh <t> - set the threshold

577:    Level: intermediate

579: .seealso: EPSEVSLGetPolParameters()
580: @*/
581: PetscErrorCode EPSEVSLSetPolParameters(EPS eps,PetscInt max_deg,PetscReal thresh)
582: {
583:   PetscFunctionBegin;
587:   PetscTryMethod(eps,"EPSEVSLSetPolParameters_C",(EPS,PetscInt,PetscReal),(eps,max_deg,thresh));
588:   PetscFunctionReturn(PETSC_SUCCESS);
589: }

591: static PetscErrorCode EPSEVSLGetPolParameters_EVSL(EPS eps,PetscInt *max_deg,PetscReal *thresh)
592: {
593:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

595:   PetscFunctionBegin;
596:   if (max_deg) *max_deg = ctx->max_deg;
597:   if (thresh)  *thresh  = ctx->thresh;
598:   PetscFunctionReturn(PETSC_SUCCESS);
599: }

601: /*@
602:    EPSEVSLGetPolParameters - Gets the parameters used for building the
603:    polynomial in the EVSL solver.

605:    Not Collective

607:    Input Parameter:
608: .  eps - the eigensolver context

610:    Output Parameters:
611: +  max_deg - the maximum degree of the polynomial
612: -  thresh  - the threshold

614:    Level: intermediate

616: .seealso: EPSEVSLSetPolParameters()
617: @*/
618: PetscErrorCode EPSEVSLGetPolParameters(EPS eps,PetscInt *max_deg,PetscReal *thresh)
619: {
620:   PetscFunctionBegin;
622:   PetscUseMethod(eps,"EPSEVSLGetPolParameters_C",(EPS,PetscInt*,PetscReal*),(eps,max_deg,thresh));
623:   PetscFunctionReturn(PETSC_SUCCESS);
624: }

626: static PetscErrorCode EPSEVSLSetDamping_EVSL(EPS eps,EPSEVSLDamping damping)
627: {
628:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

630:   PetscFunctionBegin;
631:   if (ctx->damping != damping) {
632:     ctx->damping = damping;
633:     eps->state   = EPS_STATE_INITIAL;
634:   }
635:   PetscFunctionReturn(PETSC_SUCCESS);
636: }

638: /*@
639:    EPSEVSLSetDamping - Set the type of damping to be used in EVSL.

641:    Logically Collective

643:    Input Parameters:
644: +  eps     - the eigensolver context
645: -  damping - the type of damping

647:    Options Database Key:
648: .  -eps_evsl_damping <n> - set the type of damping

650:    Notes:
651:    Damping is applied when building the polynomial to be used when solving the
652:    eigenproblem, and also during estimation of DOS with the KPM method.

654:    Level: intermediate

656: .seealso: EPSEVSLGetDamping(), EPSEVSLSetDOSParameters()
657: @*/
658: PetscErrorCode EPSEVSLSetDamping(EPS eps,EPSEVSLDamping damping)
659: {
660:   PetscFunctionBegin;
663:   PetscTryMethod(eps,"EPSEVSLSetDamping_C",(EPS,EPSEVSLDamping),(eps,damping));
664:   PetscFunctionReturn(PETSC_SUCCESS);
665: }

667: static PetscErrorCode EPSEVSLGetDamping_EVSL(EPS eps,EPSEVSLDamping *damping)
668: {
669:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

671:   PetscFunctionBegin;
672:   *damping = ctx->damping;
673:   PetscFunctionReturn(PETSC_SUCCESS);
674: }

676: /*@
677:    EPSEVSLGetDamping - Gets the type of damping.

679:    Not Collective

681:    Input Parameter:
682: .  eps - the eigensolver context

684:    Output Parameter:
685: .  damping - the type of damping

687:    Level: intermediate

689: .seealso: EPSEVSLSetDamping()
690: @*/
691: PetscErrorCode EPSEVSLGetDamping(EPS eps,EPSEVSLDamping *damping)
692: {
693:   PetscFunctionBegin;
695:   PetscAssertPointer(damping,2);
696:   PetscUseMethod(eps,"EPSEVSLGetDamping_C",(EPS,EPSEVSLDamping*),(eps,damping));
697:   PetscFunctionReturn(PETSC_SUCCESS);
698: }

700: static PetscErrorCode EPSView_EVSL(EPS eps,PetscViewer viewer)
701: {
702:   PetscBool      isascii;
703:   EPS_EVSL       *ctx = (EPS_EVSL*)eps->data;

705:   PetscFunctionBegin;
706:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
707:   if (isascii) {
708:     PetscCall(PetscViewerASCIIPrintf(viewer,"  numerical range = [%g,%g]\n",(double)ctx->lmin,(double)ctx->lmax));
709:     PetscCall(PetscViewerASCIIPrintf(viewer,"  number of slices = %" PetscInt_FMT "\n",ctx->nslices));
710:     PetscCall(PetscViewerASCIIPrintf(viewer,"  type of damping = %s\n",EPSEVSLDampings[ctx->damping]));
711:     PetscCall(PetscViewerASCIIPrintf(viewer,"  computing DOS with %s: nvec=%" PetscInt_FMT ", ",EPSEVSLDOSMethods[ctx->dos],ctx->nvec));
712:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
713:     switch (ctx->dos) {
714:       case EPS_EVSL_DOS_KPM:
715:         PetscCall(PetscViewerASCIIPrintf(viewer,"degree=%" PetscInt_FMT "\n",ctx->deg));
716:         break;
717:       case EPS_EVSL_DOS_LANCZOS:
718:         PetscCall(PetscViewerASCIIPrintf(viewer,"steps=%" PetscInt_FMT ", npoints=%" PetscInt_FMT "\n",ctx->steps,ctx->npoints));
719:         break;
720:     }
721:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
722:     PetscCall(PetscViewerASCIIPrintf(viewer,"  polynomial parameters: max degree = %" PetscInt_FMT ", threshold = %g\n",ctx->max_deg,(double)ctx->thresh));
723:   }
724:   PetscFunctionReturn(PETSC_SUCCESS);
725: }

727: static PetscErrorCode EPSSetFromOptions_EVSL(EPS eps,PetscOptionItems *PetscOptionsObject)
728: {
729:   PetscReal        array[2]={0,0},th;
730:   PetscInt         k,i1,i2,i3,i4;
731:   PetscBool        flg,flg1;
732:   EPSEVSLDOSMethod dos;
733:   EPSEVSLDamping   damping;
734:   EPS_EVSL         *ctx = (EPS_EVSL*)eps->data;

736:   PetscFunctionBegin;
737:   PetscOptionsHeadBegin(PetscOptionsObject,"EPS EVSL Options");

739:     k = 2;
740:     PetscCall(PetscOptionsRealArray("-eps_evsl_range","Interval containing all eigenvalues (two real values separated with a comma without spaces)","EPSEVSLSetRange",array,&k,&flg));
741:     if (flg) {
742:       PetscCheck(k>1,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_SIZ,"Must pass two values in -eps_evsl_range (comma-separated without spaces)");
743:       PetscCall(EPSEVSLSetRange(eps,array[0],array[1]));
744:     }

746:     PetscCall(PetscOptionsInt("-eps_evsl_slices","Number of slices","EPSEVSLSetSlices",ctx->nslices,&i1,&flg));
747:     if (flg) PetscCall(EPSEVSLSetSlices(eps,i1));

749:     PetscCall(PetscOptionsEnum("-eps_evsl_damping","Type of damping","EPSEVSLSetDamping",EPSEVSLDampings,(PetscEnum)ctx->damping,(PetscEnum*)&damping,&flg));
750:     if (flg) PetscCall(EPSEVSLSetDamping(eps,damping));

752:     PetscCall(EPSEVSLGetDOSParameters(eps,&dos,&i1,&i2,&i3,&i4));
753:     PetscCall(PetscOptionsEnum("-eps_evsl_dos_method","Method to compute the DOS","EPSEVSLSetDOSParameters",EPSEVSLDOSMethods,(PetscEnum)ctx->dos,(PetscEnum*)&dos,&flg));
754:     PetscCall(PetscOptionsInt("-eps_evsl_dos_nvec","Number of sample vectors for DOS","EPSEVSLSetDOSParameters",i1,&i1,&flg1));
755:     if (flg1) flg = PETSC_TRUE;
756:     PetscCall(PetscOptionsInt("-eps_evsl_dos_degree","Polynomial degree used for DOS","EPSEVSLSetDOSParameters",i2,&i2,&flg1));
757:     if (flg1) flg = PETSC_TRUE;
758:     PetscCall(PetscOptionsInt("-eps_evsl_dos_steps","Number of Lanczos steps in DOS","EPSEVSLSetDOSParameters",i3,&i3,&flg1));
759:     if (flg1) flg = PETSC_TRUE;
760:     PetscCall(PetscOptionsInt("-eps_evsl_dos_npoints","Number of sample points used for DOS","EPSEVSLSetDOSParameters",i4,&i4,&flg1));
761:     if (flg || flg1) PetscCall(EPSEVSLSetDOSParameters(eps,dos,i1,i2,i3,i4));

763:     PetscCall(EPSEVSLGetPolParameters(eps,&i1,&th));
764:     PetscCall(PetscOptionsInt("-eps_evsl_pol_max_deg","Maximum degree allowed for the polynomial","EPSEVSLSetPolParameters",i1,&i1,&flg));
765:     PetscCall(PetscOptionsReal("-eps_evsl_pol_threshold","Threshold for accepting polynomial","EPSEVSLSetPolParameters",th,&th,&flg1));
766:     if (flg || flg1) PetscCall(EPSEVSLSetPolParameters(eps,i1,th));

768:   PetscOptionsHeadEnd();
769:   PetscFunctionReturn(PETSC_SUCCESS);
770: }

772: static PetscErrorCode EPSDestroy_EVSL(EPS eps)
773: {
774:   EPS_EVSL       *ctx = (EPS_EVSL*)eps->data;

776:   PetscFunctionBegin;
777:   if (ctx->initialized) EVSLFinish();
778:   PetscCall(PetscLayoutDestroy(&ctx->map));
779:   PetscCall(PetscFree(ctx->sli));
780:   PetscCall(PetscFree(eps->data));
781:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetRange_C",NULL));
782:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetRange_C",NULL));
783:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetSlices_C",NULL));
784:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetSlices_C",NULL));
785:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDOSParameters_C",NULL));
786:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDOSParameters_C",NULL));
787:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetPolParameters_C",NULL));
788:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetPolParameters_C",NULL));
789:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDamping_C",NULL));
790:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDamping_C",NULL));
791:   PetscFunctionReturn(PETSC_SUCCESS);
792: }

794: static PetscErrorCode EPSReset_EVSL(EPS eps)
795: {
796:   EPS_EVSL       *ctx = (EPS_EVSL*)eps->data;

798:   PetscFunctionBegin;
799:   PetscCall(MatDestroy(&ctx->A));
800:   PetscCall(VecDestroy(&ctx->x));
801:   PetscCall(VecDestroy(&ctx->y));
802:   PetscFunctionReturn(PETSC_SUCCESS);
803: }

805: SLEPC_EXTERN PetscErrorCode EPSCreate_EVSL(EPS eps)
806: {
807:   EPS_EVSL       *ctx;

809:   PetscFunctionBegin;
810:   PetscCall(PetscNew(&ctx));
811:   eps->data = (void*)ctx;

813:   ctx->nslices = 0;
814:   ctx->lmin    = PETSC_MIN_REAL;
815:   ctx->lmax    = PETSC_MAX_REAL;
816:   ctx->dos     = EPS_EVSL_DOS_KPM;
817:   ctx->nvec    = 80;
818:   ctx->deg     = 300;
819:   ctx->steps   = 40;
820:   ctx->npoints = 200;
821:   ctx->max_deg = 10000;
822:   ctx->thresh  = 0.8;
823:   ctx->damping = EPS_EVSL_DAMPING_SIGMA;

825:   eps->categ = EPS_CATEGORY_OTHER;

827:   eps->ops->solve          = EPSSolve_EVSL;
828:   eps->ops->setup          = EPSSetUp_EVSL;
829:   eps->ops->setupsort      = EPSSetUpSort_Basic;
830:   eps->ops->setfromoptions = EPSSetFromOptions_EVSL;
831:   eps->ops->destroy        = EPSDestroy_EVSL;
832:   eps->ops->reset          = EPSReset_EVSL;
833:   eps->ops->view           = EPSView_EVSL;
834:   eps->ops->backtransform  = EPSBackTransform_Default;
835:   eps->ops->setdefaultst   = EPSSetDefaultST_NoFactor;

837:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetRange_C",EPSEVSLSetRange_EVSL));
838:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetRange_C",EPSEVSLGetRange_EVSL));
839:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetSlices_C",EPSEVSLSetSlices_EVSL));
840:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetSlices_C",EPSEVSLGetSlices_EVSL));
841:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDOSParameters_C",EPSEVSLSetDOSParameters_EVSL));
842:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDOSParameters_C",EPSEVSLGetDOSParameters_EVSL));
843:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetPolParameters_C",EPSEVSLSetPolParameters_EVSL));
844:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetPolParameters_C",EPSEVSLGetPolParameters_EVSL));
845:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDamping_C",EPSEVSLSetDamping_EVSL));
846:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDamping_C",EPSEVSLGetDamping_EVSL));
847:   PetscFunctionReturn(PETSC_SUCCESS);
848: }