Actual source code: evsl.c

slepc-3.16.2 2022-02-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: typedef struct {
 18:   PetscBool         initialized;
 19:   Mat               A;           /* problem matrix */
 20:   Vec               x,y;         /* auxiliary vectors */
 21:   PetscReal         *sli;        /* slice bounds */
 22:   PetscInt          nev;         /* approximate number of wanted eigenvalues in each slice */
 23:   PetscLayout       map;         /* used to distribute slices among MPI processes */
 24:   PetscBool         estimrange;  /* the filter range was not set by the user */
 25:   /* user parameters */
 26:   PetscInt          nslices;     /* number of slices */
 27:   PetscReal         lmin,lmax;   /* numerical range (min and max eigenvalue) */
 28:   EPSEVSLDOSMethod  dos;         /* DOS method, either KPM or Lanczos */
 29:   PetscInt          nvec;        /* number of sample vectors used for DOS */
 30:   PetscInt          deg;         /* polynomial degree used for DOS (KPM only) */
 31:   PetscInt          steps;       /* number of Lanczos steps used for DOS (Lanczos only) */
 32:   PetscInt          npoints;     /* number of sample points used for DOS (Lanczos only) */
 33:   PetscInt          max_deg;     /* maximum degree allowed for the polynomial */
 34:   PetscReal         thresh;      /* threshold for accepting polynomial */
 35:   EPSEVSLDamping    damping;     /* type of damping (for polynomial and for DOS-KPM) */
 36: } EPS_EVSL;

 38: static void AMatvec_EVSL(double *xa,double *ya,void *data)
 39: {
 41:   EPS_EVSL       *ctx = (EPS_EVSL*)data;
 42:   Vec            x = ctx->x,y = ctx->y;
 43:   Mat            A = ctx->A;

 46:   VecPlaceArray(x,(PetscScalar*)xa);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
 47:   VecPlaceArray(y,(PetscScalar*)ya);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
 48:   MatMult(A,x,y);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
 49:   VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
 50:   VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
 51:   PetscFunctionReturnVoid();
 52: }

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

 67:   EPSCheckStandard(eps);
 68:   EPSCheckHermitian(eps);
 69:   PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
 70:   if (!isshift) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations");

 72:   if (ctx->initialized) EVSLFinish();
 73:   EVSLStart();
 74:   ctx->initialized=PETSC_TRUE;

 76:   /* get number of slices per process */
 77:   MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size);
 78:   MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank);
 79:   if (!ctx->nslices) ctx->nslices = size;
 80:   PetscLayoutDestroy(&ctx->map);
 81:   PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)eps),PETSC_DECIDE,ctx->nslices,1,&ctx->map);

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

102:   if (!eps->which) eps->which=EPS_ALL;
103:   if (eps->which!=EPS_ALL || eps->inta==eps->intb) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires setting an interval with EPSSetInterval()");

