Actual source code: dvd_schm.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
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: #include davidson.h
24: #define DVD_CHECKSUM(b) \
25: ( (b)->max_size_V + (b)->max_size_auxV + (b)->max_size_auxS + \
26: (b)->own_vecs + (b)->own_scalars + (b)->max_size_oldX )
30: PetscErrorCode dvd_schm_basic_preconf(dvdDashboard *d, dvdBlackboard *b,
31: PetscInt mpd, PetscInt min_size_V, PetscInt bs,
32: PetscInt ini_size_V, PetscInt size_initV, PetscInt plusk,
33: HarmType_t harmMode, KSP ksp, InitType_t init, PetscBool allResiduals,
34: orthoV_type_t orth, PetscInt cX_proj, PetscInt cX_impr)
35: {
37: PetscInt check_sum0, check_sum1;
41: PetscMemzero(b, sizeof(dvdBlackboard));
42: b->state = DVD_STATE_PRECONF;
44: for(check_sum0=-1,check_sum1=DVD_CHECKSUM(b); check_sum0 != check_sum1;
45: check_sum0 = check_sum1, check_sum1 = DVD_CHECKSUM(b)) {
46: b->own_vecs = b->own_scalars = 0;
48: /* Setup basic management of V */
49: dvd_managementV_basic(d, b, bs, mpd, min_size_V, plusk,
50: harmMode==DVD_HARM_NONE?PETSC_FALSE:PETSC_TRUE,
51: allResiduals);
52:
53:
54: /* Setup the initial subspace for V */
55: dvd_initV(d, b, ini_size_V, size_initV,
56: init==DVD_INITV_KRYLOV?PETSC_TRUE:PETSC_FALSE);
57:
58: /* Setup the convergence in order to use the SLEPc convergence test */
59: dvd_testconv_slepc(d, b);
60:
61: /* Setup Raileigh-Ritz for selecting the best eigenpairs in V */
62: dvd_calcpairs_qz(d, b, orth, PETSC_NULL, cX_proj,
63: harmMode==DVD_HARM_NONE?PETSC_FALSE:PETSC_TRUE);
64: if (harmMode != DVD_HARM_NONE) {
65: dvd_harm_conf(d, b, harmMode, PETSC_FALSE, 0.0);
66: }
67:
68: /* Setup the method for improving the eigenvectors */
69: dvd_improvex_jd(d, b, ksp, bs, cX_impr, PETSC_FALSE);
70: dvd_improvex_jd_proj_uv(d, b, DVD_PROJ_KZX);
71: dvd_improvex_jd_lit_const(d, b, 0, 0.0, 0.0);
72: }
74: return(0);
75: }
80: PetscErrorCode dvd_schm_basic_conf(dvdDashboard *d, dvdBlackboard *b,
81: PetscInt mpd, PetscInt min_size_V, PetscInt bs,
82: PetscInt ini_size_V, PetscInt size_initV, PetscInt plusk,
83: IP ip, HarmType_t harmMode, PetscBool fixedTarget, PetscScalar t, KSP ksp,
84: PetscReal fix, InitType_t init, PetscBool allResiduals, orthoV_type_t orth,
85: PetscInt cX_proj, PetscInt cX_impr, PetscBool dynamic)
86: {
87: PetscInt check_sum0, check_sum1, maxits;
88: Vec *fv;
89: PetscScalar *fs;
90: PetscReal tol;
91: PetscErrorCode ierr;
95: b->state = DVD_STATE_CONF;
96: check_sum0 = DVD_CHECKSUM(b);
97: b->own_vecs = 0; b->own_scalars = 0;
98: fv = b->free_vecs; fs = b->free_scalars;
100: /* Setup basic management of V */
101: dvd_managementV_basic(d, b, bs, mpd, min_size_V, plusk,
102: harmMode==DVD_HARM_NONE?PETSC_FALSE:PETSC_TRUE,
103: allResiduals);
104:
106: /* Setup the initial subspace for V */
107: dvd_initV(d, b, ini_size_V, size_initV,
108: init==DVD_INITV_KRYLOV?PETSC_TRUE:PETSC_FALSE);
110: /* Setup the convergence in order to use the SLEPc convergence test */
111: dvd_testconv_slepc(d, b);
113: /* Setup Raileigh-Ritz for selecting the best eigenpairs in V */
114: dvd_calcpairs_qz(d, b, orth, ip, cX_proj,
115: harmMode==DVD_HARM_NONE?PETSC_FALSE:PETSC_TRUE);
116: if (harmMode != DVD_HARM_NONE) {
117: dvd_harm_conf(d, b, harmMode, fixedTarget, t);
118: }
120: /* Setup the method for improving the eigenvectors */
121: dvd_improvex_jd(d, b, ksp, bs, cX_impr, dynamic);
122: dvd_improvex_jd_proj_uv(d, b, DVD_PROJ_KZX);
123:
124: KSPGetTolerances(ksp, &tol, PETSC_NULL, PETSC_NULL, &maxits);
125:
126: dvd_improvex_jd_lit_const(d, b, maxits, tol, fix);
128: check_sum1 = DVD_CHECKSUM(b);
129: if ((check_sum0 != check_sum1) ||
130: (b->free_vecs - fv > b->own_vecs) ||
131: (b->free_scalars - fs > b->own_scalars))
132: SETERRQ(PETSC_COMM_WORLD,1, "Something awful happened");
133:
134: return(0);
135: }