Actual source code: epsimpl.h

slepc-3.19.2 2023-09-05
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: */

 11: #if !defined(SLEPCEPSIMPL_H)
 12: #define SLEPCEPSIMPL_H

 14: #include <slepceps.h>
 15: #include <slepc/private/slepcimpl.h>

 17: /* SUBMANSEC = EPS */

 19: SLEPC_EXTERN PetscBool EPSRegisterAllCalled;
 20: SLEPC_EXTERN PetscBool EPSMonitorRegisterAllCalled;
 21: SLEPC_EXTERN PetscErrorCode EPSRegisterAll(void);
 22: SLEPC_EXTERN PetscErrorCode EPSMonitorRegisterAll(void);
 23: SLEPC_EXTERN PetscLogEvent EPS_SetUp,EPS_Solve,EPS_CISS_SVD;

 25: typedef struct _EPSOps *EPSOps;

 27: struct _EPSOps {
 28:   PetscErrorCode (*solve)(EPS);
 29:   PetscErrorCode (*setup)(EPS);
 30:   PetscErrorCode (*setupsort)(EPS);
 31:   PetscErrorCode (*setfromoptions)(EPS,PetscOptionItems*);
 32:   PetscErrorCode (*publishoptions)(EPS);
 33:   PetscErrorCode (*destroy)(EPS);
 34:   PetscErrorCode (*reset)(EPS);
 35:   PetscErrorCode (*view)(EPS,PetscViewer);
 36:   PetscErrorCode (*backtransform)(EPS);
 37:   PetscErrorCode (*computevectors)(EPS);
 38:   PetscErrorCode (*setdefaultst)(EPS);
 39:   PetscErrorCode (*setdstype)(EPS);
 40: };

 42: /*
 43:    Maximum number of monitors you can run with a single EPS
 44: */
 45: #define MAXEPSMONITORS 5

 47: /*
 48:    The solution process goes through several states
 49: */
 50: typedef enum { EPS_STATE_INITIAL,
 51:                EPS_STATE_SETUP,
 52:                EPS_STATE_SOLVED,
 53:                EPS_STATE_EIGENVECTORS } EPSStateType;

 55: /*
 56:    To classify the different solvers into categories
 57: */
 58: typedef enum { EPS_CATEGORY_KRYLOV,      /* Krylov solver: relies on STApply and STBackTransform (same as OTHER) */
 59:                EPS_CATEGORY_PRECOND,     /* Preconditioned solver: uses ST only to manage preconditioner */
 60:                EPS_CATEGORY_CONTOUR,     /* Contour integral: ST used to solve linear systems at integration points */
 61:                EPS_CATEGORY_OTHER } EPSSolverType;

 63: /*
 64:    To check for unsupported features at EPSSetUp_XXX()
 65: */
 66: typedef enum { EPS_FEATURE_BALANCE=1,       /* balancing */
 67:                EPS_FEATURE_ARBITRARY=2,     /* arbitrary selection of eigepairs */
 68:                EPS_FEATURE_REGION=4,        /* nontrivial region for filtering */
 69:                EPS_FEATURE_EXTRACTION=8,    /* extraction technique different from Ritz */
 70:                EPS_FEATURE_CONVERGENCE=16,  /* convergence test selected by user */
 71:                EPS_FEATURE_STOPPING=32,     /* stopping test */
 72:                EPS_FEATURE_TWOSIDED=64      /* two-sided variant */
 73:              } EPSFeatureType;

 75: /*
 76:    Defines the EPS data structure
 77: */
 78: struct _p_EPS {
 79:   PETSCHEADER(struct _EPSOps);
 80:   /*------------------------- User parameters ---------------------------*/
 81:   PetscInt       max_it;           /* maximum number of iterations */
 82:   PetscInt       nev;              /* number of eigenvalues to compute */
 83:   PetscInt       ncv;              /* number of basis vectors */
 84:   PetscInt       mpd;              /* maximum dimension of projected problem */
 85:   PetscInt       nini,ninil;       /* number of initial vectors (negative means not copied yet) */
 86:   PetscInt       nds;              /* number of basis vectors of deflation space */
 87:   PetscScalar    target;           /* target value */
 88:   PetscReal      tol;              /* tolerance */
 89:   EPSConv        conv;             /* convergence test */
 90:   EPSStop        stop;             /* stopping test */
 91:   EPSWhich       which;            /* which part of the spectrum to be sought */
 92:   PetscReal      inta,intb;        /* interval [a,b] for spectrum slicing */
 93:   EPSProblemType problem_type;     /* which kind of problem to be solved */
 94:   EPSExtraction  extraction;       /* which kind of extraction to be applied */
 95:   EPSBalance     balance;          /* the balancing method */
 96:   PetscInt       balance_its;      /* number of iterations of the balancing method */
 97:   PetscReal      balance_cutoff;   /* cutoff value for balancing */
 98:   PetscBool      trueres;          /* whether the true residual norm must be computed */
 99:   PetscBool      trackall;         /* whether all the residuals must be computed */
100:   PetscBool      purify;           /* whether eigenvectors need to be purified */
101:   PetscBool      twosided;         /* whether to compute left eigenvectors (two-sided solver) */

103:   /*-------------- User-provided functions and contexts -----------------*/
104:   PetscErrorCode (*converged)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
105:   PetscErrorCode (*convergeduser)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
106:   PetscErrorCode (*convergeddestroy)(void*);
107:   PetscErrorCode (*stopping)(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);
108:   PetscErrorCode (*stoppinguser)(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);
109:   PetscErrorCode (*stoppingdestroy)(void*);
110:   PetscErrorCode (*arbitrary)(PetscScalar,PetscScalar,Vec,Vec,PetscScalar*,PetscScalar*,void*);
111:   void           *convergedctx;
112:   void           *stoppingctx;
113:   void           *arbitraryctx;
114:   PetscErrorCode (*monitor[MAXEPSMONITORS])(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
115:   PetscErrorCode (*monitordestroy[MAXEPSMONITORS])(void**);
116:   void           *monitorcontext[MAXEPSMONITORS];
117:   PetscInt       numbermonitors;

119:   /*----------------- Child objects and working data -------------------*/
120:   ST             st;               /* spectral transformation object */
121:   DS             ds;               /* direct solver object */
122:   BV             V;                /* set of basis vectors and computed eigenvectors */
123:   BV             W;                /* left basis vectors (if left eigenvectors requested) */
124:   RG             rg;               /* optional region for filtering */
125:   SlepcSC        sc;               /* sorting criterion data */
126:   Vec            D;                /* diagonal matrix for balancing */
127:   Vec            *IS,*ISL;         /* references to user-provided initial spaces */
128:   Vec            *defl;            /* references to user-provided deflation space */
129:   PetscScalar    *eigr,*eigi;      /* real and imaginary parts of eigenvalues */
130:   PetscReal      *errest;          /* error estimates */
131:   PetscScalar    *rr,*ri;          /* values computed by user's arbitrary selection function */
132:   PetscInt       *perm;            /* permutation for eigenvalue ordering */
133:   PetscInt       nwork;            /* number of work vectors */
134:   Vec            *work;            /* work vectors */
135:   void           *data;            /* placeholder for solver-specific stuff */

137:   /* ----------------------- Status variables --------------------------*/
138:   EPSStateType   state;            /* initial -> setup -> solved -> eigenvectors */
139:   EPSSolverType  categ;            /* solver category */
140:   PetscInt       nconv;            /* number of converged eigenvalues */
141:   PetscInt       its;              /* number of iterations so far computed */
142:   PetscInt       n,nloc;           /* problem dimensions (global, local) */
143:   PetscReal      nrma,nrmb;        /* computed matrix norms */
144:   PetscBool      useds;            /* whether the solver uses the DS object or not */
145:   PetscBool      isgeneralized;
146:   PetscBool      ispositive;
147:   PetscBool      ishermitian;
148:   EPSConvergedReason reason;
149: };

151: /*
152:     Macros to test valid EPS arguments
153: */
154: #if !defined(PETSC_USE_DEBUG)

156: #define EPSCheckSolved(h,arg) do {(void)(h);} while (0)

158: #else

160: #define EPSCheckSolved(h,arg) \
161:   do { \
162:     PetscCheck((h)->state>=EPS_STATE_SOLVED,PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSolve() first: Parameter #%d",arg); \
163:   } while (0)

165: #endif

167: /*
168:     Macros to check settings at EPSSetUp()
169: */

171: /* EPSCheckHermitianDefinite: the problem is HEP or GHEP */
172: #define EPSCheckHermitianDefiniteCondition(eps,condition,msg) \
173:   do { \
174:     if (condition) { \
175:       PetscCheck((eps)->ishermitian,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); \
176:       PetscCheck(!(eps)->isgeneralized || (eps)->ispositive,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); \
177:     } \
178:   } while (0)
179: #define EPSCheckHermitianDefinite(eps) EPSCheckHermitianDefiniteCondition(eps,PETSC_TRUE,"")

181: /* EPSCheckHermitian: the problem is HEP, GHEP, or GHIEP */
182: #define EPSCheckHermitianCondition(eps,condition,msg) \
183:   do { \
184:     if (condition) { \
185:       PetscCheck((eps)->ishermitian,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); \
186:     } \
187:   } while (0)
188: #define EPSCheckHermitian(eps) EPSCheckHermitianCondition(eps,PETSC_TRUE,"")

190: /* EPSCheckDefinite: the problem is not GHIEP */
191: #define EPSCheckDefiniteCondition(eps,condition,msg) \
192:   do { \
193:     if (condition) { \
194:       PetscCheck(!(eps)->isgeneralized || !(eps)->ishermitian || (eps)->ispositive,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); \
195:     } \
196:   } while (0)
197: #define EPSCheckDefinite(eps) EPSCheckDefiniteCondition(eps,PETSC_TRUE,"")

199: /* EPSCheckStandard: the problem is HEP or NHEP */
200: #define EPSCheckStandardCondition(eps,condition,msg) \
201:   do { \
202:     if (condition) { \
203:       PetscCheck(!(eps)->isgeneralized,PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s cannot be used for generalized problems",((PetscObject)(eps))->type_name,(msg)); \
204:     } \
205:   } while (0)
206: #define EPSCheckStandard(eps) EPSCheckStandardCondition(eps,PETSC_TRUE,"")

208: /* EPSCheckSinvert: shift-and-invert ST */
209: #define EPSCheckSinvertCondition(eps,condition,msg) \
210:   do { \
211:     if (condition) { \
212:       PetscBool __flg; \
213:       PetscCall(PetscObjectTypeCompare((PetscObject)(eps)->st,STSINVERT,&__flg)); \
214:       PetscCheck(__flg,PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s requires a shift-and-invert spectral transform",((PetscObject)(eps))->type_name,(msg)); \
215:     } \
216:   } while (0)
217: #define EPSCheckSinvert(eps) EPSCheckSinvertCondition(eps,PETSC_TRUE,"")

219: /* EPSCheckSinvertCayley: shift-and-invert or Cayley ST */
220: #define EPSCheckSinvertCayleyCondition(eps,condition,msg) \
221:   do { \
222:     if (condition) { \
223:       PetscBool __flg; \
224:       PetscCall(PetscObjectTypeCompareAny((PetscObject)(eps)->st,&__flg,STSINVERT,STCAYLEY,"")); \
225:       PetscCheck(__flg,PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s requires shift-and-invert or Cayley transform",((PetscObject)(eps))->type_name,(msg)); \
226:     } \
227:   } while (0)
228: #define EPSCheckSinvertCayley(eps) EPSCheckSinvertCayleyCondition(eps,PETSC_TRUE,"")

230: /* Check for unsupported features */
231: #define EPSCheckUnsupportedCondition(eps,mask,condition,msg) \
232:   do { \
233:     if (condition) { \
234:       PetscCheck(!((mask) & EPS_FEATURE_BALANCE) || (eps)->balance==EPS_BALANCE_NONE,PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s does not support balancing",((PetscObject)(eps))->type_name,(msg)); \
235:       PetscCheck(!((mask) & EPS_FEATURE_ARBITRARY) || !(eps)->arbitrary,PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s does not support arbitrary selection of eigenpairs",((PetscObject)(eps))->type_name,(msg)); \
236:       if ((mask) & EPS_FEATURE_REGION) { \
237:         PetscBool      __istrivial; \
238:         PetscCall(RGIsTrivial((eps)->rg,&__istrivial)); \
239:         PetscCheck(__istrivial,PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s does not support region filtering",((PetscObject)(eps))->type_name,(msg)); \
240:       } \
241:       PetscCheck(!((mask) & EPS_FEATURE_EXTRACTION) || (eps)->extraction==EPS_RITZ,PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s only supports Ritz extraction",((PetscObject)(eps))->type_name,(msg)); \
242:       PetscCheck(!((mask) & EPS_FEATURE_CONVERGENCE) || (eps)->converged==EPSConvergedRelative,PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default convergence test",((PetscObject)(eps))->type_name,(msg)); \
243:       PetscCheck(!((mask) & EPS_FEATURE_STOPPING) || (eps)->stopping==EPSStoppingBasic,PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default stopping test",((PetscObject)(eps))->type_name,(msg)); \
244:       PetscCheck(!((mask) & EPS_FEATURE_TWOSIDED) || !(eps)->twosided,PetscObjectComm((PetscObject)(eps)),PETSC_ERR_SUP,"The solver '%s'%s cannot compute left eigenvectors (no two-sided variant)",((PetscObject)(eps))->type_name,(msg)); \
245:     } \
246:   } while (0)
247: #define EPSCheckUnsupported(eps,mask) EPSCheckUnsupportedCondition(eps,mask,PETSC_TRUE,"")

249: /* Check for ignored features */
250: #define EPSCheckIgnoredCondition(eps,mask,condition,msg) \
251:   do { \
252:     if (condition) { \
253:       if (((mask) & EPS_FEATURE_BALANCE) && (eps)->balance!=EPS_BALANCE_NONE) PetscCall(PetscInfo((eps),"The solver '%s'%s ignores the balancing settings\n",((PetscObject)(eps))->type_name,(msg))); \
254:       if (((mask) & EPS_FEATURE_ARBITRARY) && (eps)->arbitrary) PetscCall(PetscInfo((eps),"The solver '%s'%s ignores the settings for arbitrary selection of eigenpairs\n",((PetscObject)(eps))->type_name,(msg))); \
255:       if ((mask) & EPS_FEATURE_REGION) { \
256:         PetscBool __istrivial; \
257:         PetscCall(RGIsTrivial((eps)->rg,&__istrivial)); \
258:         if (!__istrivial) PetscCall(PetscInfo((eps),"The solver '%s'%s ignores the specified region\n",((PetscObject)(eps))->type_name,(msg))); \
259:       } \
260:       if (((mask) & EPS_FEATURE_EXTRACTION) && (eps)->extraction!=EPS_RITZ) PetscCall(PetscInfo((eps),"The solver '%s'%s ignores the extraction settings\n",((PetscObject)(eps))->type_name,(msg))); \
261:       if (((mask) & EPS_FEATURE_CONVERGENCE) && (eps)->converged!=EPSConvergedRelative) PetscCall(PetscInfo((eps),"The solver '%s'%s ignores the convergence test settings\n",((PetscObject)(eps))->type_name,(msg))); \
262:       if (((mask) & EPS_FEATURE_STOPPING) && (eps)->stopping!=EPSStoppingBasic) PetscCall(PetscInfo((eps),"The solver '%s'%s ignores the stopping test settings\n",((PetscObject)(eps))->type_name,(msg))); \
263:       if (((mask) & EPS_FEATURE_TWOSIDED) && (eps)->twosided) PetscCall(PetscInfo((eps),"The solver '%s'%s ignores the two-sided flag\n",((PetscObject)(eps))->type_name,(msg))); \
264:     } \
265:   } while (0)
266: #define EPSCheckIgnored(eps,mask) EPSCheckIgnoredCondition(eps,mask,PETSC_TRUE,"")

268: /*
269:   EPS_SetInnerProduct - set B matrix for inner product if appropriate.
270: */
271: static inline PetscErrorCode EPS_SetInnerProduct(EPS eps)
272: {
273:   Mat            B;

275:   PetscFunctionBegin;
276:   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
277:   if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
278:     PetscCall(STGetBilinearForm(eps->st,&B));
279:     PetscCall(BVSetMatrix(eps->V,B,PetscNot(eps->ispositive)));
280:     PetscCall(MatDestroy(&B));
281:   } else PetscCall(BVSetMatrix(eps->V,NULL,PETSC_FALSE));
282:   PetscFunctionReturn(PETSC_SUCCESS);
283: }

285: /*
286:   EPS_Purify - purify the first k vectors in the V basis
287: */
288: static inline PetscErrorCode EPS_Purify(EPS eps,PetscInt k)
289: {
290:   PetscInt       i;
291:   Vec            v,z;

293:   PetscFunctionBegin;
294:   PetscCall(BVCreateVec(eps->V,&v));
295:   for (i=0;i<k;i++) {
296:     PetscCall(BVCopyVec(eps->V,i,v));
297:     PetscCall(BVGetColumn(eps->V,i,&z));
298:     PetscCall(STApply(eps->st,v,z));
299:     PetscCall(BVRestoreColumn(eps->V,i,&z));
300:   }
301:   PetscCall(VecDestroy(&v));
302:   PetscFunctionReturn(PETSC_SUCCESS);
303: }

305: /*
306:   EPS_KSPSetOperators - Sets the KSP matrices, see also ST_KSPSetOperators()
307: */
308: static inline PetscErrorCode EPS_KSPSetOperators(KSP ksp,Mat A,Mat B)
309: {
310:   const char     *prefix;

312:   PetscFunctionBegin;
313:   PetscCall(KSPSetOperators(ksp,A,B));
314:   PetscCall(MatGetOptionsPrefix(B,&prefix));
315:   if (!prefix) {
316:     /* set Mat prefix to be the same as KSP to enable setting command-line options (e.g. MUMPS)
317:        only applies if the Mat has no user-defined prefix */
318:     PetscCall(KSPGetOptionsPrefix(ksp,&prefix));
319:     PetscCall(MatSetOptionsPrefix(B,prefix));
320:   }
321:   PetscFunctionReturn(PETSC_SUCCESS);
322: }

324: SLEPC_INTERN PetscErrorCode EPSSetWhichEigenpairs_Default(EPS);
325: SLEPC_INTERN PetscErrorCode EPSSetDimensions_Default(EPS,PetscInt,PetscInt*,PetscInt*);
326: SLEPC_INTERN PetscErrorCode EPSBackTransform_Default(EPS);
327: SLEPC_INTERN PetscErrorCode EPSComputeVectors(EPS);
328: SLEPC_INTERN PetscErrorCode EPSComputeVectors_Hermitian(EPS);
329: SLEPC_INTERN PetscErrorCode EPSComputeVectors_Schur(EPS);
330: SLEPC_INTERN PetscErrorCode EPSComputeVectors_Indefinite(EPS);
331: SLEPC_INTERN PetscErrorCode EPSComputeVectors_Twosided(EPS);
332: SLEPC_INTERN PetscErrorCode EPSComputeVectors_Slice(EPS);
333: SLEPC_INTERN PetscErrorCode EPSComputeResidualNorm_Private(EPS,PetscBool,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*);
334: SLEPC_INTERN PetscErrorCode EPSComputeRitzVector(EPS,PetscScalar*,PetscScalar*,BV,Vec,Vec);
335: SLEPC_INTERN PetscErrorCode EPSGetStartVector(EPS,PetscInt,PetscBool*);
336: SLEPC_INTERN PetscErrorCode EPSGetLeftStartVector(EPS,PetscInt,PetscBool*);
337: SLEPC_INTERN PetscErrorCode MatEstimateSpectralRange_EPS(Mat,PetscReal*,PetscReal*);

339: /* Private functions of the solver implementations */

341: SLEPC_INTERN PetscErrorCode EPSDelayedArnoldi(EPS,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
342: SLEPC_INTERN PetscErrorCode EPSDelayedArnoldi1(EPS,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
343: SLEPC_INTERN PetscErrorCode EPSKrylovConvergence(EPS,PetscBool,PetscInt,PetscInt,PetscReal,PetscReal,PetscReal,PetscInt*);
344: SLEPC_INTERN PetscErrorCode EPSPseudoLanczos(EPS,PetscReal*,PetscReal*,PetscReal*,PetscInt,PetscInt*,PetscBool*,PetscBool*,PetscReal*,Vec);
345: SLEPC_INTERN PetscErrorCode EPSBuildBalance_Krylov(EPS);
346: SLEPC_INTERN PetscErrorCode EPSSetDefaultST(EPS);
347: SLEPC_INTERN PetscErrorCode EPSSetDefaultST_Precond(EPS);
348: SLEPC_INTERN PetscErrorCode EPSSetDefaultST_GMRES(EPS);
349: SLEPC_INTERN PetscErrorCode EPSSetDefaultST_NoFactor(EPS);
350: SLEPC_INTERN PetscErrorCode EPSSetUpSort_Basic(EPS);
351: SLEPC_INTERN PetscErrorCode EPSSetUpSort_Default(EPS);

353: #endif