Actual source code: dvd_gd2.c
1: /*
2: SLEPc eigensolver: "davidson"
4: Step: improve the eigenvectors X with GD2
6: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
7: SLEPc - Scalable Library for Eigenvalue Problem Computations
8: Copyright (c) 2002-2013, 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
27: #include <slepc-private/vecimplslepc.h> /*I "slepcvec.h" I*/
29: PetscErrorCode dvd_improvex_gd2_d(dvdDashboard *d);
30: PetscErrorCode dvd_improvex_gd2_gen(dvdDashboard *d,Vec *D,PetscInt max_size_D,PetscInt r_s,PetscInt r_e,PetscInt *size_D);
31: PetscErrorCode dvd_improvex_get_eigenvectors(dvdDashboard *d,PetscScalar *pX,PetscScalar *pY,PetscInt ld_,PetscScalar *auxS,PetscInt size_auxS);
33: #define size_Z (64*4)
35: /**** GD2 update step K*[A*X B*X] ****/
37: typedef struct {
38: PetscInt size_X;
39: void *old_improveX_data; /* old improveX_data */
40: improveX_type old_improveX; /* old improveX */
41: } dvdImprovex_gd2;
45: PetscErrorCode dvd_improvex_gd2(dvdDashboard *d,dvdBlackboard *b,KSP ksp,PetscInt max_bs)
46: {
47: PetscErrorCode ierr;
48: dvdImprovex_gd2 *data;
49: PetscBool her_probl,std_probl;
50: PC pc;
51: PetscInt s=1;
54: std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
55: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
57: /* Setting configuration constrains */
58: /* If the arithmetic is real and the problem is not Hermitian, then
59: the block size is incremented in one */
60: #if !defined(PETSC_USE_COMPLEX)
61: if (!her_probl) {
62: max_bs++;
63: b->max_size_P = PetscMax(b->max_size_P, 2);
64: s = 2;
65: } else
66: #endif
67: b->max_size_P = PetscMax(b->max_size_P, 1);
68: b->max_size_X = PetscMax(b->max_size_X, max_bs);
69: b->max_size_auxV = PetscMax(b->max_size_auxV,
70: s +
71: ((her_probl || !d->eps->trueres)?1:PetscMax(s*2,b->max_size_cX_proj+b->max_size_X))); /* testConv */
73: b->max_size_auxS = PetscMax(b->max_size_auxS,
74: (her_probl || !d->eps->trueres)?0:b->max_nev*b->max_nev+PetscMax(b->max_nev*6,(b->max_nev+b->max_size_proj)*s+b->max_nev*(b->max_size_X+b->max_size_cX_proj)*(std_probl?2:4)+64)); /* preTestConv */
76: /* Setup the preconditioner */
77: if (ksp) {
78: KSPGetPC(ksp,&pc);
79: dvd_static_precond_PC(d,b,pc);
80: } else {
81: dvd_static_precond_PC(d,b,0);
82: }
84: /* Setup the step */
85: if (b->state >= DVD_STATE_CONF) {
86: PetscMalloc(sizeof(dvdImprovex_gd2),&data);
87: PetscLogObjectMemory(d->eps,sizeof(dvdImprovex_gd2));
88: data->old_improveX_data = d->improveX_data;
89: d->improveX_data = data;
90: data->old_improveX = d->improveX;
91: data->size_X = b->max_size_X;
92: d->improveX = dvd_improvex_gd2_gen;
94: DVD_FL_ADD(d->destroyList,dvd_improvex_gd2_d);
95: }
96: return(0);
97: }
101: PetscErrorCode dvd_improvex_gd2_d(dvdDashboard *d)
102: {
103: PetscErrorCode ierr;
104: dvdImprovex_gd2 *data = (dvdImprovex_gd2*)d->improveX_data;
107: /* Restore changes in dvdDashboard */
108: d->improveX_data = data->old_improveX_data;
110: /* Free local data and objects */
111: PetscFree(data);
112: return(0);
113: }
115: #define DVD_COMPLEX_RAYLEIGH_QUOTIENT(ur,ui,Axr,Axi,Bxr,Bxi,eigr,eigi,b,ierr)\
116: { \
117: VecDot((Axr), (ur), &(b)[0]); /* r*A*r */ \
118: VecDot((Axr), (ui), &(b)[1]); /* i*A*r */ \
119: VecDot((Axi), (ur), &(b)[2]); /* r*A*i */ \
120: VecDot((Axi), (ui), &(b)[3]); /* i*A*i */ \
121: VecDot((Bxr), (ur), &(b)[4]); /* r*B*r */ \
122: VecDot((Bxr), (ui), &(b)[5]); /* i*B*r */ \
123: VecDot((Bxi), (ur), &(b)[6]); /* r*B*i */ \
124: VecDot((Bxi), (ui), &(b)[7]); /* i*B*i */ \
125: (b)[0] = (b)[0]+(b)[3]; /* rAr+iAi */ \
126: (b)[2] = (b)[2]-(b)[1]; /* rAi-iAr */ \
127: (b)[4] = (b)[4]+(b)[7]; /* rBr+iBi */ \
128: (b)[6] = (b)[6]-(b)[5]; /* rBi-iBr */ \
129: (b)[7] = (b)[4]*(b)[4] + (b)[6]*(b)[6]; /* k */ \
130: *(eigr) = ((b)[0]*(b)[4] + (b)[2]*(b)[6]) / (b)[7]; /* eig_r */ \
131: *(eigi) = ((b)[2]*(b)[4] - (b)[0]*(b)[6]) / (b)[7]; /* eig_i */ \
132: }
134: #if !defined(PETSC_USE_COMPLEX)
135: #define DVD_COMPUTE_N_RR(eps,i,i_s,n,eigr,eigi,u,Ax,Bx,b,ierr) \
136: for ((i)=0; (i)<(n); (i)++) { \
137: if ((eigi)[(i_s)+(i)] != 0.0) { \
138: /* eig_r = [(rAr+iAi)*(rBr+iBi) + (rAi-iAr)*(rBi-iBr)]/k \
139: eig_i = [(rAi-iAr)*(rBr+iBi) - (rAr+iAi)*(rBi-iBr)]/k \
140: k = (rBr+iBi)*(rBr+iBi) + (rBi-iBr)*(rBi-iBr) */ \
141: DVD_COMPLEX_RAYLEIGH_QUOTIENT((u)[(i)], (u)[(i)+1], (Ax)[(i)], \
142: (Ax)[(i)+1], (Bx)[(i)], (Bx)[(i)+1], &(b)[8], &(b)[9], (b), (ierr)); \
143: if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[8])/ \
144: PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10 || \
145: PetscAbsScalar((eigi)[(i_s)+(i)] - (b)[9])/ \
146: PetscAbsScalar((eigi)[(i_s)+(i)]) > 1e-10) { \
147: (ierr) = PetscInfo4((eps), "The eigenvalue %G+%G is far from its "\
148: "Rayleigh quotient value %G+%G\n", \
149: (eigr)[(i_s)+(i)], \
150: (eigi)[(i_s)+1], (b)[8], (b)[9]); \
151: } \
152: (i)++; \
153: } else { \
154: (ierr) = VecDot((Ax)[(i)], (u)[(i)], &(b)[0]); \
155: (ierr) = VecDot((Bx)[(i)], (u)[(i)], &(b)[1]); \
156: (b)[0] = (b)[0]/(b)[1]; \
157: if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[0])/ \
158: PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10) { \
159: (ierr) = PetscInfo3((eps), "The eigenvalue %G is far from its " \
160: "Rayleigh quotient value %G. (y'*B*x = %G)\n", \
161: (eigr)[(i_s)+(i)], \
162: (b)[0], (b)[1]); \
163: } \
164: } \
165: }
166: #else
167: #define DVD_COMPUTE_N_RR(eps,i,i_s,n,eigr,eigi,u,Ax,Bx,b,ierr) \
168: for ((i)=0; (i)<(n); (i)++) { \
169: (ierr) = VecDot((Ax)[(i)], (u)[(i)], &(b)[0]); \
170: (ierr) = VecDot((Bx)[(i)], (u)[(i)], &(b)[1]); \
171: (b)[0] = (b)[0]/(b)[1]; \
172: if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[0])/ \
173: PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10) { \
174: (ierr) = PetscInfo4((eps), "The eigenvalue %G+%G is far from its " \
175: "Rayleigh quotient value %G+%G\n", \
176: PetscRealPart((eigr)[(i_s)+(i)]), \
177: PetscImaginaryPart((eigr)[(i_s)+(i)]), PetscRealPart((b)[0]), \
178: PetscImaginaryPart((b)[0])); \
179: } \
180: }
181: #endif
185: PetscErrorCode dvd_improvex_gd2_gen(dvdDashboard *d,Vec *D,PetscInt max_size_D,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
186: {
187: dvdImprovex_gd2 *data = (dvdImprovex_gd2*)d->improveX_data;
188: PetscErrorCode ierr;
189: PetscInt i,j,n,s,ld,k;
190: PetscScalar *pX,*pY,b[10],Z[size_Z];
191: Vec *Ax,*Bx,X[4];
194: /* Compute the number of pairs to improve */
195: n = PetscMin(PetscMin(PetscMin(data->size_X*2,max_size_D),(r_e-r_s)*2),d->max_size_proj-d->size_H)/2;
196: #if !defined(PETSC_USE_COMPLEX)
197: /* If the last eigenvalue is a complex conjugate pair, n is increased by one */
198: for (i=0; i<n; i++) {
199: if (d->eigi[i] != 0.0) i++;
200: }
201: if (i > n) {
202: n = PetscMin(PetscMin(PetscMin(data->size_X*2,max_size_D),(n+1)*2),d->max_size_proj-d->size_H)/2;
203: if (i > n) n--;
204: }
205: #endif
207: /* Quick exit */
208: if (max_size_D == 0 || r_e-r_s <= 0 || n == 0) {
209: *size_D = 0;
210: /* Callback old improveX */
211: if (data->old_improveX) {
212: d->improveX_data = data->old_improveX_data;
213: data->old_improveX(d,NULL,0,0,0,NULL);
214: d->improveX_data = data;
215: }
216: return(0);
217: }
219: /* Compute the eigenvectors of the selected pairs */
220: for (i=0;i<n;) {
221: k = r_s+i+d->cX_in_H;
222: DSVectors(d->ps,DS_MAT_X,&k,NULL);
223: DSNormalize(d->ps,DS_MAT_X,r_s+i+d->cX_in_H);
224: k = r_s+i+d->cX_in_H;
225: DSVectors(d->ps,DS_MAT_Y,&k,NULL);
226: DSNormalize(d->ps,DS_MAT_Y,r_s+i+d->cX_in_H);
227: /* Jump complex conjugate pairs */
228: i = k+1;
229: }
230: DSGetArray(d->ps,DS_MAT_X,&pX);
231: DSGetArray(d->ps,DS_MAT_Y,&pY);
232: DSGetLeadingDimension(d->ps,&ld);
234: /* Bx <- B*X(i) */
235: Bx = D+n;
236: if (d->BV) {
237: /* Compute the norms of the eigenvectors */
238: if (d->correctXnorm) {
239: dvd_improvex_compute_X(d,r_s,r_s+n,Bx,pX,ld);
240: } else {
241: for (i=0; i<n; i++) d->nX[r_s+i] = 1.0;
242: }
243: SlepcUpdateVectorsZ(Bx,0.0,1.0,d->BV-d->cX_in_H,d->size_BV+d->cX_in_H,&pX[ld*r_s],ld,d->size_H,n);
244: } else if (d->B) {
245: for (i=0;i<n;i++) {
246: /* auxV(0) <- X(i) */
247: dvd_improvex_compute_X(d,r_s+i,r_s+i+1,d->auxV,pX,ld);
248: /* Bx(i) <- B*auxV(0) */
249: MatMult(d->B,d->auxV[0],Bx[i]);
250: }
251: } else {
252: /* Bx <- X */
253: dvd_improvex_compute_X(d,r_s,r_s+n,Bx,pX,ld);
254: }
256: /* Ax <- A*X(i) */
257: Ax = D;
258: SlepcUpdateVectorsZ(Ax,0.0,1.0,d->AV-d->cX_in_H,d->size_AV+d->cX_in_H,&pX[ld*r_s],ld,d->size_H,n);
260: #if !defined(PETSC_USE_COMPLEX)
261: s = d->eigi[r_s] == 0.0 ? 1 : 2;
262: /* If the available vectors allow the computation of the eigenvalue */
263: if (s <= n) {
264: #else
265: s = 1;
266: #endif
267: /* v <- Y(i) */
268: SlepcUpdateVectorsZ(d->auxV,0.0,1.0,(d->W?d->W:d->V)-d->cX_in_H,d->size_V+d->cX_in_H,&pY[ld*r_s],ld,d->size_H,s);
270: /* Recompute the eigenvalue */
271: DVD_COMPUTE_N_RR(d->eps,i,r_s,1,d->eigr,d->eigi,d->auxV,Ax,Bx,b,ierr);
272: #if !defined(PETSC_USE_COMPLEX)
273: }
274: #endif
276: DSRestoreArray(d->ps,DS_MAT_X,&pX);
277: DSRestoreArray(d->ps,DS_MAT_Y,&pY);
279: for (i=0,s=0;i<n;i+=s) {
280: #if !defined(PETSC_USE_COMPLEX)
281: if (d->eigi[r_s+i] != 0.0) {
282: /* [Ax_i Ax_i+1 Bx_i Bx_i+1]*= [ 1 0
283: 0 1
284: -eigr_i -eigi_i
285: eigi_i -eigr_i] */
286: b[0] = b[5] = 1.0/d->nX[r_s+i];
287: b[2] = b[7] = -d->eigr[r_s+i]/d->nX[r_s+i];
288: b[6] = -(b[3] = d->eigi[r_s+i]/d->nX[r_s+i]);
289: b[1] = b[4] = 0.0;
290: X[0] = Ax[i]; X[1] = Ax[i+1]; X[2] = Bx[i]; X[3] = Bx[i+1];
291: SlepcUpdateVectorsD(X,4,1.0,b,4,4,2,Z,size_Z);
292: s = 2;
293: } else
294: #endif
295: {
296: /* [Ax_i Bx_i]*= [ 1/nX_i conj(eig_i/nX_i)
297: -eig_i/nX_i 1/nX_i ] */
298: b[0] = 1.0/d->nX[r_s+i];
299: b[1] = -d->eigr[r_s+i]/d->nX[r_s+i];
300: b[2] = PetscConj(d->eigr[r_s+i]/d->nX[r_s+i]);
301: b[3] = 1.0/d->nX[r_s+i];
302: X[0] = Ax[i]; X[1] = Bx[i];
303: SlepcUpdateVectorsD(X,2,1.0,b,2,2,2,Z,size_Z);
304: s = 1;
305: }
306: for (j=0; j<s; j++) d->nX[r_s+i+j] = 1.0;
308: /* Ax = R <- P*(Ax - eig_i*Bx) */
309: d->calcpairs_proj_res(d,r_s+i,r_s+i+s,&Ax[i]);
311: /* Check if the first eigenpairs are converged */
312: if (i == 0) {
313: d->preTestConv(d,0,s,s,Ax,NULL,&d->npreconv);
314: if (d->npreconv > 0) break;
315: }
316: }
318: /* D <- K*[Ax Bx] */
319: if (d->npreconv == 0) {
320: VecCopy(D[0],d->auxV[0]);
321: for (i=0;i<2*n-1;i++) {
322: d->improvex_precond(d,r_s+(i+1)%n,D[i+1],D[i]);
323: }
324: d->improvex_precond(d,r_s,d->auxV[0],D[2*n-1]);
325: *size_D = 2*n;
326: #if !defined(PETSC_USE_COMPLEX)
327: if (d->eigi[r_s] != 0.0) {
328: s = 4;
329: } else
330: #endif
331: s = 2;
332: /* Prevent that short vectors are discarded in the orthogonalization */
333: if (d->eps->errest[d->nconv+r_s] > PETSC_MACHINE_EPSILON && d->eps->errest[d->nconv+r_s] < PETSC_MAX_REAL) {
334: for (i=0; i<s && i<*size_D; i++) {
335: VecScale(D[i],1.0/d->eps->errest[d->nconv+r_s]);
336: }
337: }
338: } else {
339: *size_D = 0;
340: }
342: /* Callback old improveX */
343: if (data->old_improveX) {
344: d->improveX_data = data->old_improveX_data;
345: data->old_improveX(d,NULL,0,0,0,NULL);
346: d->improveX_data = data;
347: }
348: return(0);
349: }