Actual source code: dvdinitv.c
slepc-3.18.0 2022-10-01
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: */
10: /*
11: SLEPc eigensolver: "davidson"
13: Step: initialize subspace V
14: */
16: #include "davidson.h"
18: typedef struct {
19: PetscInt k; /* desired initial subspace size */
20: PetscInt user; /* number of user initial vectors */
21: void *old_initV_data; /* old initV data */
22: } dvdInitV;
24: static PetscErrorCode dvd_initV_classic_0(dvdDashboard *d)
25: {
26: dvdInitV *data = (dvdInitV*)d->initV_data;
27: PetscInt i,user = PetscMin(data->user,d->eps->mpd), l,k;
29: BVGetActiveColumns(d->eps->V,&l,&k);
30: /* User vectors are added at the beginning, so no active column should be in V */
31: PetscAssert(data->user==0 || l==0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Consistency broken");
32: /* Generate a set of random initial vectors and orthonormalize them */
33: for (i=l+user;i<l+data->k && i<d->eps->ncv && i-l<d->eps->mpd;i++) BVSetRandomColumn(d->eps->V,i);
34: d->V_tra_s = 0; d->V_tra_e = 0;
35: d->V_new_s = 0; d->V_new_e = i-l;
37: /* After that the user vectors will be destroyed */
38: data->user = 0;
39: return 0;
40: }
42: static PetscErrorCode dvd_initV_krylov_0(dvdDashboard *d)
43: {
44: dvdInitV *data = (dvdInitV*)d->initV_data;
45: PetscInt i,user = PetscMin(data->user,d->eps->mpd),l,k;
46: Vec av,v1,v2;
48: BVGetActiveColumns(d->eps->V,&l,&k);
49: /* User vectors are added at the beginning, so no active column should be in V */
50: PetscAssert(data->user==0 || l==0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Consistency broken");
52: /* If needed, generate a random vector for starting the arnoldi method */
53: if (user == 0) {
54: BVSetRandomColumn(d->eps->V,l);
55: user = 1;
56: }
58: /* Perform k steps of Arnoldi with the operator K^{-1}*(t[1]*A-t[2]*B) */
59: dvd_orthV(d->eps->V,l,l+user);
60: for (i=l+user;i<l+data->k && i<d->eps->ncv && i-l<d->eps->mpd;i++) {
61: /* aux <- theta[1]A*in - theta[0]*B*in */
62: BVGetColumn(d->eps->V,i,&v1);
63: BVGetColumn(d->eps->V,i-user,&v2);
64: BVGetColumn(d->auxBV,0,&av);
65: if (d->B) {
66: MatMult(d->A,v2,v1);
67: MatMult(d->B,v2,av);
68: VecAXPBY(av,d->target[1],-d->target[0],v1);
69: } else {
70: MatMult(d->A,v2,av);
71: VecAXPBY(av,-d->target[0],d->target[1],v2);
72: }
73: d->improvex_precond(d,0,av,v1);
74: BVRestoreColumn(d->eps->V,i,&v1);
75: BVRestoreColumn(d->eps->V,i-user,&v2);
76: BVRestoreColumn(d->auxBV,0,&av);
77: dvd_orthV(d->eps->V,i,i+1);
78: }
80: d->V_tra_s = 0; d->V_tra_e = 0;
81: d->V_new_s = 0; d->V_new_e = i-l;
83: /* After that the user vectors will be destroyed */
84: data->user = 0;
85: return 0;
86: }
88: static PetscErrorCode dvd_initV_d(dvdDashboard *d)
89: {
90: dvdInitV *data = (dvdInitV*)d->initV_data;
92: /* Restore changes in dvdDashboard */
93: d->initV_data = data->old_initV_data;
95: /* Free local data */
96: PetscFree(data);
97: return 0;
98: }
100: PetscErrorCode dvd_initV(dvdDashboard *d, dvdBlackboard *b, PetscInt k,PetscInt user, PetscBool krylov)
101: {
102: dvdInitV *data;
104: /* Setting configuration constrains */
105: b->max_size_V = PetscMax(b->max_size_V, k);
107: /* Setup the step */
108: if (b->state >= DVD_STATE_CONF) {
109: PetscNew(&data);
110: data->k = k;
111: data->user = user;
112: data->old_initV_data = d->initV_data;
113: d->initV_data = data;
114: if (krylov) d->initV = dvd_initV_krylov_0;
115: else d->initV = dvd_initV_classic_0;
116: EPSDavidsonFLAdd(&d->destroyList,dvd_initV_d);
117: }
118: return 0;
119: }
121: PetscErrorCode dvd_orthV(BV V,PetscInt V_new_s,PetscInt V_new_e)
122: {
123: PetscInt i;
125: for (i=V_new_s;i<V_new_e;i++) BVOrthonormalizeColumn(V,i,PETSC_TRUE,NULL,NULL);
126: return 0;
127: }