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