Actual source code: pepimpl.h
slepc-3.19.0 2023-03-31
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(SLEPCPEPIMPL_H)
12: #define SLEPCPEPIMPL_H
14: #include <slepcpep.h>
15: #include <slepc/private/slepcimpl.h>
17: /* SUBMANSEC = PEP */
19: SLEPC_EXTERN PetscBool PEPRegisterAllCalled;
20: SLEPC_EXTERN PetscBool PEPMonitorRegisterAllCalled;
21: SLEPC_EXTERN PetscErrorCode PEPRegisterAll(void);
22: SLEPC_EXTERN PetscErrorCode PEPMonitorRegisterAll(void);
23: SLEPC_EXTERN PetscLogEvent PEP_SetUp,PEP_Solve,PEP_Refine,PEP_CISS_SVD;
25: typedef struct _PEPOps *PEPOps;
27: struct _PEPOps {
28: PetscErrorCode (*solve)(PEP);
29: PetscErrorCode (*setup)(PEP);
30: PetscErrorCode (*setfromoptions)(PEP,PetscOptionItems*);
31: PetscErrorCode (*publishoptions)(PEP);
32: PetscErrorCode (*destroy)(PEP);
33: PetscErrorCode (*reset)(PEP);
34: PetscErrorCode (*view)(PEP,PetscViewer);
35: PetscErrorCode (*backtransform)(PEP);
36: PetscErrorCode (*computevectors)(PEP);
37: PetscErrorCode (*extractvectors)(PEP);
38: PetscErrorCode (*setdefaultst)(PEP);
39: PetscErrorCode (*setdstype)(PEP);
40: };
42: /*
43: Maximum number of monitors you can run with a single PEP
44: */
45: #define MAXPEPMONITORS 5
47: typedef enum { PEP_STATE_INITIAL,
48: PEP_STATE_SETUP,
49: PEP_STATE_SOLVED,
50: PEP_STATE_EIGENVECTORS } PEPStateType;
52: /*
53: To check for unsupported features at PEPSetUp_XXX()
54: */
55: typedef enum { PEP_FEATURE_NONMONOMIAL=1, /* non-monomial bases */
56: PEP_FEATURE_REGION=4, /* nontrivial region for filtering */
57: PEP_FEATURE_EXTRACT=8, /* eigenvector extraction */
58: PEP_FEATURE_CONVERGENCE=16, /* convergence test selected by user */
59: PEP_FEATURE_STOPPING=32, /* stopping test */
60: PEP_FEATURE_SCALE=64 /* scaling */
61: } PEPFeatureType;
63: /*
64: Defines the PEP data structure.
65: */
66: struct _p_PEP {
67: PETSCHEADER(struct _PEPOps);
68: /*------------------------- User parameters ---------------------------*/
69: PetscInt max_it; /* maximum number of iterations */
70: PetscInt nev; /* number of eigenvalues to compute */
71: PetscInt ncv; /* number of basis vectors */
72: PetscInt mpd; /* maximum dimension of projected problem */
73: PetscInt nini; /* number of initial vectors (negative means not copied yet) */
74: PetscScalar target; /* target value */
75: PetscReal tol; /* tolerance */
76: PEPConv conv; /* convergence test */
77: PEPStop stop; /* stopping test */
78: PEPWhich which; /* which part of the spectrum to be sought */
79: PetscReal inta,intb; /* interval [a,b] for spectrum slicing */
80: PEPBasis basis; /* polynomial basis used to represent the problem */
81: PEPProblemType problem_type; /* which kind of problem to be solved */
82: PEPScale scale; /* scaling strategy to be used */
83: PetscReal sfactor,dsfactor; /* scaling factors */
84: PetscInt sits; /* number of iterations of the scaling method */
85: PetscReal slambda; /* norm eigenvalue approximation for scaling */
86: PEPRefine refine; /* type of refinement to be applied after solve */
87: PetscInt npart; /* number of partitions of the communicator */
88: PetscReal rtol; /* tolerance for refinement */
89: PetscInt rits; /* number of iterations of the refinement method */
90: PEPRefineScheme scheme; /* scheme for solving linear systems within refinement */
91: PEPExtract extract; /* type of extraction used */
92: PetscBool trackall; /* whether all the residuals must be computed */
94: /*-------------- User-provided functions and contexts -----------------*/
95: PetscErrorCode (*converged)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
96: PetscErrorCode (*convergeduser)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
97: PetscErrorCode (*convergeddestroy)(void*);
98: PetscErrorCode (*stopping)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*);
99: PetscErrorCode (*stoppinguser)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*);
100: PetscErrorCode (*stoppingdestroy)(void*);
101: void *convergedctx;
102: void *stoppingctx;
103: PetscErrorCode (*monitor[MAXPEPMONITORS])(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
104: PetscErrorCode (*monitordestroy[MAXPEPMONITORS])(void**);
105: void *monitorcontext[MAXPEPMONITORS];
106: PetscInt numbermonitors;
108: /*----------------- Child objects and working data -------------------*/
109: ST st; /* spectral transformation object */
110: DS ds; /* direct solver object */
111: BV V; /* set of basis vectors and computed eigenvectors */
112: RG rg; /* optional region for filtering */
113: SlepcSC sc; /* sorting criterion data */
114: Mat *A; /* coefficient matrices of the polynomial */
115: PetscInt nmat; /* number of matrices */
116: Vec Dl,Dr; /* diagonal matrices for balancing */
117: Vec *IS; /* references to user-provided initial space */
118: PetscScalar *eigr,*eigi; /* real and imaginary parts of eigenvalues */
119: PetscReal *errest; /* error estimates */
120: PetscInt *perm; /* permutation for eigenvalue ordering */
121: PetscReal *pbc; /* coefficients defining the polynomial basis */
122: PetscScalar *solvematcoeffs; /* coefficients to compute the matrix to be inverted */
123: PetscInt nwork; /* number of work vectors */
124: Vec *work; /* work vectors */
125: KSP refineksp; /* ksp used in refinement */
126: PetscSubcomm refinesubc; /* context for sub-communicators */
127: void *data; /* placeholder for solver-specific stuff */
129: /* ----------------------- Status variables --------------------------*/
130: PEPStateType state; /* initial -> setup -> solved -> eigenvectors */
131: PetscInt nconv; /* number of converged eigenvalues */
132: PetscInt its; /* number of iterations so far computed */
133: PetscInt n,nloc; /* problem dimensions (global, local) */
134: PetscReal *nrma; /* computed matrix norms */
135: PetscReal nrml[2]; /* computed matrix norms for the linearization */
136: PetscBool sfactor_set; /* flag to indicate the user gave sfactor */
137: PetscBool lineariz; /* current solver is based on linearization */
138: PEPConvergedReason reason;
139: };
141: /*
142: Macros to test valid PEP arguments
143: */
144: #if !defined(PETSC_USE_DEBUG)
146: #define PEPCheckSolved(h,arg) do {(void)(h);} while (0)
148: #else
150: #define PEPCheckSolved(h,arg) \
151: do { \
152: PetscCheck((h)->state>=PEP_STATE_SOLVED,PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"Must call PEPSolve() first: Parameter #%d",arg); \
153: } while (0)
155: #endif
157: /*
158: Macros to check settings at PEPSetUp()
159: */
161: /* PEPCheckHermitian: the problem is Hermitian or Hyperbolic */
162: #define PEPCheckHermitianCondition(pep,condition,msg) \
163: do { \
164: if (condition) { \
165: PetscCheck((pep)->problem_type==PEP_HERMITIAN || (pep)->problem_type==PEP_HYPERBOLIC,PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s can only be used for Hermitian (or hyperbolic) problems",((PetscObject)(pep))->type_name,(msg)); \
166: } \
167: } while (0)
168: #define PEPCheckHermitian(pep) PEPCheckHermitianCondition(pep,PETSC_TRUE,"")
170: /* PEPCheckQuadratic: the polynomial has degree 2 */
171: #define PEPCheckQuadraticCondition(pep,condition,msg) \
172: do { \
173: if (condition) { \
174: PetscCheck((pep)->nmat==3,PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s is only available for quadratic problems",((PetscObject)(pep))->type_name,(msg)); \
175: } \
176: } while (0)
177: #define PEPCheckQuadratic(pep) PEPCheckQuadraticCondition(pep,PETSC_TRUE,"")
179: /* PEPCheckShiftSinvert: shift or shift-and-invert ST */
180: #define PEPCheckShiftSinvertCondition(pep,condition,msg) \
181: do { \
182: if (condition) { \
183: PetscBool __flg; \
184: PetscCall(PetscObjectTypeCompareAny((PetscObject)(pep)->st,&__flg,STSINVERT,STSHIFT,"")); \
185: PetscCheck(__flg,PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s requires shift or shift-and-invert spectral transform",((PetscObject)(pep))->type_name,(msg)); \
186: } \
187: } while (0)
188: #define PEPCheckShiftSinvert(pep) PEPCheckShiftSinvertCondition(pep,PETSC_TRUE,"")
190: /* PEPCheckSinvertCayley: shift-and-invert or Cayley ST */
191: #define PEPCheckSinvertCayleyCondition(pep,condition,msg) \
192: do { \
193: if (condition) { \
194: PetscBool __flg; \
195: PetscCall(PetscObjectTypeCompareAny((PetscObject)(pep)->st,&__flg,STSINVERT,STCAYLEY,"")); \
196: PetscCheck(__flg,PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s requires shift-and-invert or Cayley transform",((PetscObject)(pep))->type_name,(msg)); \
197: } \
198: } while (0)
199: #define PEPCheckSinvertCayley(pep) PEPCheckSinvertCayleyCondition(pep,PETSC_TRUE,"")
201: /* Check for unsupported features */
202: #define PEPCheckUnsupportedCondition(pep,mask,condition,msg) \
203: do { \
204: if (condition) { \
205: PetscCheck(!((mask) & PEP_FEATURE_NONMONOMIAL) || (pep)->basis==PEP_BASIS_MONOMIAL,PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s is not implemented for non-monomial bases",((PetscObject)(pep))->type_name,(msg)); \
206: if ((mask) & PEP_FEATURE_REGION) { \
207: PetscBool __istrivial; \
208: PetscCall(RGIsTrivial((pep)->rg,&__istrivial)); \
209: PetscCheck(__istrivial,PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s does not support region filtering",((PetscObject)(pep))->type_name,(msg)); \
210: } \
211: PetscCheck(!((mask) & PEP_FEATURE_EXTRACT) || !(pep)->extract || (pep)->extract==PEP_EXTRACT_NONE,PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s does not support extraction variants",((PetscObject)(pep))->type_name,(msg)); \
212: PetscCheck(!((mask) & PEP_FEATURE_CONVERGENCE) || (pep)->converged==PEPConvergedRelative,PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default convergence test",((PetscObject)(pep))->type_name,(msg)); \
213: PetscCheck(!((mask) & PEP_FEATURE_STOPPING) || (pep)->stopping==PEPStoppingBasic,PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default stopping test",((PetscObject)(pep))->type_name,(msg)); \
214: } \
215: } while (0)
216: #define PEPCheckUnsupported(pep,mask) PEPCheckUnsupportedCondition(pep,mask,PETSC_TRUE,"")
218: /* Check for ignored features */
219: #define PEPCheckIgnoredCondition(pep,mask,condition,msg) \
220: do { \
221: if (condition) { \
222: if (((mask) & PEP_FEATURE_NONMONOMIAL) && (pep)->basis!=PEP_BASIS_MONOMIAL) PetscCall(PetscInfo((pep),"The solver '%s'%s ignores the basis settings\n",((PetscObject)(pep))->type_name,(msg))); \
223: if ((mask) & PEP_FEATURE_REGION) { \
224: PetscBool __istrivial; \
225: PetscCall(RGIsTrivial((pep)->rg,&__istrivial)); \
226: if (!__istrivial) PetscCall(PetscInfo((pep),"The solver '%s'%s ignores the specified region\n",((PetscObject)(pep))->type_name,(msg))); \
227: } \
228: if (((mask) & PEP_FEATURE_EXTRACT) && (pep)->extract && (pep)->extract!=PEP_EXTRACT_NONE) PetscCall(PetscInfo((pep),"The solver '%s'%s ignores the extract settings\n",((PetscObject)(pep))->type_name,(msg))); \
229: if (((mask) & PEP_FEATURE_CONVERGENCE) && (pep)->converged!=PEPConvergedRelative) PetscCall(PetscInfo((pep),"The solver '%s'%s ignores the convergence test settings\n",((PetscObject)(pep))->type_name,(msg))); \
230: if (((mask) & PEP_FEATURE_STOPPING) && (pep)->stopping!=PEPStoppingBasic) PetscCall(PetscInfo((pep),"The solver '%s'%s ignores the stopping test settings\n",((PetscObject)(pep))->type_name,(msg))); \
231: if (((mask) & PEP_FEATURE_SCALE) && (pep)->scale!=PEP_SCALE_NONE) PetscCall(PetscInfo((pep),"The solver '%s'%s ignores the scaling settings\n",((PetscObject)(pep))->type_name,(msg))); \
232: } \
233: } while (0)
234: #define PEPCheckIgnored(pep,mask) PEPCheckIgnoredCondition(pep,mask,PETSC_TRUE,"")
236: /*
237: PEP_KSPSetOperators - Sets the KSP matrices
238: */
239: static inline PetscErrorCode PEP_KSPSetOperators(KSP ksp,Mat A,Mat B)
240: {
241: const char *prefix;
243: PetscFunctionBegin;
244: PetscCall(KSPSetOperators(ksp,A,B));
245: PetscCall(MatGetOptionsPrefix(B,&prefix));
246: if (!prefix) {
247: /* set Mat prefix to be the same as KSP to enable setting command-line options (e.g. MUMPS)
248: only applies if the Mat has no user-defined prefix */
249: PetscCall(KSPGetOptionsPrefix(ksp,&prefix));
250: PetscCall(MatSetOptionsPrefix(B,prefix));
251: }
252: PetscFunctionReturn(PETSC_SUCCESS);
253: }
255: SLEPC_INTERN PetscErrorCode PEPSetWhichEigenpairs_Default(PEP);
256: SLEPC_INTERN PetscErrorCode PEPSetDimensions_Default(PEP,PetscInt,PetscInt*,PetscInt*);
257: SLEPC_INTERN PetscErrorCode PEPExtractVectors(PEP);
258: SLEPC_INTERN PetscErrorCode PEPBackTransform_Default(PEP);
259: SLEPC_INTERN PetscErrorCode PEPComputeVectors(PEP);
260: SLEPC_INTERN PetscErrorCode PEPComputeVectors_Default(PEP);
261: SLEPC_INTERN PetscErrorCode PEPComputeVectors_Indefinite(PEP);
262: SLEPC_INTERN PetscErrorCode PEPComputeResidualNorm_Private(PEP,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*);
263: SLEPC_INTERN PetscErrorCode PEPKrylovConvergence(PEP,PetscBool,PetscInt,PetscInt,PetscReal,PetscInt*);
264: SLEPC_INTERN PetscErrorCode PEPComputeScaleFactor(PEP);
265: SLEPC_INTERN PetscErrorCode PEPBuildDiagonalScaling(PEP);
266: SLEPC_INTERN PetscErrorCode PEPBasisCoefficients(PEP,PetscReal*);
267: SLEPC_INTERN PetscErrorCode PEPEvaluateBasis(PEP,PetscScalar,PetscScalar,PetscScalar*,PetscScalar*);
268: SLEPC_INTERN PetscErrorCode PEPEvaluateBasisDerivative(PEP,PetscScalar,PetscScalar,PetscScalar*,PetscScalar*);
269: SLEPC_INTERN PetscErrorCode PEPEvaluateBasisMat(PEP,PetscInt,PetscScalar*,PetscInt,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt);
270: SLEPC_INTERN PetscErrorCode PEPNewtonRefinement_TOAR(PEP,PetscScalar,PetscInt*,PetscReal*,PetscInt,PetscScalar*,PetscInt);
271: SLEPC_INTERN PetscErrorCode PEPNewtonRefinementSimple(PEP,PetscInt*,PetscReal,PetscInt);
272: SLEPC_INTERN PetscErrorCode PEPSetDefaultST(PEP);
273: SLEPC_INTERN PetscErrorCode PEPSetDefaultST_Transform(PEP);
275: #endif