Actual source code: evsl.c
slepc-main 2024-12-17
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: EPSCheckNotStructured(eps);
75: if (eps->nev==0) eps->nev = 1;
76: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift));
77: PetscCheck(isshift,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations");
79: if (ctx->initialized) EVSLFinish();
80: EVSLStart();
81: ctx->initialized=PETSC_TRUE;
83: /* get number of slices per process */
84: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size));
85: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank));
86: if (!ctx->nslices) ctx->nslices = size;
87: PetscCall(PetscLayoutDestroy(&ctx->map));
88: PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)eps),PETSC_DECIDE,ctx->nslices,1,&ctx->map));
90: /* get matrix and prepare auxiliary vectors */
91: PetscCall(MatDestroy(&ctx->A));
92: PetscCall(STGetMatrix(eps->st,0,&A));
93: if (size==1) {
94: PetscCall(PetscObjectReference((PetscObject)A));
95: ctx->A = A;
96: } else PetscCall(MatCreateRedundantMatrix(A,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&ctx->A));
97: SetAMatvec(eps->n,&AMatvec_EVSL,(void*)ctx);
98: if (!ctx->x) PetscCall(MatCreateVecsEmpty(ctx->A,&ctx->x,&ctx->y));
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: PetscCheck(eps->which==EPS_ALL && eps->inta!=eps->intb,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: PetscCall(MatCreateVecs(ctx->A,&v0,NULL));
108: if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
109: PetscCall(BVGetRandomContext(eps->V,&rnd));
110: PetscCall(VecSetRandom(v0,rnd));
111: PetscCall(VecGetArray(v0,&vinit));
112: PetscCallEVSL(LanTrbounds,50,200,eps->tol,vinit,1,&ctx->lmin,&ctx->lmax,NULL);
113: PetscCall(VecRestoreArray(v0,&vinit));
114: PetscCall(VecDestroy(&v0));
115: ctx->estimrange = PETSC_TRUE; /* estimate if called again with another matrix */
116: }
117: 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);
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: switch (ctx->dos) {
125: case EPS_EVSL_DOS_KPM:
126: PetscCall(PetscMalloc1(ctx->deg+1,&mu));
127: if (!rank) PetscCallEVSL(kpmdos,ctx->deg,(int)ctx->damping,ctx->nvec,xintv,mu,&ecount);
128: PetscCallMPI(MPI_Bcast(mu,ctx->deg+1,MPIU_REAL,0,PetscObjectComm((PetscObject)eps)));
129: break;
130: case EPS_EVSL_DOS_LANCZOS:
131: PetscCall(PetscMalloc2(ctx->npoints,&xdos,ctx->npoints,&ydos));
132: if (!rank) PetscCallEVSL(LanDos,ctx->nvec,PetscMin(ctx->steps,eps->n/2),ctx->npoints,xdos,ydos,&ecount,xintv);
133: PetscCallMPI(MPI_Bcast(xdos,ctx->npoints,MPIU_REAL,0,PetscObjectComm((PetscObject)eps)));
134: PetscCallMPI(MPI_Bcast(ydos,ctx->npoints,MPIU_REAL,0,PetscObjectComm((PetscObject)eps)));
135: break;
136: default:
137: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid DOS method");
138: }
139: PetscCallMPI(MPI_Bcast(&ecount,1,MPIU_REAL,0,PetscObjectComm((PetscObject)eps)));
141: PetscCall(PetscInfo(eps,"Estimated eigenvalue count in the interval: %g\n",ecount));
142: eps->ncv = (PetscInt)PetscCeilReal(1.5*ecount);
144: /* slice the spectrum */
145: PetscCall(PetscFree(ctx->sli));
146: PetscCall(PetscMalloc1(ctx->nslices+1,&ctx->sli));
147: if (ctx->dos == EPS_EVSL_DOS_KPM) {
148: PetscCallEVSL(spslicer,ctx->sli,mu,ctx->deg,xintv,ctx->nslices,10*(PetscInt)ecount);
149: PetscCall(PetscFree(mu));
150: } else if (ctx->dos == EPS_EVSL_DOS_LANCZOS) {
151: spslicer2(xdos,ydos,ctx->nslices,ctx->npoints,ctx->sli);
152: PetscCall(PetscFree2(xdos,ydos));
153: }
155: /* approximate number of eigenvalues wanted in each slice */
156: ctx->nev = (PetscInt)(1.0 + ecount/(PetscReal)ctx->nslices) + 2;
158: if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
159: if (eps->max_it==PETSC_DETERMINE) eps->max_it = 1;
160: PetscCall(EPSAllocateSolution(eps,0));
161: PetscFunctionReturn(PETSC_SUCCESS);
162: }
164: static PetscErrorCode EPSSolve_EVSL(EPS eps)
165: {
166: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
167: PetscInt i,j,k=0,sl,mlan,nevout,*ind,nevmax,rstart,rend,*nevloc,*disp,N;
168: PetscReal *res,xintv[4],*errest;
169: PetscScalar *lam,*X,*Y,*vinit,*eigr;
170: PetscMPIInt size,rank;
171: PetscRandom rnd;
172: Vec v,w,v0,x;
173: VecScatter vs;
174: IS is;
175: polparams pol;
177: PetscFunctionBegin;
178: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size));
179: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank));
180: PetscCall(PetscLayoutGetRange(ctx->map,&rstart,&rend));
181: nevmax = (rend-rstart)*ctx->nev;
182: PetscCall(MatCreateVecs(ctx->A,&v0,NULL));
183: PetscCall(BVGetRandomContext(eps->V,&rnd));
184: PetscCall(VecSetRandom(v0,rnd));
185: PetscCall(VecGetArray(v0,&vinit));
186: PetscCall(PetscMalloc5(size,&nevloc,size+1,&disp,nevmax,&eigr,nevmax,&errest,nevmax*eps->n,&X));
187: mlan = PetscMin(PetscMax(5*ctx->nev,300),eps->n);
188: for (sl=rstart; sl<rend; sl++) {
189: xintv[0] = ctx->sli[sl];
190: xintv[1] = ctx->sli[sl+1];
191: xintv[2] = ctx->lmin;
192: xintv[3] = ctx->lmax;
193: PetscCall(PetscInfo(ctx->A,"Subinterval %" PetscInt_FMT ": [%.4e, %.4e]\n",sl+1,xintv[0],xintv[1]));
194: set_pol_def(&pol);
195: pol.max_deg = ctx->max_deg;
196: pol.damping = (int)ctx->damping;
197: pol.thresh_int = ctx->thresh;
198: find_pol(xintv,&pol);
199: PetscCall(PetscInfo(ctx->A,"Polynomial [type = %" PetscInt_FMT "], deg %" PetscInt_FMT ", bar %e gam %e\n",pol.type,pol.deg,pol.bar,pol.gam));
200: PetscCallEVSL(ChebLanNr,xintv,mlan,eps->tol,vinit,&pol,&nevout,&lam,&Y,&res,NULL);
201: PetscCheck(k+nevout<=nevmax,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Too low estimation of eigenvalue count, try modifying the sampling parameters");
202: free_pol(&pol);
203: PetscCall(PetscInfo(ctx->A,"Computed %" PetscInt_FMT " eigenvalues\n",nevout));
204: PetscCall(PetscMalloc1(nevout,&ind));
205: sort_double(nevout,lam,ind);
206: for (i=0;i<nevout;i++) {
207: eigr[i+k] = lam[i];
208: errest[i+k] = res[ind[i]];
209: PetscCall(PetscArraycpy(X+(i+k)*eps->n,Y+ind[i]*eps->n,eps->n));
210: }
211: k += nevout;
212: if (lam) evsl_Free(lam);
213: if (Y) evsl_Free_device(Y);
214: if (res) evsl_Free(res);
215: PetscCall(PetscFree(ind));
216: }
217: PetscCall(VecRestoreArray(v0,&vinit));
218: PetscCall(VecDestroy(&v0));
220: /* gather eigenvalues computed by each MPI process */
221: PetscCallMPI(MPI_Allgather(&k,1,MPIU_INT,nevloc,1,MPIU_INT,PetscObjectComm((PetscObject)eps)));
222: eps->nev = nevloc[0];
223: disp[0] = 0;
224: for (i=1;i<size;i++) {
225: eps->nev += nevloc[i];
226: disp[i] = disp[i-1]+nevloc[i-1];
227: }
228: disp[size] = disp[size-1]+nevloc[size-1];
229: PetscCheck(eps->nev<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Too low estimation of eigenvalue count, try modifying the sampling parameters");
230: PetscCallMPI(MPI_Allgatherv(eigr,k,MPIU_SCALAR,eps->eigr,nevloc,disp,MPIU_SCALAR,PetscObjectComm((PetscObject)eps)));
231: PetscCallMPI(MPI_Allgatherv(errest,k,MPIU_REAL,eps->errest,nevloc,disp,MPIU_REAL,PetscObjectComm((PetscObject)eps)));
232: eps->nconv = eps->nev;
233: eps->its = 1;
234: eps->reason = EPS_CONVERGED_TOL;
236: /* scatter computed eigenvectors and store them in eps->V */
237: PetscCall(BVCreateVec(eps->V,&w));
238: for (i=0;i<size;i++) {
239: N = (rank==i)? eps->n: 0;
240: PetscCall(VecCreateSeq(PETSC_COMM_SELF,N,&x));
241: PetscCall(VecSetFromOptions(x));
242: PetscCall(ISCreateStride(PETSC_COMM_SELF,N,0,1,&is));
243: PetscCall(VecScatterCreate(x,is,w,is,&vs));
244: PetscCall(ISDestroy(&is));
245: for (j=disp[i];j<disp[i+1];j++) {
246: PetscCall(BVGetColumn(eps->V,j,&v));
247: if (rank==i) PetscCall(VecPlaceArray(x,X+(j-disp[i])*eps->n));
248: PetscCall(VecScatterBegin(vs,x,v,INSERT_VALUES,SCATTER_FORWARD));
249: PetscCall(VecScatterEnd(vs,x,v,INSERT_VALUES,SCATTER_FORWARD));
250: if (rank==i) PetscCall(VecResetArray(x));
251: PetscCall(BVRestoreColumn(eps->V,j,&v));
252: }
253: PetscCall(VecScatterDestroy(&vs));
254: PetscCall(VecDestroy(&x));
255: }
256: PetscCall(VecDestroy(&w));
257: PetscCall(PetscFree5(nevloc,disp,eigr,errest,X));
258: PetscFunctionReturn(PETSC_SUCCESS);
259: }
261: static PetscErrorCode EPSEVSLSetSlices_EVSL(EPS eps,PetscInt nslices)
262: {
263: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
265: PetscFunctionBegin;
266: if (nslices == PETSC_DECIDE || nslices == PETSC_DEFAULT) nslices = 0;
267: else PetscCheck(nslices>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Number of slices must be 1 at least");
268: if (ctx->nslices != nslices) {
269: ctx->nslices = nslices;
270: eps->state = EPS_STATE_INITIAL;
271: }
272: PetscFunctionReturn(PETSC_SUCCESS);
273: }
275: /*@
276: EPSEVSLSetSlices - Set the number of slices in which the interval must be
277: subdivided.
279: Logically Collective
281: Input Parameters:
282: + eps - the eigensolver context
283: - nslices - the number of slices
285: Options Database Key:
286: . -eps_evsl_slices <n> - set the number of slices to n
288: Notes:
289: By default, one slice per MPI process is used. Depending on the number of
290: eigenvalues, using more slices may be beneficial, but very narrow subintervals
291: imply higher polynomial degree.
293: Level: intermediate
295: .seealso: EPSEVSLGetSlices()
296: @*/
297: PetscErrorCode EPSEVSLSetSlices(EPS eps,PetscInt nslices)
298: {
299: PetscFunctionBegin;
302: PetscTryMethod(eps,"EPSEVSLSetSlices_C",(EPS,PetscInt),(eps,nslices));
303: PetscFunctionReturn(PETSC_SUCCESS);
304: }
306: static PetscErrorCode EPSEVSLGetSlices_EVSL(EPS eps,PetscInt *nslices)
307: {
308: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
310: PetscFunctionBegin;
311: *nslices = ctx->nslices;
312: PetscFunctionReturn(PETSC_SUCCESS);
313: }
315: /*@
316: EPSEVSLGetSlices - Gets the number of slices in which the interval must be
317: subdivided.
319: Not Collective
321: Input Parameter:
322: . eps - the eigensolver context
324: Output Parameter:
325: . nslices - the number of slices
327: Level: intermediate
329: .seealso: EPSEVSLSetSlices()
330: @*/
331: PetscErrorCode EPSEVSLGetSlices(EPS eps,PetscInt *nslices)
332: {
333: PetscFunctionBegin;
335: PetscAssertPointer(nslices,2);
336: PetscUseMethod(eps,"EPSEVSLGetSlices_C",(EPS,PetscInt*),(eps,nslices));
337: PetscFunctionReturn(PETSC_SUCCESS);
338: }
340: static PetscErrorCode EPSEVSLSetRange_EVSL(EPS eps,PetscReal lmin,PetscReal lmax)
341: {
342: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
344: PetscFunctionBegin;
345: PetscCheck(lmin<lmax,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: PetscFunctionReturn(PETSC_SUCCESS);
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
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: {
382: PetscFunctionBegin;
386: PetscTryMethod(eps,"EPSEVSLSetRange_C",(EPS,PetscReal,PetscReal),(eps,lmin,lmax));
387: PetscFunctionReturn(PETSC_SUCCESS);
388: }
390: static PetscErrorCode EPSEVSLGetRange_EVSL(EPS eps,PetscReal *lmin,PetscReal *lmax)
391: {
392: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
394: PetscFunctionBegin;
395: if (lmin) *lmin = ctx->lmin;
396: if (lmax) *lmax = ctx->lmax;
397: PetscFunctionReturn(PETSC_SUCCESS);
398: }
400: /*@
401: EPSEVSLGetRange - Gets the interval containing all eigenvalues.
403: Not Collective
405: Input Parameter:
406: . eps - the eigensolver context
408: Output Parameters:
409: + lmin - left end of the interval
410: - lmax - right end of the interval
412: Level: intermediate
414: .seealso: EPSEVSLSetRange()
415: @*/
416: PetscErrorCode EPSEVSLGetRange(EPS eps,PetscReal *lmin,PetscReal *lmax)
417: {
418: PetscFunctionBegin;
420: PetscUseMethod(eps,"EPSEVSLGetRange_C",(EPS,PetscReal*,PetscReal*),(eps,lmin,lmax));
421: PetscFunctionReturn(PETSC_SUCCESS);
422: }
424: static PetscErrorCode EPSEVSLSetDOSParameters_EVSL(EPS eps,EPSEVSLDOSMethod dos,PetscInt nvec,PetscInt deg,PetscInt steps,PetscInt npoints)
425: {
426: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
428: PetscFunctionBegin;
429: ctx->dos = dos;
430: if (nvec == PETSC_DETERMINE) ctx->nvec = 80;
431: else if (nvec != PETSC_CURRENT) {
432: PetscCheck(nvec>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The nvec argument must be > 0");
433: ctx->nvec = nvec;
434: }
435: switch (dos) {
436: case EPS_EVSL_DOS_KPM:
437: if (deg == PETSC_DETERMINE) ctx->deg = 300;
438: else if (deg != PETSC_CURRENT) {
439: PetscCheck(deg>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The deg argument must be > 0");
440: ctx->deg = deg;
441: }
442: break;
443: case EPS_EVSL_DOS_LANCZOS:
444: if (steps == PETSC_DETERMINE) ctx->steps = 40;
445: else if (steps != PETSC_CURRENT) {
446: PetscCheck(steps>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The steps argument must be > 0");
447: ctx->steps = steps;
448: }
449: if (npoints == PETSC_DETERMINE) ctx->npoints = 200;
450: else if (npoints != PETSC_CURRENT) {
451: PetscCheck(npoints>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The npoints argument must be > 0");
452: ctx->npoints = npoints;
453: }
454: break;
455: }
456: eps->state = EPS_STATE_INITIAL;
457: PetscFunctionReturn(PETSC_SUCCESS);
458: }
460: /*@
461: EPSEVSLSetDOSParameters - Defines the parameters used for computing the
462: density of states (DOS) in the EVSL solver.
464: Logically Collective
466: Input Parameters:
467: + eps - the eigensolver context
468: . dos - DOS method, either KPM or Lanczos
469: . nvec - number of sample vectors
470: . deg - polynomial degree (KPM only)
471: . steps - number of Lanczos steps (Lanczos only)
472: - npoints - number of sample points (Lanczos only)
474: Options Database Keys:
475: + -eps_evsl_dos_method <dos> - set the DOS method, either kpm or lanczos
476: . -eps_evsl_dos_nvec <n> - set the number of sample vectors
477: . -eps_evsl_dos_degree <n> - set the polynomial degree
478: . -eps_evsl_dos_steps <n> - set the number of Lanczos steps
479: - -eps_evsl_dos_npoints <n> - set the number of sample points
481: Notes:
482: The density of states (or spectral density) can be approximated with two
483: methods, kernel polynomial method (KPM) or Lanczos. Some parameters for
484: these methods can be set by the user with this function, with some of
485: them being relevant for one of the methods only.
487: For the integer argumens, you can use PETSC_CURRENT to keep the current
488: value, and PETSC_DETERMINE to set them to a reasonable default.
490: Level: intermediate
492: .seealso: EPSEVSLGetDOSParameters()
493: @*/
494: PetscErrorCode EPSEVSLSetDOSParameters(EPS eps,EPSEVSLDOSMethod dos,PetscInt nvec,PetscInt deg,PetscInt steps,PetscInt npoints)
495: {
496: PetscFunctionBegin;
503: PetscTryMethod(eps,"EPSEVSLSetDOSParameters_C",(EPS,EPSEVSLDOSMethod,PetscInt,PetscInt,PetscInt,PetscInt),(eps,dos,nvec,deg,steps,npoints));
504: PetscFunctionReturn(PETSC_SUCCESS);
505: }
507: static PetscErrorCode EPSEVSLGetDOSParameters_EVSL(EPS eps,EPSEVSLDOSMethod *dos,PetscInt *nvec,PetscInt *deg,PetscInt *steps,PetscInt *npoints)
508: {
509: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
511: PetscFunctionBegin;
512: if (dos) *dos = ctx->dos;
513: if (nvec) *nvec = ctx->nvec;
514: if (deg) *deg = ctx->deg;
515: if (steps) *steps = ctx->steps;
516: if (npoints) *npoints = ctx->npoints;
517: PetscFunctionReturn(PETSC_SUCCESS);
518: }
520: /*@
521: EPSEVSLGetDOSParameters - Gets the parameters used for computing the
522: density of states (DOS) in the EVSL solver.
524: Not Collective
526: Input Parameter:
527: . eps - the eigensolver context
529: Output Parameters:
530: + dos - DOS method, either KPM or Lanczos
531: . nvec - number of sample vectors
532: . deg - polynomial degree (KPM only)
533: . steps - number of Lanczos steps (Lanczos only)
534: - npoints - number of sample points (Lanczos only)
536: Level: intermediate
538: .seealso: EPSEVSLSetDOSParameters()
539: @*/
540: PetscErrorCode EPSEVSLGetDOSParameters(EPS eps,EPSEVSLDOSMethod *dos,PetscInt *nvec,PetscInt *deg,PetscInt *steps,PetscInt *npoints)
541: {
542: PetscFunctionBegin;
544: PetscUseMethod(eps,"EPSEVSLGetDOSParameters_C",(EPS,EPSEVSLDOSMethod*,PetscInt*,PetscInt*,PetscInt*,PetscInt*),(eps,dos,nvec,deg,steps,npoints));
545: PetscFunctionReturn(PETSC_SUCCESS);
546: }
548: static PetscErrorCode EPSEVSLSetPolParameters_EVSL(EPS eps,PetscInt max_deg,PetscReal thresh)
549: {
550: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
552: PetscFunctionBegin;
553: if (max_deg == PETSC_DETERMINE) ctx->max_deg = 10000;
554: else if (max_deg != PETSC_CURRENT) {
555: PetscCheck(max_deg>2,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The max_deg argument must be > 2");
556: ctx->max_deg = max_deg;
557: }
558: if (thresh == (PetscReal)PETSC_DETERMINE) ctx->thresh = 0.8;
559: else if (thresh != (PetscReal)PETSC_CURRENT) {
560: PetscCheck(thresh>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The thresh argument must be > 0.0");
561: ctx->thresh = thresh;
562: }
563: eps->state = EPS_STATE_INITIAL;
564: PetscFunctionReturn(PETSC_SUCCESS);
565: }
567: /*@
568: EPSEVSLSetPolParameters - Defines the parameters used for building the
569: building the polynomial in the EVSL solver.
571: Logically Collective
573: Input Parameters:
574: + eps - the eigensolver context
575: . max_deg - maximum degree allowed for the polynomial
576: - thresh - threshold for accepting polynomial
578: Options Database Keys:
579: + -eps_evsl_pol_max_deg <d> - set maximum polynomial degree
580: - -eps_evsl_pol_thresh <t> - set the threshold
582: Note:
583: PETSC_CURRENT can be used to preserve the current value of any of the
584: arguments, and PETSC_DETERMINE to set them to a default value.
586: Level: intermediate
588: .seealso: EPSEVSLGetPolParameters()
589: @*/
590: PetscErrorCode EPSEVSLSetPolParameters(EPS eps,PetscInt max_deg,PetscReal thresh)
591: {
592: PetscFunctionBegin;
596: PetscTryMethod(eps,"EPSEVSLSetPolParameters_C",(EPS,PetscInt,PetscReal),(eps,max_deg,thresh));
597: PetscFunctionReturn(PETSC_SUCCESS);
598: }
600: static PetscErrorCode EPSEVSLGetPolParameters_EVSL(EPS eps,PetscInt *max_deg,PetscReal *thresh)
601: {
602: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
604: PetscFunctionBegin;
605: if (max_deg) *max_deg = ctx->max_deg;
606: if (thresh) *thresh = ctx->thresh;
607: PetscFunctionReturn(PETSC_SUCCESS);
608: }
610: /*@
611: EPSEVSLGetPolParameters - Gets the parameters used for building the
612: polynomial in the EVSL solver.
614: Not Collective
616: Input Parameter:
617: . eps - the eigensolver context
619: Output Parameters:
620: + max_deg - the maximum degree of the polynomial
621: - thresh - the threshold
623: Level: intermediate
625: .seealso: EPSEVSLSetPolParameters()
626: @*/
627: PetscErrorCode EPSEVSLGetPolParameters(EPS eps,PetscInt *max_deg,PetscReal *thresh)
628: {
629: PetscFunctionBegin;
631: PetscUseMethod(eps,"EPSEVSLGetPolParameters_C",(EPS,PetscInt*,PetscReal*),(eps,max_deg,thresh));
632: PetscFunctionReturn(PETSC_SUCCESS);
633: }
635: static PetscErrorCode EPSEVSLSetDamping_EVSL(EPS eps,EPSEVSLDamping damping)
636: {
637: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
639: PetscFunctionBegin;
640: if (ctx->damping != damping) {
641: ctx->damping = damping;
642: eps->state = EPS_STATE_INITIAL;
643: }
644: PetscFunctionReturn(PETSC_SUCCESS);
645: }
647: /*@
648: EPSEVSLSetDamping - Set the type of damping to be used in EVSL.
650: Logically Collective
652: Input Parameters:
653: + eps - the eigensolver context
654: - damping - the type of damping
656: Options Database Key:
657: . -eps_evsl_damping <n> - set the type of damping
659: Notes:
660: Damping is applied when building the polynomial to be used when solving the
661: eigenproblem, and also during estimation of DOS with the KPM method.
663: Level: intermediate
665: .seealso: EPSEVSLGetDamping(), EPSEVSLSetDOSParameters()
666: @*/
667: PetscErrorCode EPSEVSLSetDamping(EPS eps,EPSEVSLDamping damping)
668: {
669: PetscFunctionBegin;
672: PetscTryMethod(eps,"EPSEVSLSetDamping_C",(EPS,EPSEVSLDamping),(eps,damping));
673: PetscFunctionReturn(PETSC_SUCCESS);
674: }
676: static PetscErrorCode EPSEVSLGetDamping_EVSL(EPS eps,EPSEVSLDamping *damping)
677: {
678: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
680: PetscFunctionBegin;
681: *damping = ctx->damping;
682: PetscFunctionReturn(PETSC_SUCCESS);
683: }
685: /*@
686: EPSEVSLGetDamping - Gets the type of damping.
688: Not Collective
690: Input Parameter:
691: . eps - the eigensolver context
693: Output Parameter:
694: . damping - the type of damping
696: Level: intermediate
698: .seealso: EPSEVSLSetDamping()
699: @*/
700: PetscErrorCode EPSEVSLGetDamping(EPS eps,EPSEVSLDamping *damping)
701: {
702: PetscFunctionBegin;
704: PetscAssertPointer(damping,2);
705: PetscUseMethod(eps,"EPSEVSLGetDamping_C",(EPS,EPSEVSLDamping*),(eps,damping));
706: PetscFunctionReturn(PETSC_SUCCESS);
707: }
709: static PetscErrorCode EPSView_EVSL(EPS eps,PetscViewer viewer)
710: {
711: PetscBool isascii;
712: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
714: PetscFunctionBegin;
715: PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
716: if (isascii) {
717: PetscCall(PetscViewerASCIIPrintf(viewer," numerical range = [%g,%g]\n",(double)ctx->lmin,(double)ctx->lmax));
718: PetscCall(PetscViewerASCIIPrintf(viewer," number of slices = %" PetscInt_FMT "\n",ctx->nslices));
719: PetscCall(PetscViewerASCIIPrintf(viewer," type of damping = %s\n",EPSEVSLDampings[ctx->damping]));
720: PetscCall(PetscViewerASCIIPrintf(viewer," computing DOS with %s: nvec=%" PetscInt_FMT ", ",EPSEVSLDOSMethods[ctx->dos],ctx->nvec));
721: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
722: switch (ctx->dos) {
723: case EPS_EVSL_DOS_KPM:
724: PetscCall(PetscViewerASCIIPrintf(viewer,"degree=%" PetscInt_FMT "\n",ctx->deg));
725: break;
726: case EPS_EVSL_DOS_LANCZOS:
727: PetscCall(PetscViewerASCIIPrintf(viewer,"steps=%" PetscInt_FMT ", npoints=%" PetscInt_FMT "\n",ctx->steps,ctx->npoints));
728: break;
729: }
730: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
731: PetscCall(PetscViewerASCIIPrintf(viewer," polynomial parameters: max degree = %" PetscInt_FMT ", threshold = %g\n",ctx->max_deg,(double)ctx->thresh));
732: }
733: PetscFunctionReturn(PETSC_SUCCESS);
734: }
736: static PetscErrorCode EPSSetFromOptions_EVSL(EPS eps,PetscOptionItems *PetscOptionsObject)
737: {
738: PetscReal array[2]={0,0},th;
739: PetscInt k,i1,i2,i3,i4;
740: PetscBool flg,flg1;
741: EPSEVSLDOSMethod dos;
742: EPSEVSLDamping damping;
743: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
745: PetscFunctionBegin;
746: PetscOptionsHeadBegin(PetscOptionsObject,"EPS EVSL Options");
748: k = 2;
749: PetscCall(PetscOptionsRealArray("-eps_evsl_range","Interval containing all eigenvalues (two real values separated with a comma without spaces)","EPSEVSLSetRange",array,&k,&flg));
750: if (flg) {
751: PetscCheck(k>1,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_SIZ,"Must pass two values in -eps_evsl_range (comma-separated without spaces)");
752: PetscCall(EPSEVSLSetRange(eps,array[0],array[1]));
753: }
755: PetscCall(PetscOptionsInt("-eps_evsl_slices","Number of slices","EPSEVSLSetSlices",ctx->nslices,&i1,&flg));
756: if (flg) PetscCall(EPSEVSLSetSlices(eps,i1));
758: PetscCall(PetscOptionsEnum("-eps_evsl_damping","Type of damping","EPSEVSLSetDamping",EPSEVSLDampings,(PetscEnum)ctx->damping,(PetscEnum*)&damping,&flg));
759: if (flg) PetscCall(EPSEVSLSetDamping(eps,damping));
761: PetscCall(EPSEVSLGetDOSParameters(eps,&dos,&i1,&i2,&i3,&i4));
762: PetscCall(PetscOptionsEnum("-eps_evsl_dos_method","Method to compute the DOS","EPSEVSLSetDOSParameters",EPSEVSLDOSMethods,(PetscEnum)ctx->dos,(PetscEnum*)&dos,&flg));
763: PetscCall(PetscOptionsInt("-eps_evsl_dos_nvec","Number of sample vectors for DOS","EPSEVSLSetDOSParameters",i1,&i1,&flg1));
764: if (flg1) flg = PETSC_TRUE;
765: PetscCall(PetscOptionsInt("-eps_evsl_dos_degree","Polynomial degree used for DOS","EPSEVSLSetDOSParameters",i2,&i2,&flg1));
766: if (flg1) flg = PETSC_TRUE;
767: PetscCall(PetscOptionsInt("-eps_evsl_dos_steps","Number of Lanczos steps in DOS","EPSEVSLSetDOSParameters",i3,&i3,&flg1));
768: if (flg1) flg = PETSC_TRUE;
769: PetscCall(PetscOptionsInt("-eps_evsl_dos_npoints","Number of sample points used for DOS","EPSEVSLSetDOSParameters",i4,&i4,&flg1));
770: if (flg || flg1) PetscCall(EPSEVSLSetDOSParameters(eps,dos,i1,i2,i3,i4));
772: PetscCall(EPSEVSLGetPolParameters(eps,&i1,&th));
773: PetscCall(PetscOptionsInt("-eps_evsl_pol_max_deg","Maximum degree allowed for the polynomial","EPSEVSLSetPolParameters",i1,&i1,&flg));
774: PetscCall(PetscOptionsReal("-eps_evsl_pol_threshold","Threshold for accepting polynomial","EPSEVSLSetPolParameters",th,&th,&flg1));
775: if (flg || flg1) PetscCall(EPSEVSLSetPolParameters(eps,i1,th));
777: PetscOptionsHeadEnd();
778: PetscFunctionReturn(PETSC_SUCCESS);
779: }
781: static PetscErrorCode EPSDestroy_EVSL(EPS eps)
782: {
783: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
785: PetscFunctionBegin;
786: if (ctx->initialized) EVSLFinish();
787: PetscCall(PetscLayoutDestroy(&ctx->map));
788: PetscCall(PetscFree(ctx->sli));
789: PetscCall(PetscFree(eps->data));
790: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetRange_C",NULL));
791: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetRange_C",NULL));
792: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetSlices_C",NULL));
793: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetSlices_C",NULL));
794: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDOSParameters_C",NULL));
795: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDOSParameters_C",NULL));
796: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetPolParameters_C",NULL));
797: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetPolParameters_C",NULL));
798: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDamping_C",NULL));
799: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDamping_C",NULL));
800: PetscFunctionReturn(PETSC_SUCCESS);
801: }
803: static PetscErrorCode EPSReset_EVSL(EPS eps)
804: {
805: EPS_EVSL *ctx = (EPS_EVSL*)eps->data;
807: PetscFunctionBegin;
808: PetscCall(MatDestroy(&ctx->A));
809: PetscCall(VecDestroy(&ctx->x));
810: PetscCall(VecDestroy(&ctx->y));
811: PetscFunctionReturn(PETSC_SUCCESS);
812: }
814: SLEPC_EXTERN PetscErrorCode EPSCreate_EVSL(EPS eps)
815: {
816: EPS_EVSL *ctx;
818: PetscFunctionBegin;
819: PetscCall(PetscNew(&ctx));
820: eps->data = (void*)ctx;
822: ctx->nslices = 0;
823: ctx->lmin = PETSC_MIN_REAL;
824: ctx->lmax = PETSC_MAX_REAL;
825: ctx->dos = EPS_EVSL_DOS_KPM;
826: ctx->nvec = 80;
827: ctx->deg = 300;
828: ctx->steps = 40;
829: ctx->npoints = 200;
830: ctx->max_deg = 10000;
831: ctx->thresh = 0.8;
832: ctx->damping = EPS_EVSL_DAMPING_SIGMA;
834: eps->categ = EPS_CATEGORY_OTHER;
836: eps->ops->solve = EPSSolve_EVSL;
837: eps->ops->setup = EPSSetUp_EVSL;
838: eps->ops->setupsort = EPSSetUpSort_Basic;
839: eps->ops->setfromoptions = EPSSetFromOptions_EVSL;
840: eps->ops->destroy = EPSDestroy_EVSL;
841: eps->ops->reset = EPSReset_EVSL;
842: eps->ops->view = EPSView_EVSL;
843: eps->ops->backtransform = EPSBackTransform_Default;
844: eps->ops->setdefaultst = EPSSetDefaultST_NoFactor;
846: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetRange_C",EPSEVSLSetRange_EVSL));
847: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetRange_C",EPSEVSLGetRange_EVSL));
848: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetSlices_C",EPSEVSLSetSlices_EVSL));
849: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetSlices_C",EPSEVSLGetSlices_EVSL));
850: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDOSParameters_C",EPSEVSLSetDOSParameters_EVSL));
851: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDOSParameters_C",EPSEVSLGetDOSParameters_EVSL));
852: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetPolParameters_C",EPSEVSLSetPolParameters_EVSL));
853: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetPolParameters_C",EPSEVSLGetPolParameters_EVSL));
854: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLSetDamping_C",EPSEVSLSetDamping_EVSL));
855: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSEVSLGetDamping_C",EPSEVSLGetDamping_EVSL));
856: PetscFunctionReturn(PETSC_SUCCESS);
857: }