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