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