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