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-2012, Universitat Politecnica de Valencia, Spain

 10:    This file is part of SLEPc.
 11:       
 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 <slepcvec.h>
 28: #include <slepcblaslapack.h>

 30: PetscErrorCode dvd_improvex_gd2_d(dvdDashboard *d);
 31: PetscErrorCode dvd_improvex_gd2_gen(dvdDashboard *d, Vec *D,
 32:                                    PetscInt max_size_D, PetscInt r_s,
 33:                                    PetscInt r_e, PetscInt *size_D);
 34: PetscErrorCode dvd_improvex_get_eigenvectors(dvdDashboard *d, PetscScalar *pX,
 35:   PetscScalar *pY, PetscInt ld_,
 36:   PetscScalar *auxS, PetscInt size_auxS);

 38: #define size_Z (64*4)

 40: /**** GD2 update step K*[A*X B*X]  ****/

 42: typedef struct {
 43:   PetscInt size_X;
 44:   void
 45:     *old_improveX_data;   /* old improveX_data */
 46:   improveX_type
 47:     old_improveX;         /* old improveX */
 48: } dvdImprovex_gd2;

 52: PetscErrorCode dvd_improvex_gd2(dvdDashboard *d,dvdBlackboard *b,KSP ksp,PetscInt max_bs)
 53: {
 54:   PetscErrorCode  ierr;
 55:   dvdImprovex_gd2 *data;
 56:   PetscBool       her_probl,std_probl;
 57:   PC              pc;
 58:   PetscInt        s=1;

 61:   std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
 62:   her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;

 64:   /* Setting configuration constrains */
 65:   /* If the arithmetic is real and the problem is not Hermitian, then
 66:      the block size is incremented in one */
 67: #if !defined(PETSC_USE_COMPLEX)
 68:   if (!her_probl) {
 69:     max_bs++;
 70:     b->max_size_P = PetscMax(b->max_size_P, 2);
 71:     s = 2;
 72:   } else
 73: #endif
 74:     b->max_size_P = PetscMax(b->max_size_P, 1);
 75:   b->max_size_X = PetscMax(b->max_size_X, max_bs);
 76:   b->max_size_auxV = PetscMax(b->max_size_auxV,
 77:      s +
 78:      ((her_probl || !d->eps->trueres)?1:PetscMax(s*2,b->max_size_cX_proj+b->max_size_X))); /* testConv */
 79: 
 80:   b->max_size_auxS = PetscMax(b->max_size_auxS,
 81:       (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 */

 83:   /* Setup the preconditioner */
 84:   if (ksp) {
 85:     KSPGetPC(ksp,&pc);
 86:     dvd_static_precond_PC(d,b,pc);
 87:   } else {
 88:     dvd_static_precond_PC(d,b,0);
 89:   }

 91:   /* Setup the step */
 92:   if (b->state >= DVD_STATE_CONF) {
 93:     PetscMalloc(sizeof(dvdImprovex_gd2),&data);
 94:     data->old_improveX_data = d->improveX_data;
 95:     d->improveX_data = data;
 96:     data->old_improveX = d->improveX;
 97:     data->size_X = b->max_size_X;
 98:     d->improveX = dvd_improvex_gd2_gen;

100:     DVD_FL_ADD(d->destroyList,dvd_improvex_gd2_d);
101:   }
102:   return(0);
103: }


108: PetscErrorCode dvd_improvex_gd2_d(dvdDashboard *d)
109: {
110:   PetscErrorCode  ierr;
111:   dvdImprovex_gd2 *data = (dvdImprovex_gd2*)d->improveX_data;

114: 
115:   /* Restore changes in dvdDashboard */
116:   d->improveX_data = data->old_improveX_data;

118:   /* Free local data and objects */
119:   PetscFree(data);

121:   return(0);
122: }

124: #define DVD_COMPLEX_RAYLEIGH_QUOTIENT(ur,ui,Axr,Axi,Bxr,Bxi,eigr,eigi,b,ierr)\
125: { \
126:   VecDot((Axr), (ur), &(b)[0]);  /* r*A*r */ \
127:   VecDot((Axr), (ui), &(b)[1]);  /* i*A*r */ \
128:   VecDot((Axi), (ur), &(b)[2]);  /* r*A*i */ \
129:   VecDot((Axi), (ui), &(b)[3]);  /* i*A*i */ \
130:   VecDot((Bxr), (ur), &(b)[4]);  /* r*B*r */ \
131:   VecDot((Bxr), (ui), &(b)[5]);  /* i*B*r */ \
132:   VecDot((Bxi), (ur), &(b)[6]);  /* r*B*i */ \
133:   VecDot((Bxi), (ui), &(b)[7]);  /* i*B*i */ \
134:   (b)[0]  = (b)[0]+(b)[3]; /* rAr+iAi */ \
135:   (b)[2] =  (b)[2]-(b)[1]; /* rAi-iAr */ \
136:   (b)[4] = (b)[4]+(b)[7]; /* rBr+iBi */ \
137:   (b)[6] = (b)[6]-(b)[5]; /* rBi-iBr */ \
138:   (b)[7] = (b)[4]*(b)[4] + (b)[6]*(b)[6]; /* k */ \
139:   *(eigr) = ((b)[0]*(b)[4] + (b)[2]*(b)[6]) / (b)[7]; /* eig_r */ \
140:   *(eigi) = ((b)[2]*(b)[4] - (b)[0]*(b)[6]) / (b)[7]; /* eig_i */ \
141: }

143: #if !defined(PETSC_USE_COMPLEX)
144: #define DVD_COMPUTE_N_RR(eps,i,i_s,n,eigr,eigi,u,Ax,Bx,b,ierr) \
145:   for((i)=0; (i)<(n); (i)++) { \
146:     if ((eigi)[(i_s)+(i)] != 0.0) { \
147:       /* eig_r = [(rAr+iAi)*(rBr+iBi) + (rAi-iAr)*(rBi-iBr)]/k \
148:          eig_i = [(rAi-iAr)*(rBr+iBi) - (rAr+iAi)*(rBi-iBr)]/k \
149:          k     =  (rBr+iBi)*(rBr+iBi) + (rBi-iBr)*(rBi-iBr)    */ \
150:       DVD_COMPLEX_RAYLEIGH_QUOTIENT((u)[(i)], (u)[(i)+1], (Ax)[(i)], \
151:         (Ax)[(i)+1], (Bx)[(i)], (Bx)[(i)+1], &(b)[8], &(b)[9], (b), (ierr)); \
152:       if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[8])/ \
153:             PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10    || \
154:           PetscAbsScalar((eigi)[(i_s)+(i)] - (b)[9])/ \
155:             PetscAbsScalar((eigi)[(i_s)+(i)]) > 1e-10         ) { \
156:         (ierr) = PetscInfo4((eps), "The eigenvalue %G+%G is far from its "\
157:                             "Rayleigh quotient value %G+%G\n", \
158:                             (eigr)[(i_s)+(i)], \
159:                             (eigi)[(i_s)+1], (b)[8], (b)[9]); \
160:       } \
161:       (i)++; \
162:     } else { \
163:       (ierr) = VecDot((Ax)[(i)], (u)[(i)], &(b)[0]);  \
164:       (ierr) = VecDot((Bx)[(i)], (u)[(i)], &(b)[1]);  \
165:       (b)[0] = (b)[0]/(b)[1]; \
166:       if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[0])/ \
167:             PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10     ) { \
168:         (ierr) = PetscInfo3((eps), "The eigenvalue %G is far from its " \
169:                "Rayleigh quotient value %G. (y'*B*x = %G)\n", \
170:                (eigr)[(i_s)+(i)], \
171:                (b)[0], (b)[1]); \
172:       } \
173:     } \
174:   }
175: #else
176: #define DVD_COMPUTE_N_RR(eps,i,i_s,n,eigr,eigi,u,Ax,Bx,b,ierr) \
177:   for((i)=0; (i)<(n); (i)++) { \
178:       (ierr) = VecDot((Ax)[(i)], (u)[(i)], &(b)[0]);  \
179:       (ierr) = VecDot((Bx)[(i)], (u)[(i)], &(b)[1]);  \
180:       (b)[0] = (b)[0]/(b)[1]; \
181:       if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[0])/ \
182:             PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10     ) { \
183:         (ierr) = PetscInfo4((eps), "The eigenvalue %G+%G is far from its " \
184:                "Rayleigh quotient value %G+%G\n", \
185:                PetscRealPart((eigr)[(i_s)+(i)]), \
186:                PetscImaginaryPart((eigr)[(i_s)+(i)]), PetscRealPart((b)[0]), \
187:                PetscImaginaryPart((b)[0])); \
188:       } \
189:     }
190: #endif



