Actual source code: dvdschm.c
slepc-3.22.2 2024-12-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: PetscFunctionBegin;
20: PetscCall(PetscMemzero(b,sizeof(dvdBlackboard)));
21: b->state = DVD_STATE_PRECONF;
23: for (check_sum0=-1,check_sum1=DVD_CHECKSUM(b); check_sum0 != check_sum1; check_sum0 = check_sum1, check_sum1 = DVD_CHECKSUM(b)) {
25: /* Setup basic management of V */
26: PetscCall(dvd_managementV_basic(d,b,bs,mpd,min_size_V,plusk,PetscNot(harmMode==DVD_HARM_NONE),allResiduals));
28: /* Setup the initial subspace for V */
29: PetscCall(dvd_initV(d,b,ini_size_V,size_initV,(init==DVD_INITV_KRYLOV)?PETSC_TRUE:PETSC_FALSE));
31: /* Setup the convergence in order to use the SLEPc convergence test */
32: PetscCall(dvd_testconv_slepc(d,b));
34: /* Setup Raileigh-Ritz for selecting the best eigenpairs in V */
35: PetscCall(dvd_calcpairs_qz(d,b,orth,PetscNot(harmMode==DVD_HARM_NONE)));
36: if (harmMode != DVD_HARM_NONE) PetscCall(dvd_harm_conf(d,b,harmMode,PETSC_FALSE,0.0));
38: /* Setup the method for improving the eigenvectors */
39: if (doubleexp) PetscCall(dvd_improvex_gd2(d,b,ksp,bs));
40: else {
41: PetscCall(dvd_improvex_jd(d,b,ksp,bs,PETSC_FALSE));
42: PetscCall(dvd_improvex_jd_proj_uv(d,b));
43: PetscCall(dvd_improvex_jd_lit_const(d,b,0,0.0,0.0));
44: }
45: }
46: PetscFunctionReturn(PETSC_SUCCESS);
47: }
49: 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)
50: {
51: PetscInt check_sum0,check_sum1,maxits;
52: PetscReal tol;
54: PetscFunctionBegin;
55: b->state = DVD_STATE_CONF;
56: check_sum0 = DVD_CHECKSUM(b);
58: /* Setup basic management of V */
59: PetscCall(dvd_managementV_basic(d,b,bs,mpd,min_size_V,plusk,PetscNot(harmMode==DVD_HARM_NONE),allResiduals));
61: /* Setup the initial subspace for V */
62: PetscCall(dvd_initV(d,b,ini_size_V,size_initV,(init==DVD_INITV_KRYLOV)?PETSC_TRUE:PETSC_FALSE));
64: /* Setup the convergence in order to use the SLEPc convergence test */
65: PetscCall(dvd_testconv_slepc(d,b));
67: /* Setup Raileigh-Ritz for selecting the best eigenpairs in V */
68: PetscCall(dvd_calcpairs_qz(d,b,orth,PetscNot(harmMode==DVD_HARM_NONE)));
69: if (harmMode != DVD_HARM_NONE) PetscCall(dvd_harm_conf(d,b,harmMode,fixedTarget,t));
71: /* Setup the method for improving the eigenvectors */
72: if (doubleexp) PetscCall(dvd_improvex_gd2(d,b,ksp,bs));
73: else {
74: PetscCall(dvd_improvex_jd(d,b,ksp,bs,dynamic));
75: PetscCall(dvd_improvex_jd_proj_uv(d,b));
76: PetscCall(KSPGetTolerances(ksp,&tol,NULL,NULL,&maxits));
77: PetscCall(dvd_improvex_jd_lit_const(d,b,maxits,tol,fix));
78: }
80: check_sum1 = DVD_CHECKSUM(b);
81: PetscAssert(check_sum0==check_sum1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Something awful happened");
82: PetscFunctionReturn(PETSC_SUCCESS);
83: }