Actual source code: nepimpl.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(SLEPCNEPIMPL_H)
12: #define SLEPCNEPIMPL_H
14: #include <slepcnep.h>
15: #include <slepc/private/slepcimpl.h>
17: /* SUBMANSEC = NEP */
19: SLEPC_EXTERN PetscBool NEPRegisterAllCalled;
20: SLEPC_EXTERN PetscBool NEPMonitorRegisterAllCalled;
21: SLEPC_EXTERN PetscErrorCode NEPRegisterAll(void);
22: SLEPC_EXTERN PetscErrorCode NEPMonitorRegisterAll(void);
23: SLEPC_EXTERN PetscLogEvent NEP_SetUp,NEP_Solve,NEP_Refine,NEP_FunctionEval,NEP_JacobianEval,NEP_Resolvent,NEP_CISS_SVD;
25: typedef struct _NEPOps *NEPOps;
27: struct _NEPOps {
28: PetscErrorCode (*solve)(NEP);
29: PetscErrorCode (*setup)(NEP);
30: PetscErrorCode (*setfromoptions)(NEP,PetscOptionItems*);
31: PetscErrorCode (*publishoptions)(NEP);
32: PetscErrorCode (*destroy)(NEP);
33: PetscErrorCode (*reset)(NEP);
34: PetscErrorCode (*view)(NEP,PetscViewer);
35: PetscErrorCode (*computevectors)(NEP);
36: PetscErrorCode (*setdstype)(NEP);
37: };
39: /*
40: Maximum number of monitors you can run with a single NEP
41: */
42: #define MAXNEPMONITORS 5
44: typedef enum { NEP_STATE_INITIAL,
45: NEP_STATE_SETUP,
46: NEP_STATE_SOLVED,
47: NEP_STATE_EIGENVECTORS } NEPStateType;
49: /*
50: How the problem function T(lambda) has been defined by the user
51: - Callback: one callback to build the function matrix, another one for the Jacobian
52: - Split: in split form sum_j(A_j*f_j(lambda))
53: */
54: typedef enum { NEP_USER_INTERFACE_CALLBACK=1,
55: NEP_USER_INTERFACE_SPLIT } NEPUserInterface;
57: /*
58: To check for unsupported features at NEPSetUp_XXX()
59: */
60: typedef enum { NEP_FEATURE_CALLBACK=1, /* callback user interface */
61: NEP_FEATURE_REGION=4, /* nontrivial region for filtering */
62: NEP_FEATURE_CONVERGENCE=16, /* convergence test selected by user */
63: NEP_FEATURE_STOPPING=32, /* stopping test */
64: NEP_FEATURE_TWOSIDED=64 /* two-sided variant */
65: } NEPFeatureType;
67: /*
68: Defines the NEP data structure.
69: */
70: struct _p_NEP {
71: PETSCHEADER(struct _NEPOps);
72: /*------------------------- User parameters ---------------------------*/
73: PetscInt max_it; /* maximum number of iterations */
74: PetscInt nev; /* number of eigenvalues to compute */
75: PetscInt ncv; /* number of basis vectors */
76: PetscInt mpd; /* maximum dimension of projected problem */
77: PetscInt nini; /* number of initial vectors (negative means not copied yet) */
78: PetscScalar target; /* target value */
79: PetscReal tol; /* tolerance */
80: NEPConv conv; /* convergence test */
81: NEPStop stop; /* stopping test */
82: NEPWhich which; /* which part of the spectrum to be sought */
83: NEPProblemType problem_type; /* which kind of problem to be solved */
84: NEPRefine refine; /* type of refinement to be applied after solve */
85: PetscInt npart; /* number of partitions of the communicator */
86: PetscReal rtol; /* tolerance for refinement */
87: PetscInt rits; /* number of iterations of the refinement method */
88: NEPRefineScheme scheme; /* scheme for solving linear systems within refinement */
89: PetscBool trackall; /* whether all the residuals must be computed */
90: PetscBool twosided; /* whether to compute left eigenvectors (two-sided solver) */
92: /*-------------- User-provided functions and contexts -----------------*/
93: PetscErrorCode (*computefunction)(NEP,PetscScalar,Mat,Mat,void*);
94: PetscErrorCode (*computejacobian)(NEP,PetscScalar,Mat,void*);
95: void *functionctx;
96: void *jacobianctx;
97: PetscErrorCode (*converged)(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
98: PetscErrorCode (*convergeduser)(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
99: PetscErrorCode (*convergeddestroy)(void*);
100: PetscErrorCode (*stopping)(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*);
101: PetscErrorCode (*stoppinguser)(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*);
102: PetscErrorCode (*stoppingdestroy)(void*);
103: void *convergedctx;
104: void *stoppingctx;
105: PetscErrorCode (*monitor[MAXNEPMONITORS])(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
106: PetscErrorCode (*monitordestroy[MAXNEPMONITORS])(void**);
107: void *monitorcontext[MAXNEPMONITORS];
108: PetscInt numbermonitors;
110: /*----------------- Child objects and working data -------------------*/
111: DS ds; /* direct solver object */
112: BV V; /* set of basis vectors and computed eigenvectors */
113: BV W; /* left basis vectors (if left eigenvectors requested) */
114: RG rg; /* optional region for filtering */
115: SlepcSC sc; /* sorting criterion data */
116: Mat function; /* function matrix */
117: Mat function_pre; /* function matrix (preconditioner) */
118: Mat jacobian; /* Jacobian matrix */
119: Mat *A; /* matrix coefficients of split form */
120: FN *f; /* matrix functions of split form */
121: PetscInt nt; /* number of terms in split form */
122: MatStructure mstr; /* pattern of split matrices */
123: Mat *P; /* matrix coefficients of split form (preconditioner) */
124: MatStructure mstrp; /* pattern of split matrices (preconditioner) */
125: Vec *IS; /* references to user-provided initial space */
126: PetscScalar *eigr,*eigi; /* real and imaginary parts of eigenvalues */
127: PetscReal *errest; /* error estimates */
128: PetscInt *perm; /* permutation for eigenvalue ordering */
129: PetscInt nwork; /* number of work vectors */
130: Vec *work; /* work vectors */
131: KSP refineksp; /* ksp used in refinement */
132: PetscSubcomm refinesubc; /* context for sub-communicators */
133: void *data; /* placeholder for solver-specific stuff */
135: /* ----------------------- Status variables --------------------------*/
136: NEPStateType state; /* initial -> setup -> solved -> eigenvectors */
137: PetscInt nconv; /* number of converged eigenvalues */
138: PetscInt its; /* number of iterations so far computed */
139: PetscInt n,nloc; /* problem dimensions (global, local) */
140: PetscReal *nrma; /* computed matrix norms */
141: NEPUserInterface fui; /* how the user has defined the nonlinear operator */
142: PetscBool useds; /* whether the solver uses the DS object or not */
143: Mat resolvent; /* shell matrix to be used in NEPApplyResolvent */
144: NEPConvergedReason reason;
145: };
147: /*
148: Macros to test valid NEP arguments
149: */
150: #if !defined(PETSC_USE_DEBUG)
152: #define NEPCheckProblem(h,arg) do {(void)(h);} while (0)
153: #define NEPCheckCallback(h,arg) do {(void)(h);} while (0)
154: #define NEPCheckSplit(h,arg) do {(void)(h);} while (0)
155: #define NEPCheckDerivatives(h,arg) do {(void)(h);} while (0)
156: #define NEPCheckSolved(h,arg) do {(void)(h);} while (0)
158: #else
160: #define NEPCheckProblem(h,arg) \
161: do { \
162: PetscCheck(((h)->fui),PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"The nonlinear eigenproblem has not been specified yet. Parameter #%d",arg); \
163: } while (0)
165: #define NEPCheckCallback(h,arg) \
166: do { \
167: PetscCheck((h)->fui==NEP_USER_INTERFACE_CALLBACK,PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"This operation requires the nonlinear eigenproblem specified with callbacks. Parameter #%d",arg); \
168: } while (0)
170: #define NEPCheckSplit(h,arg) \
171: do { \
172: PetscCheck((h)->fui==NEP_USER_INTERFACE_SPLIT,PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"This operation requires the nonlinear eigenproblem in split form. Parameter #%d",arg); \
173: } while (0)
175: #define NEPCheckSolved(h,arg) \
176: do { \
177: PetscCheck((h)->state>=NEP_STATE_SOLVED,PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"Must call NEPSolve() first: Parameter #%d",arg); \
178: } while (0)
180: #endif
182: /* Check for unsupported features */
183: #define NEPCheckUnsupportedCondition(nep,mask,condition,msg) \
184: do { \
185: if (condition) { \
186: PetscCheck(!((mask) & NEP_FEATURE_CALLBACK) || (nep)->fui!=NEP_USER_INTERFACE_CALLBACK,PetscObjectComm((PetscObject)(nep)),PETSC_ERR_SUP,"The solver '%s'%s cannot be used with callback functions (use the split operator)",((PetscObject)(nep))->type_name,(msg)); \
187: if ((mask) & NEP_FEATURE_REGION) { \
188: PetscBool __istrivial; \
189: PetscCall(RGIsTrivial((nep)->rg,&__istrivial)); \
190: PetscCheck(__istrivial,PetscObjectComm((PetscObject)(nep)),PETSC_ERR_SUP,"The solver '%s'%s does not support region filtering",((PetscObject)(nep))->type_name,(msg)); \
191: } \
192: PetscCheck(!((mask) & NEP_FEATURE_CONVERGENCE) || (nep)->converged==NEPConvergedRelative,PetscObjectComm((PetscObject)(nep)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default convergence test",((PetscObject)(nep))->type_name,(msg)); \
193: PetscCheck(!((mask) & NEP_FEATURE_STOPPING) || (nep)->stopping==NEPStoppingBasic,PetscObjectComm((PetscObject)(nep)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default stopping test",((PetscObject)(nep))->type_name,(msg)); \
194: PetscCheck(!((mask) & NEP_FEATURE_TWOSIDED) || !(nep)->twosided,PetscObjectComm((PetscObject)(nep)),PETSC_ERR_SUP,"The solver '%s'%s cannot compute left eigenvectors (no two-sided variant)",((PetscObject)(nep))->type_name,(msg)); \
195: } \
196: } while (0)
197: #define NEPCheckUnsupported(nep,mask) NEPCheckUnsupportedCondition(nep,mask,PETSC_TRUE,"")
199: /* Check for ignored features */
200: #define NEPCheckIgnoredCondition(nep,mask,condition,msg) \
201: do { \
202: if (condition) { \
203: if (((mask) & NEP_FEATURE_CALLBACK) && (nep)->fui==NEP_USER_INTERFACE_CALLBACK) PetscCall(PetscInfo((nep),"The solver '%s'%s ignores the user interface settings\n",((PetscObject)(nep))->type_name,(msg))); \
204: if ((mask) & NEP_FEATURE_REGION) { \
205: PetscBool __istrivial; \
206: PetscCall(RGIsTrivial((nep)->rg,&__istrivial)); \
207: if (!__istrivial) PetscCall(PetscInfo((nep),"The solver '%s'%s ignores the specified region\n",((PetscObject)(nep))->type_name,(msg))); \
208: } \
209: if (((mask) & NEP_FEATURE_CONVERGENCE) && (nep)->converged!=NEPConvergedRelative) PetscCall(PetscInfo((nep),"The solver '%s'%s ignores the convergence test settings\n",((PetscObject)(nep))->type_name,(msg))); \
210: if (((mask) & NEP_FEATURE_STOPPING) && (nep)->stopping!=NEPStoppingBasic) PetscCall(PetscInfo((nep),"The solver '%s'%s ignores the stopping test settings\n",((PetscObject)(nep))->type_name,(msg))); \
211: if (((mask) & NEP_FEATURE_TWOSIDED) && (nep)->twosided) PetscCall(PetscInfo((nep),"The solver '%s'%s ignores the two-sided flag\n",((PetscObject)(nep))->type_name,(msg))); \
212: } \
213: } while (0)
214: #define NEPCheckIgnored(nep,mask) NEPCheckIgnoredCondition(nep,mask,PETSC_TRUE,"")
216: /*
217: NEP_KSPSetOperators - Sets the KSP matrices
218: */
219: static inline PetscErrorCode NEP_KSPSetOperators(KSP ksp,Mat A,Mat B)
220: {
221: const char *prefix;
223: PetscFunctionBegin;
224: PetscCall(KSPSetOperators(ksp,A,B));
225: PetscCall(MatGetOptionsPrefix(B,&prefix));
226: if (!prefix) {
227: /* set Mat prefix to be the same as KSP to enable setting command-line options (e.g. MUMPS)
228: only applies if the Mat has no user-defined prefix */
229: PetscCall(KSPGetOptionsPrefix(ksp,&prefix));
230: PetscCall(MatSetOptionsPrefix(B,prefix));
231: }
232: PetscFunctionReturn(PETSC_SUCCESS);
233: }
235: SLEPC_INTERN PetscErrorCode NEPSetDimensions_Default(NEP,PetscInt,PetscInt*,PetscInt*);
236: SLEPC_INTERN PetscErrorCode NEPComputeVectors(NEP);
237: SLEPC_INTERN PetscErrorCode NEPReset_Problem(NEP);
238: SLEPC_INTERN PetscErrorCode NEPGetDefaultShift(NEP,PetscScalar*);
239: SLEPC_INTERN PetscErrorCode NEPComputeVectors_Schur(NEP);
240: SLEPC_INTERN PetscErrorCode NEPComputeResidualNorm_Private(NEP,PetscBool,PetscScalar,Vec,Vec*,PetscReal*);
241: SLEPC_INTERN PetscErrorCode NEPNewtonRefinementSimple(NEP,PetscInt*,PetscReal,PetscInt);
243: #endif