Actual source code: epsimpl.h
slepc-3.6.1 2015-09-03
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: #if !defined(_EPSIMPL)
23: #define _EPSIMPL
25: #include <slepceps.h>
26: #include <slepc/private/slepcimpl.h>
28: PETSC_EXTERN PetscBool EPSRegisterAllCalled;
29: PETSC_EXTERN PetscErrorCode EPSRegisterAll(void);
30: PETSC_EXTERN PetscLogEvent EPS_SetUp,EPS_Solve;
32: typedef struct _EPSOps *EPSOps;
34: struct _EPSOps {
35: PetscErrorCode (*solve)(EPS);
36: PetscErrorCode (*setup)(EPS);
37: PetscErrorCode (*setfromoptions)(PetscOptions*,EPS);
38: PetscErrorCode (*publishoptions)(EPS);
39: PetscErrorCode (*destroy)(EPS);
40: PetscErrorCode (*reset)(EPS);
41: PetscErrorCode (*view)(EPS,PetscViewer);
42: PetscErrorCode (*backtransform)(EPS);
43: PetscErrorCode (*computevectors)(EPS);
44: };
46: /*
47: Maximum number of monitors you can run with a single EPS
48: */
49: #define MAXEPSMONITORS 5
51: typedef enum { EPS_STATE_INITIAL,
52: EPS_STATE_SETUP,
53: EPS_STATE_SOLVED,
54: EPS_STATE_EIGENVECTORS } EPSStateType;
56: /*
57: Defines the EPS data structure.
58: */
59: struct _p_EPS {
60: PETSCHEADER(struct _EPSOps);
61: /*------------------------- User parameters ---------------------------*/
62: PetscInt max_it; /* maximum number of iterations */
63: PetscInt nev; /* number of eigenvalues to compute */
64: PetscInt ncv; /* number of basis vectors */
65: PetscInt mpd; /* maximum dimension of projected problem */
66: PetscInt nini; /* number of initial vectors (negative means not copied yet) */
67: PetscInt nds; /* number of basis vectors of deflation space */
68: PetscScalar target; /* target value */
69: PetscReal tol; /* tolerance */
70: EPSConv conv; /* convergence test */
71: EPSWhich which; /* which part of the spectrum to be sought */
72: PetscReal inta,intb; /* interval [a,b] for spectrum slicing */
73: EPSProblemType problem_type; /* which kind of problem to be solved */
74: EPSExtraction extraction; /* which kind of extraction to be applied */
75: EPSBalance balance; /* the balancing method */
76: PetscInt balance_its; /* number of iterations of the balancing method */
77: PetscReal balance_cutoff; /* cutoff value for balancing */
78: PetscBool trueres; /* whether the true residual norm must be computed */
79: PetscBool trackall; /* whether all the residuals must be computed */
80: PetscBool purify; /* whether eigenvectors need to be purified */
82: /*-------------- User-provided functions and contexts -----------------*/
83: PetscErrorCode (*converged)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
84: PetscErrorCode (*convergeddestroy)(void*);
85: PetscErrorCode (*arbitrary)(PetscScalar,PetscScalar,Vec,Vec,PetscScalar*,PetscScalar*,void*);
86: void *convergedctx;
87: void *arbitraryctx;
88: PetscErrorCode (*monitor[MAXEPSMONITORS])(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
89: PetscErrorCode (*monitordestroy[MAXEPSMONITORS])(void**);
90: void *monitorcontext[MAXEPSMONITORS];
91: PetscInt numbermonitors;
93: /*----------------- Child objects and working data -------------------*/
94: ST st; /* spectral transformation object */
95: DS ds; /* direct solver object */
96: BV V; /* set of basis vectors and computed eigenvectors */
97: RG rg; /* optional region for filtering */
98: PetscRandom rand; /* random number generator */
99: SlepcSC sc; /* sorting criterion data */
100: Vec D; /* diagonal matrix for balancing */
101: Vec *IS; /* references to user-provided initial space */
102: Vec *defl; /* references to user-provided deflation space */
103: PetscScalar *eigr,*eigi; /* real and imaginary parts of eigenvalues */
104: PetscReal *errest; /* error estimates */
105: PetscScalar *rr,*ri; /* values computed by user's arbitrary selection function */
106: PetscInt *perm; /* permutation for eigenvalue ordering */
107: PetscInt nwork; /* number of work vectors */
108: Vec *work; /* work vectors */
109: void *data; /* placeholder for solver-specific stuff */
111: /* ----------------------- Status variables --------------------------*/
112: EPSStateType state; /* initial -> setup -> solved -> eigenvectors */
113: PetscInt nconv; /* number of converged eigenvalues */
114: PetscInt its; /* number of iterations so far computed */
115: PetscInt n,nloc; /* problem dimensions (global, local) */
116: PetscReal nrma,nrmb; /* computed matrix norms */
117: PetscBool isgeneralized;
118: PetscBool ispositive;
119: PetscBool ishermitian;
120: EPSConvergedReason reason;
121: };
123: /*
124: Macros to test valid EPS arguments
125: */
126: #if !defined(PETSC_USE_DEBUG)
128: #define EPSCheckSolved(h,arg) do {} while (0)
130: #else
132: #define EPSCheckSolved(h,arg) \
133: do { \
134: if (h->state<EPS_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSolve() first: Parameter #%d",arg); \
135: } while (0)
137: #endif
141: /*
142: EPS_SetInnerProduct - set B matrix for inner product if appropriate.
143: */
144: PETSC_STATIC_INLINE PetscErrorCode EPS_SetInnerProduct(EPS eps)
145: {
147: Mat B;
150: if (!eps->V) { EPSGetBV(eps,&eps->V); }
151: if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
152: STGetBilinearForm(eps->st,&B);
153: BVSetMatrix(eps->V,B,eps->ispositive?PETSC_FALSE:PETSC_TRUE);
154: MatDestroy(&B);
155: } else {
156: BVSetMatrix(eps->V,NULL,PETSC_FALSE);
157: }
158: return(0);
159: }
161: PETSC_INTERN PetscErrorCode EPSSetWhichEigenpairs_Default(EPS);
162: PETSC_INTERN PetscErrorCode EPSSetDimensions_Default(EPS,PetscInt,PetscInt*,PetscInt*);
163: PETSC_INTERN PetscErrorCode EPSBackTransform_Default(EPS);
164: PETSC_INTERN PetscErrorCode EPSComputeVectors(EPS);
165: PETSC_INTERN PetscErrorCode EPSComputeVectors_Hermitian(EPS);
166: PETSC_INTERN PetscErrorCode EPSComputeVectors_Schur(EPS);
167: PETSC_INTERN PetscErrorCode EPSComputeVectors_Indefinite(EPS);
168: PETSC_INTERN PetscErrorCode EPSComputeVectors_Slice(EPS);
169: PETSC_INTERN PetscErrorCode EPSComputeResidualNorm_Private(EPS,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*);
170: PETSC_INTERN PetscErrorCode EPSComputeRitzVector(EPS,PetscScalar*,PetscScalar*,BV,Vec,Vec);
171: PETSC_INTERN PetscErrorCode EPSGetStartVector(EPS,PetscInt,PetscBool*);
173: /* Private functions of the solver implementations */
175: PETSC_INTERN PetscErrorCode EPSBasicArnoldi(EPS,PetscBool,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
176: PETSC_INTERN PetscErrorCode EPSDelayedArnoldi(EPS,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
177: PETSC_INTERN PetscErrorCode EPSDelayedArnoldi1(EPS,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
178: PETSC_INTERN PetscErrorCode EPSKrylovConvergence(EPS,PetscBool,PetscInt,PetscInt,PetscReal,PetscReal,PetscInt*);
179: PETSC_INTERN PetscErrorCode EPSFullLanczos(EPS,PetscReal*,PetscReal*,PetscInt,PetscInt*,PetscBool*);
180: PETSC_INTERN PetscErrorCode EPSPseudoLanczos(EPS,PetscReal*,PetscReal*,PetscReal*,PetscInt,PetscInt*,PetscBool*,PetscBool*,PetscReal*,Vec);
181: PETSC_INTERN PetscErrorCode EPSBuildBalance_Krylov(EPS);
183: #endif