105:   /* estimate numerical range */
106:   if (ctx->estimrange || ctx->lmin == PETSC_MIN_REAL || ctx->lmax == PETSC_MAX_REAL) {
107:     MatCreateVecs(ctx->A,&v0,NULL);
108:     if (!eps->V) { EPSGetBV(eps,&eps->V); }
109:     BVGetRandomContext(eps->V,&rnd);
110:     VecSetRandom(v0,rnd);
111:     VecGetArray(v0,&vinit);
112:     LanTrbounds(50,200,eps->tol,vinit,1,&ctx->lmin,&ctx->lmax,NULL);
113:     VecRestoreArray(v0,&vinit);
114:     VecDestroy(&v0);
115:     ctx->estimrange = PETSC_TRUE;   /* estimate if called again with another matrix */
116:   }
117:   if (ctx->lmin > eps->inta || ctx->lmax < eps->intb) SETERRQ4(PetscObjectComm((PetscObject)eps),1,"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);
118:   xintv[0] = eps->inta;
119:   xintv[1] = eps->intb;
120:   xintv[2] = ctx->lmin;
121:   xintv[3] = ctx->lmax;

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

136:   PetscInfo1(eps,"Estimated eigenvalue count in the interval: %g\n",ecount);
137:   eps->ncv = (PetscInt)PetscCeilReal(1.5*ecount);

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

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

153:   if (eps->mpd!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
154:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
155:   EPSAllocateSolution(eps,0);
156:   return(0);
157: }

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

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

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

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

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

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

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

275:    Logically Collective on eps

277:    Input Parameters:
278: +  eps     - the eigensolver context
279: -  nslices - the number of slices

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

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

289:    Level: intermediate

291: .seealso: EPSEVSLGetSlices()
292: @*/
293: PetscErrorCode EPSEVSLSetSlices(EPS eps,PetscInt nslices)
294: {

300:   PetscTryMethod(eps,"EPSEVSLSetSlices_C",(EPS,PetscInt),(eps,nslices));
301:   return(0);
302: }

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

309:   *nslices = ctx->nslices;
310:   return(0);
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: {

336:   PetscUseMethod(eps,"EPSEVSLGetSlices_C",(EPS,PetscInt*),(eps,nslices));
337:   return(0);
338: }

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

345:   if (lmin>lmax) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be lmin<lmax");
346:   if (ctx->lmin != lmin || ctx->lmax != lmax) {
347:     ctx->lmin  = lmin;
348:     ctx->lmax  = lmax;
349:     eps->state = EPS_STATE_INITIAL;
350:   }
351:   return(0);
352: }

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

358:    Logically Collective on eps

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

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

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

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

376:    Level: intermediate

378: .seealso: EPSEVSLGetRange(), EPSSetInterval()
379: @*/
380: PetscErrorCode EPSEVSLSetRange(EPS eps,PetscReal lmin,PetscReal lmax)
381: {

388:   PetscTryMethod(eps,"EPSEVSLSetRange_C",(EPS,PetscReal,PetscReal),(eps,lmin,lmax));
389:   return(0);
390: }

392: static PetscErrorCode EPSEVSLGetRange_EVSL(EPS eps,PetscReal *lmin,PetscReal *lmax)
393: {
394:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

397:   if (lmin) *lmin = ctx->lmin;
398:   if (lmax) *lmax = ctx->lmax;
399:   return(0);
400: }

402: /*@
403:    EPSEVSLGetRange - Gets the interval containing all eigenvalues.

405:    Not Collective

407:    Input Parameter:
408: .  eps - the eigensolver context

410:    Output Parameters:
411: +  lmin - left end of the interval
412: -  lmax - right end of the interval

414:    Level: intermediate

416: .seealso: EPSEVSLSetRange()
417: @*/
418: PetscErrorCode EPSEVSLGetRange(EPS eps,PetscReal *lmin,PetscReal *lmax)
419: {

424:   PetscUseMethod(eps,"EPSEVSLGetRange_C",(EPS,PetscReal*,PetscReal*),(eps,lmin,lmax));
425:   return(0);
426: }

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

