Actual source code: dvdschm.c
slepc-3.18.1 2022-11-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: #include "davidson.h"
13: #define DVD_CHECKSUM(b) ((b)->max_size_V + (b)->max_size_oldX)
15: 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,PetscBool orth,PetscBool doubleexp)
16: {
17: PetscInt check_sum0,check_sum1;
19: PetscMemzero(b,sizeof(dvdBlackboard));
20: b->state = DVD_STATE_PRECONF;
22: for (check_sum0=-1,check_sum1=DVD_CHECKSUM(b); check_sum0 != check_sum1; check_sum0 = check_sum1, check_sum1 = DVD_CHECKSUM(b)) {
24: /* Setup basic management of V */
25: dvd_managementV_basic(d,b,bs,mpd,min_size_V,plusk,PetscNot(harmMode==DVD_HARM_NONE),allResiduals);
27: /* Setup the initial subspace for V */
28: dvd_initV(d,b,ini_size_V,size_initV,(init==DVD_INITV_KRYLOV)?PETSC_TRUE:PETSC_FALSE);
30: /* Setup the convergence in order to use the SLEPc convergence test */
31: dvd_testconv_slepc(d,b);
33: /* Setup Raileigh-Ritz for selecting the best eigenpairs in V */
34: dvd_calcpairs_qz(d,b,orth,PetscNot(harmMode==DVD_HARM_NONE));
35: if (harmMode != DVD_HARM_NONE) dvd_harm_conf(d,b,harmMode,PETSC_FALSE,0.0);
37: /* Setup the method for improving the eigenvectors */
38: if (doubleexp) dvd_improvex_gd2(d,b,ksp,bs);
39: else {
40: dvd_improvex_jd(d,b,ksp,bs,PETSC_FALSE);
41: dvd_improvex_jd_proj_uv(d,b);
42: dvd_improvex_jd_lit_const(d,b,0,0.0,0.0);
43: }
44: }
45: return 0;
46: }
48: 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,HarmType_t harmMode,PetscBool fixedTarget,PetscScalar t,KSP ksp,PetscReal fix,InitType_t init,PetscBool allResiduals,PetscBool orth,PetscBool dynamic,PetscBool doubleexp)
49: {
50: PetscInt check_sum0,check_sum1,maxits;
51: PetscReal tol;
53: b->state = DVD_STATE_CONF;
54: check_sum0 = DVD_CHECKSUM(b);
56: /* Setup basic management of V */
57: dvd_managementV_basic(d,b,bs,mpd,min_size_V,plusk,PetscNot(harmMode==DVD_HARM_NONE),allResiduals);
59: /* Setup the initial subspace for V */
60: dvd_initV(d,b,ini_size_V,size_initV,(init==DVD_INITV_KRYLOV)?PETSC_TRUE:PETSC_FALSE);
62: /* Setup the convergence in order to use the SLEPc convergence test */
63: dvd_testconv_slepc(d,b);
65: /* Setup Raileigh-Ritz for selecting the best eigenpairs in V */
66: dvd_calcpairs_qz(d,b,orth,PetscNot(harmMode==DVD_HARM_NONE));
67: if (harmMode != DVD_HARM_NONE) dvd_harm_conf(d,b,harmMode,fixedTarget,t);
69: /* Setup the method for improving the eigenvectors */
70: if (doubleexp) dvd_improvex_gd2(d,b,ksp,bs);
71: else {
72: dvd_improvex_jd(d,b,ksp,bs,dynamic);
73: dvd_improvex_jd_proj_uv(d,b);
74: KSPGetTolerances(ksp,&tol,NULL,NULL,&maxits);
75: dvd_improvex_jd_lit_const(d,b,maxits,tol,fix);
76: }
78: check_sum1 = DVD_CHECKSUM(b);
79: PetscAssert(check_sum0==check_sum1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Something awful happened");
80: return 0;
81: }