Actual source code: epsimpl.h
slepc-3.14.2 2021-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: */
11: #if !defined(SLEPCEPSIMPL_H)
12: #define SLEPCEPSIMPL_H
14: #include <slepceps.h>
15: #include <slepc/private/slepcimpl.h>
17: SLEPC_EXTERN PetscBool EPSRegisterAllCalled;
18: SLEPC_EXTERN PetscErrorCode EPSRegisterAll(void);
19: SLEPC_EXTERN PetscLogEvent EPS_SetUp,EPS_Solve;
21: typedef struct _EPSOps *EPSOps;
23: struct _EPSOps {
24: PetscErrorCode (*solve)(EPS);
25: PetscErrorCode (*setup)(EPS);
26: PetscErrorCode (*setupsort)(EPS);
27: PetscErrorCode (*setfromoptions)(PetscOptionItems*,EPS);
28: PetscErrorCode (*publishoptions)(EPS);
29: PetscErrorCode (*destroy)(EPS);
30: PetscErrorCode (*reset)(EPS);
31: PetscErrorCode (*view)(EPS,PetscViewer);
32: PetscErrorCode (*backtransform)(EPS);
33: PetscErrorCode (*computevectors)(EPS);
34: PetscErrorCode (*setdefaultst)(EPS);
35: };
37: /*
38: Maximum number of monitors you can run with a single EPS
39: */
40: #define MAXEPSMONITORS 5
42: /*
43: The solution process goes through several states
44: */
45: typedef enum { EPS_STATE_INITIAL,
46: EPS_STATE_SETUP,
47: EPS_STATE_SOLVED,
48: EPS_STATE_EIGENVECTORS } EPSStateType;
50: /*
51: To classify the different solvers into categories
52: */
53: typedef enum { EPS_CATEGORY_KRYLOV, /* Krylov solver: relies on STApply and STBackTransform (same as OTHER) */
54: EPS_CATEGORY_PRECOND, /* Preconditioned solver: uses ST only to manage preconditioner */
55: EPS_CATEGORY_CONTOUR, /* Contour integral: ST used to solve linear systems at integration points */
56: EPS_CATEGORY_OTHER } EPSSolverType;
58: /*
59: To check for unsupported features at EPSSetUp_XXX()
60: */
61: typedef enum { EPS_FEATURE_BALANCE=1, /* balancing */
62: EPS_FEATURE_ARBITRARY=2, /* arbitrary selection of eigepairs */
63: EPS_FEATURE_REGION=4, /* nontrivial region for filtering */
64: EPS_FEATURE_EXTRACTION=8, /* extraction technique different from Ritz */
65: EPS_FEATURE_CONVERGENCE=16, /* convergence test selected by user */
66: EPS_FEATURE_STOPPING=32, /* stopping test */
67: EPS_FEATURE_TWOSIDED=64 /* two-sided variant */
68: } EPSFeatureType;
70: /*
71: Defines the EPS data structure
72: */
73: struct _p_EPS {
74: PETSCHEADER(struct _EPSOps);
75: /*------------------------- User parameters ---------------------------*/
76: PetscInt max_it; /* maximum number of iterations */
77: PetscInt nev; /* number of eigenvalues to compute */
78: PetscInt ncv; /* number of basis vectors */
79: PetscInt mpd; /* maximum dimension of projected problem */
80: PetscInt nini,ninil; /* number of initial vectors (negative means not copied yet) */
81: PetscInt nds; /* number of basis vectors of deflation space */
82: PetscScalar target; /* target value */
83: PetscReal tol; /* tolerance */
84: EPSConv conv; /* convergence test */
85: EPSStop stop; /* stopping test */
86: EPSWhich which; /* which part of the spectrum to be sought */
87: PetscReal inta,intb; /* interval [a,b] for spectrum slicing */
88: EPSProblemType problem_type; /* which kind of problem to be solved */
89: EPSExtraction extraction; /* which kind of extraction to be applied */
90: EPSBalance balance; /* the balancing method */
91: PetscInt balance_its; /* number of iterations of the balancing method */
92: PetscReal balance_cutoff; /* cutoff value for balancing */
93: PetscBool trueres; /* whether the true residual norm must be computed */
94: PetscBool trackall; /* whether all the residuals must be computed */
95: PetscBool purify; /* whether eigenvectors need to be purified */
96: PetscBool twosided; /* whether to compute left eigenvectors (two-sided solver) */
98: /*-------------- User-provided functions and contexts -----------------*/
99: PetscErrorCode (*converged)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
100: PetscErrorCode (*convergeduser)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
101: PetscErrorCode (*convergeddestroy)(void*);
102: PetscErrorCode (*stopping)(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);
103: PetscErrorCode (*stoppinguser)(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);
104: PetscErrorCode (*stoppingdestroy)(void*);
105: PetscErrorCode (*arbitrary)(PetscScalar,PetscScalar,Vec,Vec,PetscScalar*,PetscScalar*,void*);
106: void *convergedctx;
107: void *stoppingctx;
108: void *arbitraryctx;
109: PetscErrorCode (*monitor[MAXEPSMONITORS])(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
110: PetscErrorCode (*monitordestroy[MAXEPSMONITORS])(void**);
111: void *monitorcontext[MAXEPSMONITORS];
112: PetscInt numbermonitors;
114: /*----------------- Child objects and working data -------------------*/
115: ST st; /* spectral transformation object */
116: DS ds; /* direct solver object */
117: DS dsts; /* auxiliary direct solver object used in two-sided case */
118: BV V; /* set of basis vectors and computed eigenvectors */
119: BV W; /* left basis vectors (if left eigenvectors requested) */
120: RG rg; /* optional region for filtering */
121: SlepcSC sc; /* sorting criterion data */
122: Vec D; /* diagonal matrix for balancing */
123: Vec *IS,*ISL; /* references to user-provided initial spaces */
124: Vec *defl; /* references to user-provided deflation space */
125: PetscScalar *eigr,*eigi; /* real and imaginary parts of eigenvalues */
126: PetscReal *errest; /* error estimates */
127: PetscScalar *rr,*ri; /* values computed by user's arbitrary selection function */
128: PetscInt *perm; /* permutation for eigenvalue ordering */
129: PetscInt nwork; /* number of work vectors */
130: Vec *work; /* work vectors */
131: void *data; /* placeholder for solver-specific stuff */
133: /* ----------------------- Status variables --------------------------*/
134: EPSStateType state; /* initial -> setup -> solved -> eigenvectors */
135: EPSSolverType categ; /* solver category */
136: PetscInt nconv; /* number of converged eigenvalues */
137: PetscInt its; /* number of iterations so far computed */
138: PetscInt n,nloc; /* problem dimensions (global, local) */
139: PetscReal nrma,nrmb; /* computed matrix norms */
140: PetscBool useds; /* whether the solver uses the DS object or not */
141: PetscBool isgeneralized;
142: PetscBool ispositive;
143: PetscBool ishermitian;
144: EPSConvergedReason reason;
145: };
147: /*
148: Macros to test valid EPS arguments
149: */
150: #if !defined(PETSC_USE_DEBUG)
152: #define EPSCheckSolved(h,arg) do {} while (0)
154: #else
156: #define EPSCheckSolved(h,arg) \
157: do { \
158: if ((h)->state<EPS_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSolve() first: Parameter #%d",arg); \
159: } while (0)
161: #endif
163: /*
164: Macros to check settings at EPSSetUp()
165: */
167: /* EPSCheckHermitianDefinite: the problem is HEP or GHEP */
168: #define EPSCheckHermitianDefiniteCondition(eps,condition,msg) \
169: do { \
170: if (condition) { \
171: if (!(eps)->ishermitian) SETERRQ3(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s cannot be used for non-%s problems",((PetscObject)(eps))->type_name,(msg),SLEPC_STRING_HERMITIAN); \
172: else if ((eps)->isgeneralized && !(eps)->ispositive) SETERRQ3(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s requires that the problem is %s-definite",((PetscObject)(eps))->type_name,(msg),SLEPC_STRING_HERMITIAN); \
173: } \
174: } while (0)
175: #define EPSCheckHermitianDefinite(eps) EPSCheckHermitianDefiniteCondition(eps,PETSC_TRUE,"")
177: /* EPSCheckHermitian: the problem is HEP, GHEP, or GHIEP */
178: #define EPSCheckHermitianCondition(eps,condition,msg) \
179: do { \
180: if (condition) { \
181: if (!(eps)->ishermitian) SETERRQ3(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s cannot be used for non-%s problems",((PetscObject)(eps))->type_name,(msg),SLEPC_STRING_HERMITIAN); \
182: } \
183: } while (0)
184: #define EPSCheckHermitian(eps) EPSCheckHermitianCondition(eps,PETSC_TRUE,"")
186: /* EPSCheckDefinite: the problem is not GHIEP */
187: #define EPSCheckDefiniteCondition(eps,condition,msg) \
188: do { \
189: if (condition) { \
190: if ((eps)->isgeneralized && (eps)->ishermitian && !(eps)->ispositive) SETERRQ3(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s cannot be used for %s-indefinite problems",((PetscObject)(eps))->type_name,(msg),SLEPC_STRING_HERMITIAN); \
191: } \
192: } while (0)
193: #define EPSCheckDefinite(eps) EPSCheckDefiniteCondition(eps,PETSC_TRUE,"")
195: /* EPSCheckStandard: the problem is HEP or NHEP */
196: #define EPSCheckStandardCondition(eps,condition,msg) \
197: do { \
198: if (condition) { \
199: if ((eps)->isgeneralized) SETERRQ2(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s cannot be used for generalized problems",((PetscObject)(eps))->type_name,(msg)); \
200: } \
201: } while (0)
202: #define EPSCheckStandard(eps) EPSCheckStandardCondition(eps,PETSC_TRUE,"")
204: /* EPSCheckSinvert: shift-and-invert ST */
205: #define EPSCheckSinvertCondition(eps,condition,msg) \
206: do { \
207: if (condition) { \
208: PetscBool __flg; \
209: PetscObjectTypeCompare((PetscObject)(eps)->st,STSINVERT,&__flg); \
210: if (!__flg) SETERRQ2(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s requires a shift-and-invert spectral transform",((PetscObject)(eps))->type_name,(msg)); \
211: } \
212: } while (0)
213: #define EPSCheckSinvert(eps) EPSCheckSinvertCondition(eps,PETSC_TRUE,"")
215: /* EPSCheckSinvertCayley: shift-and-invert or Cayley ST */
216: #define EPSCheckSinvertCayleyCondition(eps,condition,msg) \
217: do { \
218: if (condition) { \
219: PetscBool __flg; \
220: PetscObjectTypeCompareAny((PetscObject)(eps)->st,&__flg,STSINVERT,STCAYLEY,""); \
221: if (!__flg) SETERRQ2(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s requires shift-and-invert or Cayley transform",((PetscObject)(eps))->type_name,(msg)); \
222: } \
223: } while (0)
224: #define EPSCheckSinvertCayley(eps) EPSCheckSinvertCayleyCondition(eps,PETSC_TRUE,"")
226: /* Check for unsupported features */
227: #define EPSCheckUnsupportedCondition(eps,mask,condition,msg) \
228: do { \
229: if (condition) { \
230: if (((mask) & EPS_FEATURE_BALANCE) && (eps)->balance!=EPS_BALANCE_NONE) SETERRQ2(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s does not support balancing",((PetscObject)(eps))->type_name,(msg)); \
231: if (((mask) & EPS_FEATURE_ARBITRARY) && (eps)->arbitrary) SETERRQ2(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s does not support arbitrary selection of eigenpairs",((PetscObject)(eps))->type_name,(msg)); \
232: if ((mask) & EPS_FEATURE_REGION) { \
233: PetscBool __istrivial; \
234: PetscErrorCode __RGIsTrivial((eps)->rg,&__istrivial);CHKERRQ(__ierr); \
235: if (!__istrivial) SETERRQ2(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s does not support region filtering",((PetscObject)(eps))->type_name,(msg)); \
236: } \
237: if (((mask) & EPS_FEATURE_EXTRACTION) && (eps)->extraction!=EPS_RITZ) SETERRQ2(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s only supports Ritz extraction",((PetscObject)(eps))->type_name,(msg)); \
238: if (((mask) & EPS_FEATURE_CONVERGENCE) && (eps)->converged!=EPSConvergedRelative) SETERRQ2(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default convergence test",((PetscObject)(eps))->type_name,(msg)); \
239: if (((mask) & EPS_FEATURE_STOPPING) && (eps)->stopping!=EPSStoppingBasic) SETERRQ2(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default stopping test",((PetscObject)(eps))->type_name,(msg)); \
240: if (((mask) & EPS_FEATURE_TWOSIDED) && (eps)->twosided) SETERRQ2(PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s cannot compute left eigenvectors (no two-sided variant)",((PetscObject)(eps))->type_name,(msg)); \
241: } \
242: } while (0)
243: #define EPSCheckUnsupported(eps,mask) EPSCheckUnsupportedCondition(eps,mask,PETSC_TRUE,"")
245: /* Check for ignored features */
246: #define EPSCheckIgnoredCondition(eps,mask,condition,msg) \
247: do { \
248: PetscErrorCode __ierr; \
249: if (condition) { \
250: if (((mask) & EPS_FEATURE_BALANCE) && (eps)->balance!=EPS_BALANCE_NONE) { __PetscInfo2((eps),"The solver '%s'%s ignores the balancing settings\n",((PetscObject)(eps))->type_name,(msg)); } \
251: if (((mask) & EPS_FEATURE_ARBITRARY) && (eps)->arbitrary) { __PetscInfo2((eps),"The solver '%s'%s ignores the settings for arbitrary selection of eigenpairs\n",((PetscObject)(eps))->type_name,(msg)); } \
252: if ((mask) & EPS_FEATURE_REGION) { \
253: PetscBool __istrivial; \
254: __RGIsTrivial((eps)->rg,&__istrivial);CHKERRQ(__ierr); \
255: if (!__istrivial) { __PetscInfo2((eps),"The solver '%s'%s ignores the specified region\n",((PetscObject)(eps))->type_name,(msg)); } \
256: } \
257: if (((mask) & EPS_FEATURE_EXTRACTION) && (eps)->extraction!=EPS_RITZ) { __PetscInfo2((eps),"The solver '%s'%s ignores the extraction settings\n",((PetscObject)(eps))->type_name,(msg)); } \
258: if (((mask) & EPS_FEATURE_CONVERGENCE) && (eps)->converged!=EPSConvergedRelative) { __PetscInfo2((eps),"The solver '%s'%s ignores the convergence test settings\n",((PetscObject)(eps))->type_name,(msg)); } \
259: if (((mask) & EPS_FEATURE_STOPPING) && (eps)->stopping!=EPSStoppingBasic) { __PetscInfo2((eps),"The solver '%s'%s ignores the stopping test settings\n",((PetscObject)(eps))->type_name,(msg)); } \
260: if (((mask) & EPS_FEATURE_TWOSIDED) && (eps)->twosided) { __PetscInfo2((eps),"The solver '%s'%s ignores the two-sided flag\n",((PetscObject)(eps))->type_name,(msg)); } \
261: } \
262: } while (0)
263: #define EPSCheckIgnored(eps,mask) EPSCheckIgnoredCondition(eps,mask,PETSC_TRUE,"")
265: /*
266: EPS_SetInnerProduct - set B matrix for inner product if appropriate.
267: */
268: PETSC_STATIC_INLINE PetscErrorCode EPS_SetInnerProduct(EPS eps)
269: {
271: Mat B;
274: if (!eps->V) { EPSGetBV(eps,&eps->V); }
275: if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
276: STGetBilinearForm(eps->st,&B);
277: BVSetMatrix(eps->V,B,PetscNot(eps->ispositive));
278: MatDestroy(&B);
279: } else {
280: BVSetMatrix(eps->V,NULL,PETSC_FALSE);
281: }
282: return(0);
283: }
285: /*
286: EPS_Purify - purify the first k vectors in the V basis
287: */
288: PETSC_STATIC_INLINE PetscErrorCode EPS_Purify(EPS eps,PetscInt k)
289: {
291: PetscInt i;
292: Vec v,z;
295: BVCreateVec(eps->V,&v);
296: for (i=0;i<k;i++) {
297: BVCopyVec(eps->V,i,v);
298: BVGetColumn(eps->V,i,&z);
299: STApply(eps->st,v,z);
300: BVRestoreColumn(eps->V,i,&z);
301: }
302: VecDestroy(&v);
303: return(0);
304: }
306: SLEPC_INTERN PetscErrorCode EPSSetWhichEigenpairs_Default(EPS);
307: SLEPC_INTERN PetscErrorCode EPSSetDimensions_Default(EPS,PetscInt,PetscInt*,PetscInt*);
308: SLEPC_INTERN PetscErrorCode EPSBackTransform_Default(EPS);
309: SLEPC_INTERN PetscErrorCode EPSComputeVectors(EPS);
310: SLEPC_INTERN PetscErrorCode EPSComputeVectors_Hermitian(EPS);
311: SLEPC_INTERN PetscErrorCode EPSComputeVectors_Schur(EPS);
312: SLEPC_INTERN PetscErrorCode EPSComputeVectors_Indefinite(EPS);
313: SLEPC_INTERN PetscErrorCode EPSComputeVectors_Twosided(EPS);
314: SLEPC_INTERN PetscErrorCode EPSComputeVectors_Slice(EPS);
315: SLEPC_INTERN PetscErrorCode EPSComputeResidualNorm_Private(EPS,PetscBool,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*);
316: SLEPC_INTERN PetscErrorCode EPSComputeRitzVector(EPS,PetscScalar*,PetscScalar*,BV,Vec,Vec);
317: SLEPC_INTERN PetscErrorCode EPSGetStartVector(EPS,PetscInt,PetscBool*);
318: SLEPC_INTERN PetscErrorCode EPSGetLeftStartVector(EPS,PetscInt,PetscBool*);
320: /* Private functions of the solver implementations */
322: SLEPC_INTERN PetscErrorCode EPSDelayedArnoldi(EPS,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
323: SLEPC_INTERN PetscErrorCode EPSDelayedArnoldi1(EPS,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
324: SLEPC_INTERN PetscErrorCode EPSKrylovConvergence(EPS,PetscBool,PetscInt,PetscInt,PetscReal,PetscReal,PetscReal,PetscInt*);
325: SLEPC_INTERN PetscErrorCode EPSPseudoLanczos(EPS,PetscReal*,PetscReal*,PetscReal*,PetscInt,PetscInt*,PetscBool*,PetscBool*,PetscReal*,Vec);
326: SLEPC_INTERN PetscErrorCode EPSBuildBalance_Krylov(EPS);
327: SLEPC_INTERN PetscErrorCode EPSSetDefaultST(EPS);
328: SLEPC_INTERN PetscErrorCode EPSSetDefaultST_Precond(EPS);
329: SLEPC_INTERN PetscErrorCode EPSSetDefaultST_GMRES(EPS);
330: SLEPC_INTERN PetscErrorCode EPSSetDefaultST_NoFactor(EPS);
331: SLEPC_INTERN PetscErrorCode EPSSetUpSort_Basic(EPS);
332: SLEPC_INTERN PetscErrorCode EPSSetUpSort_Default(EPS);
334: #endif