433:   ctx->dos = dos;
434:   if (nvec == PETSC_DECIDE || nvec == PETSC_DEFAULT) ctx->nvec = 80;
435:   else if (nvec<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The nvec argument must be > 0");
436:   else ctx->nvec = nvec;
437:   switch (dos) {
438:     case EPS_EVSL_DOS_KPM:
439:       if (deg == PETSC_DECIDE || deg == PETSC_DEFAULT) ctx->deg = 300;
440:       else if (deg<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The deg argument must be > 0");
441:       else ctx->deg = deg;
442:       break;
443:     case EPS_EVSL_DOS_LANCZOS:
444:       if (steps == PETSC_DECIDE || steps == PETSC_DEFAULT) ctx->steps = 40;
445:       else if (steps<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The steps argument must be > 0");
446:       else ctx->steps = steps;
447:       if (npoints == PETSC_DECIDE || npoints == PETSC_DEFAULT) ctx->npoints = 200;
448:       else if (npoints<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The npoints argument must be > 0");
449:       else ctx->npoints = npoints;
450:       break;
451:   }
452:   eps->state = EPS_STATE_INITIAL;
453:   return(0);
454: }

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

460:    Logically Collective on eps

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

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

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

483:    Level: intermediate

485: .seealso: EPSEVSLGetDOSParameters()
486: @*/
487: PetscErrorCode EPSEVSLSetDOSParameters(EPS eps,EPSEVSLDOSMethod dos,PetscInt nvec,PetscInt deg,PetscInt steps,PetscInt npoints)
488: {

498:   PetscTryMethod(eps,"EPSEVSLSetDOSParameters_C",(EPS,EPSEVSLDOSMethod,PetscInt,PetscInt,PetscInt,PetscInt),(eps,dos,nvec,deg,steps,npoints));
499:   return(0);
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;

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:   return(0);
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: {

541:   PetscUseMethod(eps,"EPSEVSLGetDOSParameters_C",(EPS,EPSEVSLDOSMethod*,PetscInt*,PetscInt*,PetscInt*,PetscInt*),(eps,dos,nvec,deg,steps,npoints));
542:   return(0);
543: }

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

550:   if (max_deg == PETSC_DECIDE || max_deg == PETSC_DEFAULT) ctx->max_deg = 10000;
551:   else if (max_deg<3) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The max_deg argument must be > 2");
552:   else ctx->max_deg = max_deg;
553:   if (thresh == PETSC_DECIDE || thresh == PETSC_DEFAULT) ctx->thresh = 0.8;
554:   else if (thresh<0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The thresh argument must be > 0.0");
555:   else ctx->thresh = thresh;
556:   eps->state = EPS_STATE_INITIAL;
557:   return(0);
558: }

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

564:    Logically Collective on eps

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

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

575:    Level: intermediate

577: .seealso: EPSEVSLGetPolParameters()
578: @*/
579: PetscErrorCode EPSEVSLSetPolParameters(EPS eps,PetscInt max_deg,PetscReal thresh)
580: {

587:   PetscTryMethod(eps,"EPSEVSLSetPolParameters_C",(EPS,PetscInt,PetscReal),(eps,max_deg,thresh));
588:   return(0);
589: }

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

596:   if (max_deg) *max_deg = ctx->max_deg;
597:   if (thresh)  *thresh  = ctx->thresh;
598:   return(0);
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: {

624:   PetscUseMethod(eps,"EPSEVSLGetPolParameters_C",(EPS,PetscInt*,PetscReal*),(eps,max_deg,thresh));
625:   return(0);
626: }

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

633:   if (ctx->damping != damping) {
634:     ctx->damping = damping;
635:     eps->state   = EPS_STATE_INITIAL;
636:   }
637:   return(0);
638: }

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

643:    Logically Collective on eps

645:    Input Parameters:
646: +  eps     - the eigensolver context
647: -  damping - the type of damping

649:    Options Database Key:
650: .  -eps_evsl_damping <n> - set the type of damping

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

656:    Level: intermediate

658: .seealso: EPSEVSLGetDamping(), EPSEVSLSetDOSParameters()
659: @*/
660: PetscErrorCode EPSEVSLSetDamping(EPS eps,EPSEVSLDamping damping)
661: {

667:   PetscTryMethod(eps,"EPSEVSLSetDamping_C",(EPS,EPSEVSLDamping),(eps,damping));
668:   return(0);
669: }

671: static PetscErrorCode EPSEVSLGetDamping_EVSL(EPS eps,EPSEVSLDamping *damping)
672: {
673:   EPS_EVSL *ctx = (EPS_EVSL*)eps->data;

676:   *damping = ctx->damping;
677:   return(0);
678: }

680: /*@
681:    EPSEVSLGetDamping - Gets the type of damping.

683:    Not Collective

685:    Input Parameter:
686: .  eps - the eigensolver context

688:    Output Parameter:
689: .  damping - the type of damping

691:    Level: intermediate

693: .seealso: EPSEVSLSetDamping()
694: @*/
695: PetscErrorCode EPSEVSLGetDamping(EPS eps,EPSEVSLDamping *damping)
696: {

702:   PetscUseMethod(eps,"EPSEVSLGetDamping_C",(EPS,EPSEVSLDamping*),(eps,damping));
703:   return(0);
704: }

706: PetscErrorCode EPSView_EVSL(EPS eps,PetscViewer viewer)
707: {
709:   PetscBool      isascii;
710:   EPS_EVSL       *ctx = (EPS_EVSL*)eps->data;

713:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
714:   if (isascii) {
715:     PetscViewerASCIIPrintf(viewer,"  numerical range = [%g,%g]\n",(double)ctx->lmin,(double)ctx->lmax);
716:     PetscViewerASCIIPrintf(viewer,"  number of slices = %D\n",ctx->nslices);
717:     PetscViewerASCIIPrintf(viewer,"  type of damping = %s\n",EPSEVSLDampings[ctx->damping]);
718:     PetscViewerASCIIPrintf(viewer,"  computing DOS with %s: nvec=%D, ",EPSEVSLDOSMethods[ctx->dos],ctx->nvec);
719:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
720:     switch (ctx->dos) {
721:       case EPS_EVSL_DOS_KPM:
722:         PetscViewerASCIIPrintf(viewer,"degree=%D\n",ctx->deg);
723:         break;
724:       case EPS_EVSL_DOS_LANCZOS:
725:         PetscViewerASCIIPrintf(viewer,"steps=%D, npoints=%D\n",ctx->steps,ctx->npoints);
726:         break;
727:     }
728:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
729:     PetscViewerASCIIPrintf(viewer,"  polynomial parameters: max degree = %D, threshold = %g\n",ctx->max_deg,(double)ctx->thresh);
730:   }
731:   return(0);
732: }

734: PetscErrorCode EPSSetFromOptions_EVSL(PetscOptionItems *PetscOptionsObject,EPS eps)
735: {
736:   PetscErrorCode   ierr;
737:   PetscReal        array[2]={0,0},th;
738:   PetscInt         k,i1,i2,i3,i4;
739:   PetscBool        flg,flg1;
740:   EPSEVSLDOSMethod dos;
741:   EPSEVSLDamping   damping;
742:   EPS_EVSL         *ctx = (EPS_EVSL*)eps->data;

745:   PetscOptionsHead(PetscOptionsObject,"EPS EVSL Options");

747:     k = 2;
748:     PetscOptionsRealArray("-eps_evsl_range","Interval containing all eigenvalues (two real values separated with a comma without spaces)","EPSEVSLSetRange",array,&k,&flg);
749:     if (flg) {
750:       if (k<2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_SIZ,"Must pass two values in -eps_evsl_range (comma-separated without spaces)");
751:       EPSEVSLSetRange(eps,array[0],array[1]);
752:     }

754:     PetscOptionsInt("-eps_evsl_slices","Number of slices","EPSEVSLSetSlices",ctx->nslices,&i1,&flg);
755:     if (flg) { EPSEVSLSetSlices(eps,i1); }

757:     PetscOptionsEnum("-eps_evsl_damping","Type of damping","EPSEVSLSetDamping",EPSEVSLDampings,(PetscEnum)ctx->damping,(PetscEnum*)&damping,&flg);
758:     if (flg) { EPSEVSLSetDamping(eps,damping); }

760:     EPSEVSLGetDOSParameters(eps,&dos,&i1,&i2,&i3,&i4);
761:     PetscOptionsEnum("-eps_evsl_dos_method","Method to compute the DOS","EPSEVSLSetDOSParameters",EPSEVSLDOSMethods,(PetscEnum)ctx->dos,(PetscEnum*)&dos,&flg);
762:     PetscOptionsInt("-eps_evsl_dos_nvec","Number of sample vectors for DOS","EPSEVSLSetDOSParameters",i1,&i1,&flg1);
763:     if (flg1) flg = PETSC_TRUE;
764:     PetscOptionsInt("-eps_evsl_dos_degree","Polynomial degree used for DOS","EPSEVSLSetDOSParameters",i2,&i2,&flg1);
765:     if (flg1) flg = PETSC_TRUE;
766:     PetscOptionsInt("-eps_evsl_dos_steps","Number of Lanczos steps in DOS","EPSEVSLSetDOSParameters",i3,&i3,&flg1);
767:     if (flg1) flg = PETSC_TRUE;
768:     PetscOptionsInt("-eps_evsl_dos_npoints","Number of sample points used for DOS","EPSEVSLSetDOSParameters",i4,&i4,&flg1);
769:     if (flg || flg1) { EPSEVSLSetDOSParameters(eps,dos,i1,i2,i3,i4); }

771:     EPSEVSLGetPolParameters(eps,&i1,&th);
772:     PetscOptionsInt("-eps_evsl_pol_max_deg","Maximum degree allowed for the polynomial","EPSEVSLSetPolParameters",i1,&i1,&flg);
773:     PetscOptionsReal("-eps_evsl_pol_threshold","Threshold for accepting polynomial","EPSEVSLSetPolParameters",th,&th,&flg1);
774:     if (flg || flg1) { EPSEVSLSetPolParameters(eps,i1,th); }

776:   PetscOptionsTail();
777:   return(0);
778: }

780: PetscErrorCode EPSDestroy_EVSL(EPS eps)
781: {
783:   EPS_EVSL       *ctx = (EPS_EVSL*)eps->data;

786:   if (ctx->initialized) EVSLFinish();
787:   PetscLayoutDestroy(&ctx->map);
788:   PetscFree(ctx->sli);
789:   PetscFree(eps->data);
790:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetRange_C",NULL);
791:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetRange_C",NULL);
792:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetSlices_C",NULL);
793:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetSlices_C",NULL);
794:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDOSParameters_C",NULL);
795:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDOSParameters_C",NULL);
796:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetPolParameters_C",NULL);
797:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetPolParameters_C",NULL);
798:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDamping_C",NULL);
799:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDamping_C",NULL);
800:   return(0);
801: }

803: PetscErrorCode EPSReset_EVSL(EPS eps)
804: {
806:   EPS_EVSL       *ctx = (EPS_EVSL*)eps->data;

809:   MatDestroy(&ctx->A);
810:   VecDestroy(&ctx->x);
811:   VecDestroy(&ctx->y);
812:   return(0);
813: }

815: SLEPC_EXTERN PetscErrorCode EPSCreate_EVSL(EPS eps)
816: {
817:   EPS_EVSL       *ctx;

821:   PetscNewLog(eps,&ctx);
822:   eps->data = (void*)ctx;

824:   ctx->nslices = 0;
825:   ctx->lmin    = PETSC_MIN_REAL;
826:   ctx->lmax    = PETSC_MAX_REAL;
827:   ctx->dos     = EPS_EVSL_DOS_KPM;
828:   ctx->nvec    = 80;
829:   ctx->deg     = 300;
830:   ctx->steps   = 40;
831:   ctx->npoints = 200;
832:   ctx->max_deg = 10000;
833:   ctx->thresh  = 0.8;
834:   ctx->damping = EPS_EVSL_DAMPING_SIGMA;

836:   eps->categ = EPS_CATEGORY_OTHER;

838:   eps->ops->solve          = EPSSolve_EVSL;
839:   eps->ops->setup          = EPSSetUp_EVSL;
840:   eps->ops->setupsort      = EPSSetUpSort_Basic;
841:   eps->ops->setfromoptions = EPSSetFromOptions_EVSL;
842:   eps->ops->destroy        = EPSDestroy_EVSL;
843:   eps->ops->reset          = EPSReset_EVSL;
844:   eps->ops->view           = EPSView_EVSL;
845:   eps->ops->backtransform  = EPSBackTransform_Default;
846:   eps->ops->setdefaultst   = EPSSetDefaultST_NoFactor;

848:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetRange_C",EPSEVSLSetRange_EVSL);
849:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetRange_C",EPSEVSLGetRange_EVSL);
850:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetSlices_C",EPSEVSLSetSlices_EVSL);
851:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetSlices_C",EPSEVSLGetSlices_EVSL);
852:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDOSParameters_C",EPSEVSLSetDOSParameters_EVSL);
853:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDOSParameters_C",EPSEVSLGetDOSParameters_EVSL);
854:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetPolParameters_C",EPSEVSLSetPolParameters_EVSL);
855:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetPolParameters_C",EPSEVSLGetPolParameters_EVSL);
856:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDamping_C",EPSEVSLSetDamping_EVSL);
857:   PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDamping_C",EPSEVSLGetDamping_EVSL);
858:   return(0);
859: }