196: PetscErrorCode dvd_improvex_gd2_gen(dvdDashboard *d,Vec *D,PetscInt max_size_D,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
197: {
198:   dvdImprovex_gd2 *data = (dvdImprovex_gd2*)d->improveX_data;
199:   PetscErrorCode  ierr;
200:   PetscInt        i,j,n,s,ld,k;
201:   PetscScalar     *pX,*pY,b[10],Z[size_Z];
202:   Vec             *Ax,*Bx,X[4];


206:   /* Compute the number of pairs to improve */
207:   n = PetscMin(PetscMin(PetscMin(data->size_X*2,max_size_D),(r_e-r_s)*2),d->max_size_proj-d->size_H)/2;
208: #if !defined(PETSC_USE_COMPLEX)
209:   /* If the last eigenvalue is a complex conjugate pair, n is increased by one */
210:   for (i=0; i<n; i++) {
211:     if (d->eigi[i] != 0.0) i++;
212:   }
213:   if (i > n) {
214:     n = PetscMin(PetscMin(PetscMin(data->size_X*2,max_size_D),(n+1)*2),d->max_size_proj-d->size_H)/2;
215:     if (i > n) n--;
216:   }
217: #endif

219:   /* Quick exit */
220:   if (max_size_D == 0 || r_e-r_s <= 0 || n == 0) {
221:    *size_D = 0;
222:    /* Callback old improveX */
223:     if (data->old_improveX) {
224:       d->improveX_data = data->old_improveX_data;
225:       data->old_improveX(d,PETSC_NULL,0,0,0,PETSC_NULL);
226:       d->improveX_data = data;
227:     }
228:     return(0);
229:   }

231:   /* Compute the eigenvectors of the selected pairs */
232:   for (i=0; i<n; ) {
233:     k = r_s+i+d->cX_in_H;
234:     DSVectors(d->ps,DS_MAT_X,&k,PETSC_NULL);
235:     DSNormalize(d->ps,DS_MAT_X,r_s+i+d->cX_in_H);
236:     k = r_s+i+d->cX_in_H;
237:     DSVectors(d->ps,DS_MAT_Y,&k,PETSC_NULL);
238:     DSNormalize(d->ps,DS_MAT_Y,r_s+i+d->cX_in_H);
239:     /* Jump complex conjugate pairs */
240:     i = k+1;
241:   }
242:   DSGetArray(d->ps,DS_MAT_X,&pX);
243:   DSGetArray(d->ps,DS_MAT_Y,&pY);
244:   DSGetLeadingDimension(d->ps,&ld);

246:   /* Bx <- B*X(i) */
247:   Bx = D+n;
248:   if (d->BV) {
249:     /* Compute the norms of the eigenvectors */
250:     if (d->correctXnorm) {
251:       dvd_improvex_compute_X(d,r_s,r_s+n,Bx,pX,ld);
252:     } else {
253:       for (i=0; i<n; i++) d->nX[r_s+i] = 1.0;
254:     }
255:     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);
256:   } else if (d->B) {
257:     for(i=0; i<n; i++) {
258:       /* auxV(0) <- X(i) */
259:       dvd_improvex_compute_X(d,r_s+i,r_s+i+1,d->auxV,pX,ld);
260:       /* Bx(i) <- B*auxV(0) */
261:       MatMult(d->B,d->auxV[0],Bx[i]);
262:     }
263:   } else {
264:     /* Bx <- X */
265:     dvd_improvex_compute_X(d,r_s,r_s+n,Bx,pX,ld);
266:   }

268:   /* Ax <- A*X(i) */
269:   Ax = D;
270:   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);

