Actual source code: davidson.h
slepc-3.7.2 2016-07-19
1: /*
2: Method: General Davidson Method (includes GD and JD)
4: References:
5: - Ernest R. Davidson. Super-matrix methods. Computer Physics Communications,
6: 53:49-60, May 1989.
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: SLEPc - Scalable Library for Eigenvalue Problem Computations
10: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
12: This file is part of SLEPc.
14: SLEPc is free software: you can redistribute it and/or modify it under the
15: terms of version 3 of the GNU Lesser General Public License as published by
16: the Free Software Foundation.
18: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
19: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
21: more details.
23: You should have received a copy of the GNU Lesser General Public License
24: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
25: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26: */
28: #include <slepc/private/epsimpl.h>
29: #include <slepc/private/vecimplslepc.h>
31: struct _dvdDashboard;
32: typedef PetscErrorCode (*dvdCallback)(struct _dvdDashboard*);
33: typedef struct _dvdFunctionList {
34: dvdCallback f;
35: struct _dvdFunctionList *next;
36: } dvdFunctionList;
38: typedef enum {
39: DVD_HARM_NONE,
40: DVD_HARM_RR,
41: DVD_HARM_RRR,
42: DVD_HARM_REIGS,
43: DVD_HARM_LEIGS
44: } HarmType_t;
46: typedef enum {
47: DVD_INITV_CLASSIC,
48: DVD_INITV_KRYLOV
49: } InitType_t;
51: typedef enum {
52: DVD_PROJ_KXX,
53: DVD_PROJ_KZX
54: } ProjType_t;
56: /*
57: Dashboard struct: contains the methods that will be employed and the tuning
58: options.
59: */
61: typedef struct _dvdDashboard {
62: /**** Function steps ****/
63: /* Initialize V */
64: PetscErrorCode (*initV)(struct _dvdDashboard*);
65: void *initV_data;
67: /* Find the approximate eigenpairs from V */
68: PetscErrorCode (*calcPairs)(struct _dvdDashboard*);
69: void *calcPairs_data;
71: /* Eigenpair test for convergence */
72: PetscBool (*testConv)(struct _dvdDashboard*,PetscScalar,PetscScalar,PetscReal,PetscReal*);
73: void *testConv_data;
75: /* Improve the selected eigenpairs */
76: PetscErrorCode (*improveX)(struct _dvdDashboard*,PetscInt,PetscInt,PetscInt*);
77: void *improveX_data;
79: /* Check for restarting */
80: PetscErrorCode (*isRestarting)(struct _dvdDashboard*,PetscBool*);
81: void *isRestarting_data;
83: /* Perform restarting */
84: PetscErrorCode (*restartV)(struct _dvdDashboard*);
85: void *restartV_data;
87: /* Update V */
88: PetscErrorCode (*updateV)(struct _dvdDashboard*);
89: void *updateV_data;
91: /**** Problem specification ****/
92: Mat A,B; /* problem matrices */
93: MatType_t sA,sB; /* matrix specifications */
94: EPType_t sEP; /* problem specifications */
95: PetscInt nev; /* number of eigenpairs */
96: EPSWhich which; /* spectrum selection */
97: PetscBool withTarget; /* if there is a target */
98: PetscScalar target[2]; /* target value */
99: PetscReal tol; /* tolerance */
100: PetscBool correctXnorm; /* if true, norm of X are computed */
102: /**** Subspaces specification ****/
103: PetscInt nconv; /* number of converged eigenpairs */
104: PetscInt npreconv; /* number of pairs ready to converge */
106: BV W; /* left basis for harmonic case */
107: BV AX; /* A*V */
108: BV BX; /* B*V */
109: PetscInt size_D; /* active vectors */
110: PetscInt max_size_proj; /* max size projected problem */
111: PetscInt max_cX_in_proj; /* max vectors from cX in the projected problem */
112: PetscInt max_cX_in_impr; /* max vectros from cX in the projector */
113: PetscInt max_size_P; /* max unconverged vectors in the projector */
114: PetscInt bs; /* max vectors that expands the subspace every iteration */
115: EPS eps; /* connection to SLEPc */
117: /**** Auxiliary space ****/
118: VecPool auxV; /* auxiliary vectors */
119: BV auxBV; /* auxiliary vectors */
121: /**** Eigenvalues and errors ****/
122: PetscScalar *ceigr,*ceigi; /* converged eigenvalues */
123: PetscScalar *eigr,*eigi; /* current eigenvalues */
124: PetscReal *nR; /* residual norm */
125: PetscReal *real_nR; /* original nR */
126: PetscReal *nX; /* X norm */
127: PetscReal *real_nX; /* original nX */
128: PetscReal *errest; /* relative error eigenpairs */
129: PetscReal *nBds; /* B-norms of projected problem */
131: /**** Shared function and variables ****/
132: PetscErrorCode (*e_Vchanged)(struct _dvdDashboard*,PetscInt,PetscInt,PetscInt,PetscInt);
133: void *e_Vchanged_data;
134: PetscErrorCode (*calcpairs_residual)(struct _dvdDashboard*,PetscInt,PetscInt);
135: PetscErrorCode (*calcpairs_selectPairs)(struct _dvdDashboard*,PetscInt);
136: void *calcpairs_residual_data;
137: PetscErrorCode (*improvex_precond)(struct _dvdDashboard*,PetscInt,Vec,Vec);
138: void *improvex_precond_data;
139: PetscErrorCode (*improvex_jd_proj_uv)(struct _dvdDashboard*,PetscInt,PetscInt,Vec*,Vec*,Vec*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt);
140: PetscErrorCode (*improvex_jd_lit)(struct _dvdDashboard*,PetscInt,PetscScalar*,PetscScalar*,PetscInt*,PetscReal*);
141: PetscErrorCode (*calcpairs_W)(struct _dvdDashboard*);
142: void *calcpairs_W_data;
143: PetscErrorCode (*calcpairs_proj_trans)(struct _dvdDashboard*);
144: PetscErrorCode (*calcpairs_eigs_trans)(struct _dvdDashboard*);
145: PetscErrorCode (*calcpairs_eig_backtrans)(struct _dvdDashboard*,PetscScalar,PetscScalar,PetscScalar*,PetscScalar*);
146: PetscErrorCode (*calcpairs_proj_res)(struct _dvdDashboard*,PetscInt,PetscInt,Vec*);
147: PetscErrorCode (*preTestConv)(struct _dvdDashboard*,PetscInt,PetscInt,PetscInt,PetscInt*);
148: PetscErrorCode (*e_newIteration)(struct _dvdDashboard*);
149: void *e_newIteration_data;
151: dvdFunctionList *startList; /* starting list */
152: dvdFunctionList *endList; /* ending list */
153: dvdFunctionList *destroyList;/* destructor list */
155: Mat H,G; /* projected problem matrices */
156: Mat auxM; /* auxiliary dense matrix */
157: PetscInt size_MT; /* rows in MT */
159: PetscInt V_tra_s;
160: PetscInt V_tra_e; /* cX <- [cX V*MT(0:V_tra_s-1)], V <- V*MT(V_tra_s:V_tra_e) */
161: PetscInt V_new_s;
162: PetscInt V_new_e; /* added to V the columns V_new_s:V_new_e */
163: PetscBool BV_shift; /* if true BV is shifted when vectors converge */
164: PetscBool W_shift; /* if true W is shifted when vectors converge */
165: } dvdDashboard;
167: typedef struct {
168: /*------------------------- User parameters ---------------------------*/
169: PetscInt blocksize; /* block size */
170: PetscInt initialsize; /* initial size of V */
171: PetscInt minv; /* size of V after restarting */
172: PetscInt plusk; /* keep plusk eigenvectors from the last iteration */
173: PetscBool ipB; /* true if B-ortho is used */
174: PetscReal fix; /* the fix parameter */
175: PetscBool krylovstart; /* true if the starting subspace is a Krylov basis */
176: PetscBool dynamic; /* true if dynamic stopping criterion is used */
177: PetscInt cX_in_proj; /* converged vectors in the projected problem */
178: PetscInt cX_in_impr; /* converged vectors in the projector */
179: PetscBool doubleexp; /* double expansion in GD (GD2) */
181: /*----------------- Child objects and working data -------------------*/
182: dvdDashboard ddb;
183: } EPS_DAVIDSON;
187: PETSC_STATIC_INLINE PetscErrorCode EPSDavidsonFLAdd(dvdFunctionList **fl,dvdCallback f)
188: {
190: dvdFunctionList *l;
193: PetscNew(&l);
194: l->f = f;
195: l->next = *fl;
196: *fl = l;
197: return(0);
198: }
199:
202: PETSC_STATIC_INLINE PetscErrorCode EPSDavidsonFLCall(dvdFunctionList *fl,dvdDashboard *d)
203: {
205: dvdFunctionList *l;
208: for (l=fl;l;l=l->next) { (l->f)(d); }
209: return(0);
210: }
214: PETSC_STATIC_INLINE PetscErrorCode EPSDavidsonFLDestroy(dvdFunctionList **fl)
215: {
216: PetscErrorCode ierr;
217: dvdFunctionList *l,*l0;
220: for (l=*fl;l;l=l0) {
221: l0 = l->next;
222: PetscFree(l);
223: }
224: *fl = NULL;
225: return(0);
226: }
228: /*
229: The blackboard configuration structure: saves information about the memory
230: and other requirements.
232: The starting memory structure:
234: V W? AV BV? tKZ
235: |-----------|-----------|-----------|------------|-------|
236: nev+mpd nev+mpd scP+mpd nev?+mpd sP+scP
237: scP+mpd scP+mpd
239: The final memory structure considering W_shift and BV_shift:
241: cX V cY? W? cAV AV BcX? BV? KZ tKZ
242: |---|-------|----|------|---|-------|----|-------|---|---|
243: nev mpd nev mpd scP mpd nev mpd scP sP <- shift
244: scP scP <- !shift
245: */
246: typedef struct {
247: PetscInt max_size_V; /* max size of the searching subspace (mpd) */
248: PetscInt max_size_X; /* max size of X (bs) */
249: PetscInt size_V; /* real size of V (nev+size_P+mpd) */
250: PetscInt max_size_oldX; /* max size of oldX */
251: PetscInt max_nev; /* max number of converged pairs */
252: PetscInt max_size_P; /* number of computed vectors for the projector */
253: PetscInt max_size_cP; /* number of converged vectors in the projectors */
254: PetscInt max_size_proj; /* max size projected problem */
255: PetscInt max_size_cX_proj; /* max converged vectors in the projected problem */
256: PetscInt state; /* method states:
257: 0: preconfiguring
258: 1: configuring
259: 2: running */
260: } dvdBlackboard;
262: #define DVD_STATE_PRECONF 0
263: #define DVD_STATE_CONF 1
264: #define DVD_STATE_RUN 2
266: /* Prototypes of non-static auxiliary functions */
267: PETSC_INTERN PetscErrorCode dvd_calcpairs_qz(dvdDashboard*,dvdBlackboard*,PetscBool,PetscInt,PetscBool);
268: PETSC_INTERN PetscErrorCode dvd_improvex_gd2(dvdDashboard*,dvdBlackboard*,KSP,PetscInt);
269: PETSC_INTERN PetscErrorCode dvd_improvex_jd(dvdDashboard*,dvdBlackboard*,KSP,PetscInt,PetscInt,PetscBool);
270: PETSC_INTERN PetscErrorCode dvd_improvex_jd_proj_uv(dvdDashboard*,dvdBlackboard*,ProjType_t);
271: PETSC_INTERN PetscErrorCode dvd_improvex_jd_lit_const(dvdDashboard*,dvdBlackboard*,PetscInt,PetscReal,PetscReal);
272: PETSC_INTERN PetscErrorCode dvd_improvex_compute_X(dvdDashboard*,PetscInt,PetscInt,Vec*,PetscScalar*,PetscInt);
273: PETSC_INTERN PetscErrorCode dvd_initV(dvdDashboard*,dvdBlackboard*,PetscInt,PetscInt,PetscBool);
274: PETSC_INTERN PetscErrorCode dvd_orthV(BV,PetscInt,PetscInt);
275: PETSC_INTERN PetscErrorCode dvd_schm_basic_preconf(dvdDashboard*,dvdBlackboard*,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,HarmType_t,KSP,InitType_t,PetscBool,PetscBool,PetscInt,PetscInt,PetscBool);
276: PETSC_INTERN PetscErrorCode dvd_schm_basic_conf(dvdDashboard*,dvdBlackboard*,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,HarmType_t,PetscBool,PetscScalar,KSP,PetscReal,InitType_t,PetscBool,PetscBool,PetscInt,PetscInt,PetscBool,PetscBool);
277: PETSC_INTERN PetscErrorCode dvd_testconv_basic(dvdDashboard*,dvdBlackboard*);
278: PETSC_INTERN PetscErrorCode dvd_testconv_slepc(dvdDashboard*,dvdBlackboard*);
279: PETSC_INTERN PetscErrorCode dvd_managementV_basic(dvdDashboard*,dvdBlackboard*,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool,PetscBool);
280: PETSC_INTERN PetscErrorCode dvd_static_precond_PC(dvdDashboard*,dvdBlackboard*,PC);
281: PETSC_INTERN PetscErrorCode dvd_jacobi_precond(dvdDashboard*,dvdBlackboard*);
282: PETSC_INTERN PetscErrorCode dvd_harm_updateproj(dvdDashboard*);
283: PETSC_INTERN PetscErrorCode dvd_harm_conf(dvdDashboard*,dvdBlackboard*,HarmType_t,PetscBool,PetscScalar);
285: /* Internal interface routines */
286: PETSC_INTERN PetscErrorCode EPSReset_XD(EPS);
287: PETSC_INTERN PetscErrorCode EPSSetUp_XD(EPS);
288: PETSC_INTERN PetscErrorCode EPSSolve_XD(EPS);
289: PETSC_INTERN PetscErrorCode EPSComputeVectors_XD(EPS);
290: PETSC_INTERN PetscErrorCode EPSXDSetKrylovStart_XD(EPS,PetscBool);
291: PETSC_INTERN PetscErrorCode EPSXDGetKrylovStart_XD(EPS,PetscBool*);
292: PETSC_INTERN PetscErrorCode EPSXDSetBlockSize_XD(EPS,PetscInt);
293: PETSC_INTERN PetscErrorCode EPSXDGetBlockSize_XD(EPS,PetscInt*);
294: PETSC_INTERN PetscErrorCode EPSXDSetRestart_XD(EPS,PetscInt,PetscInt);
295: PETSC_INTERN PetscErrorCode EPSXDGetRestart_XD(EPS,PetscInt*,PetscInt*);
296: PETSC_INTERN PetscErrorCode EPSXDGetInitialSize_XD(EPS,PetscInt*);
297: PETSC_INTERN PetscErrorCode EPSXDSetInitialSize_XD(EPS,PetscInt);
298: PETSC_INTERN PetscErrorCode EPSXDGetFix_XD(EPS,PetscReal*);
299: PETSC_INTERN PetscErrorCode EPSJDSetFix_JD(EPS,PetscReal);
300: PETSC_INTERN PetscErrorCode EPSXDSetBOrth_XD(EPS,PetscBool);
301: PETSC_INTERN PetscErrorCode EPSXDGetBOrth_XD(EPS,PetscBool*);
302: PETSC_INTERN PetscErrorCode EPSJDSetConstCorrectionTol_JD(EPS,PetscBool);
303: PETSC_INTERN PetscErrorCode EPSJDGetConstCorrectionTol_JD(EPS,PetscBool*);
304: PETSC_INTERN PetscErrorCode EPSXDSetWindowSizes_XD(EPS,PetscInt,PetscInt);
305: PETSC_INTERN PetscErrorCode EPSXDGetWindowSizes_XD(EPS,PetscInt*,PetscInt*);