Actual source code: dvd_schm.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2013, 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:  #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,PetscInt mpd,PetscInt min_size_V,PetscInt bs,PetscInt ini_size_V,PetscInt size_initV,PetscInt plusk,HarmType_t harmMode,KSP ksp,InitType_t init,PetscBool allResiduals,EPSOrthType orth,PetscInt cX_proj,PetscInt cX_impr,Method_t method)
 31: {
 33:   PetscInt       check_sum0, check_sum1;

 36:   PetscMemzero(b, sizeof(dvdBlackboard));
 37:   b->state = DVD_STATE_PRECONF;

 39:   for (check_sum0=-1,check_sum1=DVD_CHECKSUM(b); check_sum0 != check_sum1;
 40:        check_sum0 = check_sum1, check_sum1 = DVD_CHECKSUM(b)) {
 41:     b->own_vecs = b->own_scalars = 0;

 43:     /* Setup basic management of V */
 44:     dvd_managementV_basic(d, b, bs, mpd, min_size_V, plusk,
 45:                                harmMode==DVD_HARM_NONE?PETSC_FALSE:PETSC_TRUE,
 46:                                allResiduals);

 48:     /* Setup the initial subspace for V */
 49:     dvd_initV(d, b, ini_size_V, size_initV,
 50:                      init==DVD_INITV_KRYLOV?PETSC_TRUE:PETSC_FALSE);

 52:     /* Setup the convergence in order to use the SLEPc convergence test */
 53:     dvd_testconv_slepc(d, b);

 55:     /* Setup Raileigh-Ritz for selecting the best eigenpairs in V */
 56:     dvd_calcpairs_qz(d, b, orth, NULL, cX_proj,
 57:                 harmMode==DVD_HARM_NONE?PETSC_FALSE:PETSC_TRUE);
 58:     if (harmMode != DVD_HARM_NONE) {
 59:       dvd_harm_conf(d, b, harmMode, PETSC_FALSE, 0.0);
 60:     }

 62:     /* Setup the method for improving the eigenvectors */
 63:     switch (method) {
 64:       case DVD_METH_GD:
 65:       case DVD_METH_JD:
 66:       dvd_improvex_jd(d, b, ksp, bs, cX_impr, PETSC_FALSE);
 67:       dvd_improvex_jd_proj_uv(d, b, DVD_PROJ_KZX);
 68:       dvd_improvex_jd_lit_const(d, b, 0, 0.0, 0.0);
 69:       break;
 70:       case DVD_METH_GD2:
 71:       dvd_improvex_gd2(d, b, ksp, bs);
 72:       break;
 73:     }
 74:   }
 75:   return(0);
 76: }

 80: PetscErrorCode dvd_schm_basic_conf(dvdDashboard *d,dvdBlackboard *b,PetscInt mpd,PetscInt min_size_V,PetscInt bs,PetscInt ini_size_V,PetscInt size_initV,PetscInt plusk,IP ip,HarmType_t harmMode,PetscBool fixedTarget,PetscScalar t,KSP ksp,PetscReal fix,InitType_t init,PetscBool allResiduals,EPSOrthType orth,PetscInt cX_proj,PetscInt cX_impr,PetscBool dynamic,Method_t method)
 81: {
 82:   PetscInt        check_sum0, check_sum1, maxits;
 83:   Vec             *fv;
 84:   PetscScalar     *fs;
 85:   PetscReal       tol;
 86:   PetscErrorCode  ierr;

 89:   b->state = DVD_STATE_CONF;
 90:   check_sum0 = DVD_CHECKSUM(b);
 91:   b->own_vecs = 0; b->own_scalars = 0;
 92:   fv = b->free_vecs; fs = b->free_scalars;

 94:   /* Setup basic management of V */
 95:   dvd_managementV_basic(d, b, bs, mpd, min_size_V, plusk,
 96:                         harmMode==DVD_HARM_NONE?PETSC_FALSE:PETSC_TRUE,
 97:                         allResiduals);

 99:   /* Setup the initial subspace for V */
100:   dvd_initV(d, b, ini_size_V, size_initV,
101:                    init==DVD_INITV_KRYLOV?PETSC_TRUE:PETSC_FALSE);

103:   /* Setup the convergence in order to use the SLEPc convergence test */
104:   dvd_testconv_slepc(d, b);

106:   /* Setup Raileigh-Ritz for selecting the best eigenpairs in V */
107:   dvd_calcpairs_qz(d, b, orth, ip, cX_proj,
108:                 harmMode==DVD_HARM_NONE?PETSC_FALSE:PETSC_TRUE);
109:   if (harmMode != DVD_HARM_NONE) {
110:     dvd_harm_conf(d, b, harmMode, fixedTarget, t);
111:   }

113:   /* Setup the method for improving the eigenvectors */
114:   switch (method) {
115:     case DVD_METH_GD:
116:     case DVD_METH_JD:
117:       dvd_improvex_jd(d, b, ksp, bs, cX_impr, dynamic);
118:       dvd_improvex_jd_proj_uv(d, b, DVD_PROJ_KZX);
119:       KSPGetTolerances(ksp, &tol, NULL, NULL, &maxits);
120:       dvd_improvex_jd_lit_const(d, b, maxits, tol, fix);
121:       break;
122:     case DVD_METH_GD2:
123:       dvd_improvex_gd2(d, b, ksp, bs);
124:       break;
125:   }

127:   check_sum1 = DVD_CHECKSUM(b);
128:   if ((check_sum0 != check_sum1) ||
129:       (b->free_vecs - fv > b->own_vecs) ||
130:       (b->free_scalars - fs > b->own_scalars))
131:     SETERRQ(PETSC_COMM_SELF,1, "Something awful happened");
132:   return(0);
133: }