272: #if !defined(PETSC_USE_COMPLEX)
273:   s = d->eigi[r_s] == 0.0 ? 1 : 2;
274:   /* If the available vectors allow the computation of the eigenvalue */
275:   if (s <= n) {
276: #else
277:   s = 1;
278: #endif
279:   /* v <- Y(i) */
280:   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);

282:   /* Recompute the eigenvalue */
283:   DVD_COMPUTE_N_RR(d->eps,i,r_s,1,d->eigr,d->eigi,d->auxV,Ax,Bx,b,ierr);
284: #if !defined(PETSC_USE_COMPLEX)
285:   }
286: #endif

288:   DSRestoreArray(d->ps,DS_MAT_X,&pX);
289:   DSRestoreArray(d->ps,DS_MAT_Y,&pY);

291:   for(i=0,s=0; i<n; i+=s) {
292: #if !defined(PETSC_USE_COMPLEX)
293:     if (d->eigi[r_s+i] != 0.0) {
294:        /* [Ax_i Ax_i+1 Bx_i Bx_i+1]*= [   1        0 
295:                                           0        1
296:                                        -eigr_i -eigi_i
297:                                         eigi_i -eigr_i] */
298:       b[0] = b[5] = 1.0/d->nX[r_s+i];
299:       b[2] = b[7] = -d->eigr[r_s+i]/d->nX[r_s+i];
300:       b[6] = -(b[3] = d->eigi[r_s+i]/d->nX[r_s+i]);
301:       b[1] = b[4] = 0.0;
302:       X[0] = Ax[i]; X[1] = Ax[i+1]; X[2] = Bx[i]; X[3] = Bx[i+1];
303:       SlepcUpdateVectorsD(X,4,1.0,b,4,4,2,Z,size_Z);
304:       s = 2;
305:     } else
306: #endif
307:     {
308:       /* [Ax_i Bx_i]*= [ 1/nX_i    conj(eig_i/nX_i)
309:                        -eig_i/nX_i     1/nX_i       ] */
310:       b[0] = 1.0/d->nX[r_s+i];
311:       b[1] = -d->eigr[r_s+i]/d->nX[r_s+i];
312:       b[2] = PetscConj(d->eigr[r_s+i]/d->nX[r_s+i]);
313:       b[3] = 1.0/d->nX[r_s+i];
314:       X[0] = Ax[i]; X[1] = Bx[i];
315:       SlepcUpdateVectorsD(X,2,1.0,b,2,2,2,Z,size_Z);
316:       s = 1;
317:     }
318:     for (j=0; j<s; j++) d->nX[r_s+i+j] = 1.0;

320:     /* Ax = R <- P*(Ax - eig_i*Bx) */
321:     d->calcpairs_proj_res(d,r_s+i,r_s+i+s,&Ax[i]);

323:     /* Check if the first eigenpairs are converged */
324:     if (i == 0) {
325:       d->preTestConv(d,0,s,s,Ax,PETSC_NULL,&d->npreconv);
326:       if (d->npreconv > 0) break;
327:     }
328:   }
329: 
330:   /* D <- K*[Ax Bx] */
331:   if (d->npreconv == 0) {
332:     VecCopy(D[0],d->auxV[0]);
333:     for(i=0; i<2*n-1; i++) {
334:       d->improvex_precond(d,r_s+(i+1)%n,D[i+1],D[i]);
335:     }
336:     d->improvex_precond(d,r_s,d->auxV[0],D[2*n-1]);
337:     *size_D = 2*n;
338: #if !defined(PETSC_USE_COMPLEX)
339:     if (d->eigi[r_s] != 0.0) {
340:       s = 4;
341:     } else
342: #endif
343:     s = 2;
344:     /* Prevent that short vectors are discarded in the orthogonalization */
345:     if (d->eps->errest[d->nconv+r_s] > PETSC_MACHINE_EPSILON && d->eps->errest[d->nconv+r_s] < PETSC_MAX_REAL) {
346:       for (i=0; i<s && i<*size_D; i++) {
347:         VecScale(D[i],1.0/d->eps->errest[d->nconv+r_s]);
348:       }
349:     }
350:   } else {
351:     *size_D = 0;
352:   }
353: 
354:   /* Callback old improveX */
355:   if (data->old_improveX) {
356:     d->improveX_data = data->old_improveX_data;
357:     data->old_improveX(d,PETSC_NULL,0,0,0,PETSC_NULL);
358:     d->improveX_data = data;
359:   }

361:   return(0);
362: }