Actual source code: dvdinitv.c
slepc-3.7.4 2017-05-17
1: /*
2: SLEPc eigensolver: "davidson"
4: Step: init subspace V
6: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
7: SLEPc - Scalable Library for Eigenvalue Problem Computations
8: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
10: This file is part of SLEPc.
12: SLEPc is free software: you can redistribute it and/or modify it under the
13: terms of version 3 of the GNU Lesser General Public License as published by
14: the Free Software Foundation.
16: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
17: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
19: more details.
21: You should have received a copy of the GNU Lesser General Public License
22: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
23: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24: */
26: #include davidson.h
28: typedef struct {
29: PetscInt k; /* desired initial subspace size */
30: PetscInt user; /* number of user initial vectors */
31: void *old_initV_data; /* old initV data */
32: } dvdInitV;
36: static PetscErrorCode dvd_initV_classic_0(dvdDashboard *d)
37: {
39: dvdInitV *data = (dvdInitV*)d->initV_data;
40: PetscInt i,user = PetscMin(data->user,d->eps->mpd), l,k;
43: BVGetActiveColumns(d->eps->V,&l,&k);
44: /* User vectors are added at the beginning, so no active column should be in V */
45: if (data->user>0&&l>0) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
46: /* Generate a set of random initial vectors and orthonormalize them */
47: for (i=l+user;i<l+data->k && i<d->eps->ncv && i-l<d->eps->mpd;i++) {
48: BVSetRandomColumn(d->eps->V,i);
49: }
50: d->V_tra_s = 0; d->V_tra_e = 0;
51: d->V_new_s = 0; d->V_new_e = i-l;
53: /* After that the user vectors will be destroyed */
54: data->user = 0;
55: return(0);
56: }
60: static PetscErrorCode dvd_initV_krylov_0(dvdDashboard *d)
61: {
63: dvdInitV *data = (dvdInitV*)d->initV_data;
64: PetscInt i,user = PetscMin(data->user,d->eps->mpd),l,k;
65: Vec av,v1,v2;
68: BVGetActiveColumns(d->eps->V,&l,&k);
69: /* User vectors are added at the beginning, so no active column should be in V */
70: if (data->user>0&&l>0) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
72: /* If needed, generate a random vector for starting the arnoldi method */
73: if (user == 0) {
74: BVSetRandomColumn(d->eps->V,l);
75: user = 1;
76: }
78: /* Perform k steps of Arnoldi with the operator K^{-1}*(t[1]*A-t[2]*B) */
79: dvd_orthV(d->eps->V,l,l+user);
80: for (i=l+user;i<l+data->k && i<d->eps->ncv && i-l<d->eps->mpd;i++) {
81: /* aux <- theta[1]A*in - theta[0]*B*in */
82: BVGetColumn(d->eps->V,i,&v1);
83: BVGetColumn(d->eps->V,i-user,&v2);
84: BVGetColumn(d->auxBV,0,&av);
85: if (d->B) {
86: MatMult(d->A,v2,v1);
87: MatMult(d->B,v2,av);
88: VecAXPBY(av,d->target[1],-d->target[0],v1);
89: } else {
90: MatMult(d->A,v2,av);
91: VecAXPBY(av,-d->target[0],d->target[1],v2);
92: }
93: d->improvex_precond(d,0,av,v1);
94: BVRestoreColumn(d->eps->V,i,&v1);
95: BVRestoreColumn(d->eps->V,i-user,&v2);
96: BVRestoreColumn(d->auxBV,0,&av);
97: dvd_orthV(d->eps->V,i,i+1);
98: }
100: d->V_tra_s = 0; d->V_tra_e = 0;
101: d->V_new_s = 0; d->V_new_e = i-l;
103: /* After that the user vectors will be destroyed */
104: data->user = 0;
105: return(0);
106: }
110: static PetscErrorCode dvd_initV_d(dvdDashboard *d)
111: {
113: dvdInitV *data = (dvdInitV*)d->initV_data;
116: /* Restore changes in dvdDashboard */
117: d->initV_data = data->old_initV_data;
119: /* Free local data */
120: PetscFree(data);
121: return(0);
122: }
126: PetscErrorCode dvd_initV(dvdDashboard *d, dvdBlackboard *b, PetscInt k,PetscInt user, PetscBool krylov)
127: {
129: dvdInitV *data;
132: /* Setting configuration constrains */
133: b->max_size_V = PetscMax(b->max_size_V, k);
135: /* Setup the step */
136: if (b->state >= DVD_STATE_CONF) {
137: PetscNewLog(d->eps,&data);
138: data->k = k;
139: data->user = user;
140: data->old_initV_data = d->initV_data;
141: d->initV_data = data;
142: if (krylov) d->initV = dvd_initV_krylov_0;
143: else d->initV = dvd_initV_classic_0;
144: EPSDavidsonFLAdd(&d->destroyList,dvd_initV_d);
145: }
146: return(0);
147: }
151: PetscErrorCode dvd_orthV(BV V,PetscInt V_new_s,PetscInt V_new_e)
152: {
154: PetscInt i,j,l,k;
155: PetscBool lindep;
156: PetscReal norm;
159: BVGetActiveColumns(V,&l,&k);
160: for (i=V_new_s;i<V_new_e;i++) {
161: BVSetActiveColumns(V,0,i);
162: for (j=0;j<3;j++) {
163: if (j>0) {
164: BVSetRandomColumn(V,i);
165: PetscInfo1(V,"Orthonormalization problems adding the vector %D to the searching subspace\n",i);
166: }
167: BVOrthogonalizeColumn(V,i,NULL,&norm,&lindep);
168: if (!lindep && (PetscAbsReal(norm) > PETSC_SQRT_MACHINE_EPSILON)) break;
169: }
170: if (lindep || (PetscAbsReal(norm) < PETSC_SQRT_MACHINE_EPSILON)) SETERRQ(PetscObjectComm((PetscObject)V),1, "Error during the orthonormalization of the vectors");
171: BVScaleColumn(V,i,1.0/norm);
172: }
173: BVSetActiveColumns(V,l,k);
174: return(0);
175: }