Actual source code: dvd_gd2.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  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-2014, 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: static PetscErrorCode dvd_improvex_gd2_d(dvdDashboard *d);
 29: static PetscErrorCode dvd_improvex_gd2_gen(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D);

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

 33: typedef struct {
 34:   PetscInt size_X;
 35: } dvdImprovex_gd2;

 39: PetscErrorCode dvd_improvex_gd2(dvdDashboard *d,dvdBlackboard *b,KSP ksp,PetscInt max_bs)
 40: {
 41:   PetscErrorCode  ierr;
 42:   dvdImprovex_gd2 *data;
 43:   PC              pc;


 47:   /* Setting configuration constrains */
 48:   /* If the arithmetic is real and the problem is not Hermitian, then
 49:      the block size is incremented in one */
 50: #if !defined(PETSC_USE_COMPLEX)
 51:   if (!DVD_IS(d->sEP, DVD_EP_HERMITIAN)) {
 52:     max_bs++;
 53:     b->max_size_P = PetscMax(b->max_size_P, 2);
 54:   } else
 55: #endif
 56:   {
 57:     b->max_size_P = PetscMax(b->max_size_P, 1);
 58:   }
 59:   b->max_size_X = PetscMax(b->max_size_X, max_bs);

 61:   /* Setup the preconditioner */
 62:   if (ksp) {
 63:     KSPGetPC(ksp,&pc);
 64:     dvd_static_precond_PC(d,b,pc);
 65:   } else {
 66:     dvd_static_precond_PC(d,b,0);
 67:   }

 69:   /* Setup the step */
 70:   if (b->state >= DVD_STATE_CONF) {
 71:     PetscMalloc(sizeof(dvdImprovex_gd2),&data);
 72:     PetscLogObjectMemory((PetscObject)d->eps,sizeof(dvdImprovex_gd2));
 73:     d->improveX_data = data;
 74:     data->size_X = b->max_size_X;
 75:     d->improveX = dvd_improvex_gd2_gen;

 77:     EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_gd2_d);
 78:   }
 79:   return(0);
 80: }

 84: static PetscErrorCode dvd_improvex_gd2_d(dvdDashboard *d)
 85: {
 86:   PetscErrorCode  ierr;
 87:   dvdImprovex_gd2 *data = (dvdImprovex_gd2*)d->improveX_data;

 90:   /* Free local data and objects */
 91:   PetscFree(data);
 92:   return(0);
 93: }

 97: static PetscErrorCode dvd_improvex_gd2_gen(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
 98: {
 99:   dvdImprovex_gd2 *data = (dvdImprovex_gd2*)d->improveX_data;
100:   PetscErrorCode  ierr;
101:   PetscInt        i,j,n,s,ld,k,lv,kv,max_size_D;
102:   PetscScalar     *pX,*b;
103:   Vec             *Ax,*Bx,v,*x;
104:   Mat             M;
105:   BV              X;

108:   /* Compute the number of pairs to improve */
109:   BVGetActiveColumns(d->eps->V,&lv,&kv);
110:   max_size_D = d->eps->ncv-kv;
111:   n = PetscMin(PetscMin(data->size_X*2,max_size_D),(r_e-r_s)*2)/2;
112: #if !defined(PETSC_USE_COMPLEX)
113:   /* If the last eigenvalue is a complex conjugate pair, n is increased by one */
114:   for (i=0; i<n; i++) {
115:     if (d->eigi[i] != 0.0) i++;
116:   }
117:   if (i > n) {
118:     n = PetscMin(PetscMin(data->size_X*2,max_size_D),(n+1)*2)/2;
119:     if (i > n) n--;
120:   }
121: #endif

123:   /* Quick exit */
124:   if (max_size_D == 0 || r_e-r_s <= 0 || n == 0) {
125:    *size_D = 0;
126:     return(0);
127:   }

129:   BVDuplicateResize(d->eps->V,4,&X);
130:   MatCreateSeqDense(PETSC_COMM_SELF,4,2,NULL,&M);

132:   /* Compute the eigenvectors of the selected pairs */
133:   for (i=0;i<n;) {
134:     k = r_s+i;
135:     DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
136:     DSNormalize(d->eps->ds,DS_MAT_X,r_s+i);
137:     /* Jump complex conjugate pairs */
138:     i = k+1;
139:   }
140:   DSGetArray(d->eps->ds,DS_MAT_X,&pX);
141:   DSGetLeadingDimension(d->eps->ds,&ld);

143:   SlepcVecPoolGetVecs(d->auxV,n,&Ax);
144:   SlepcVecPoolGetVecs(d->auxV,n,&Bx);

146:   /* Bx <- B*X(i) */
147:   if (d->BX) {
148:     /* Compute the norms of the eigenvectors */
149:     if (d->correctXnorm) {
150:       dvd_improvex_compute_X(d,r_s,r_s+n,Bx,pX,ld);
151:     } else {
152:       for (i=0; i<n; i++) d->nX[r_s+i] = 1.0;
153:     }
154:     for (i=0; i<n; ++i) {
155:       BVMultVec(d->BX,1.0,0.0,Bx[i],&pX[ld*(r_s+i)]);
156:     }
157:   } else if (d->B) {
158:     SlepcVecPoolGetVecs(d->auxV,1,&x);
159:     for (i=0;i<n;i++) {
160:       /* auxV(0) <- X(i) */
161:       dvd_improvex_compute_X(d,r_s+i,r_s+i+1,x,pX,ld);
162:       /* Bx(i) <- B*auxV(0) */
163:       MatMult(d->B,x[0],Bx[i]);
164:     }
165:     SlepcVecPoolRestoreVecs(d->auxV,1,&x);
166:   } else {
167:     /* Bx <- X */
168:     dvd_improvex_compute_X(d,r_s,r_s+n,Bx,pX,ld);
169:   }

171:   /* Ax <- A*X(i) */
172:   for (i=0; i<n; ++i) {
173:     BVMultVec(d->AX,1.0,0.0,Ax[i],&pX[ld*(i+r_s)]);
174:   }

176:   DSRestoreArray(d->eps->ds,DS_MAT_X,&pX);

178:   for (i=0,s=0;i<n;i+=s) {
179: #if !defined(PETSC_USE_COMPLEX)
180:     if (d->eigi[r_s+i] != 0.0) {
181:        /* [Ax_i Ax_i+1 Bx_i Bx_i+1]*= [   1        0
182:                                           0        1
183:                                        -eigr_i -eigi_i
184:                                         eigi_i -eigr_i] */
185:       MatDenseGetArray(M,&b);
186:       b[0] = b[5] = 1.0/d->nX[r_s+i];
187:       b[2] = b[7] = -d->eigr[r_s+i]/d->nX[r_s+i];
188:       b[6] = -(b[3] = d->eigi[r_s+i]/d->nX[r_s+i]);
189:       b[1] = b[4] = 0.0;
190:       MatDenseRestoreArray(M,&b);
191:       BVInsertVec(X,0,Ax[i]);
192:       BVInsertVec(X,1,Ax[i+1]);
193:       BVInsertVec(X,2,Bx[i]);
194:       BVInsertVec(X,3,Bx[i+1]);
195:       BVSetActiveColumns(X,0,4);
196:       BVMultInPlace(X,M,0,2);
197:       BVCopyVec(X,0,Ax[i]);
198:       BVCopyVec(X,1,Ax[i+1]);
199:       s = 2;
200:     } else
201: #endif
202:     {
203:       /* [Ax_i Bx_i]*= [ 1/nX_i    conj(eig_i/nX_i)
204:                        -eig_i/nX_i     1/nX_i       ] */
205:       MatDenseGetArray(M,&b);
206:       b[0] = 1.0/d->nX[r_s+i];
207:       b[1] = -d->eigr[r_s+i]/d->nX[r_s+i];
208:       b[4] = PetscConj(d->eigr[r_s+i]/d->nX[r_s+i]);
209:       b[5] = 1.0/d->nX[r_s+i];
210:       MatDenseRestoreArray(M,&b);
211:       BVInsertVec(X,0,Ax[i]);
212:       BVInsertVec(X,1,Bx[i]);
213:       BVSetActiveColumns(X,0,2);
214:       BVMultInPlace(X,M,0,2);
215:       BVCopyVec(X,0,Ax[i]);
216:       BVCopyVec(X,1,Bx[i]);
217:       s = 1;
218:     }
219:     for (j=0; j<s; j++) d->nX[r_s+i+j] = 1.0;

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

224:     /* Check if the first eigenpairs are converged */
225:     if (i == 0) {
226:       d->preTestConv(d,0,s,s,&d->npreconv);
227:       if (d->npreconv > 0) break;
228:     }
229:   }

231:   /* D <- K*[Ax Bx] */
232:   if (d->npreconv == 0) {
233:     for (i=0;i<n;i++) {
234:       BVGetColumn(d->eps->V,kv+i,&v);
235:       d->improvex_precond(d,r_s+i,Ax[i],v);
236:       BVRestoreColumn(d->eps->V,kv+i,&v);
237:     }
238:     for (i=n;i<n*2;i++) {
239:       BVGetColumn(d->eps->V,kv+i,&v);
240:       d->improvex_precond(d,r_s+i-n,Bx[i-n],v);
241:       BVRestoreColumn(d->eps->V,kv+i,&v);
242:     }
243:     *size_D = 2*n;
244: #if !defined(PETSC_USE_COMPLEX)
245:     if (d->eigi[r_s] != 0.0) {
246:       s = 4;
247:     } else
248: #endif
249:     {
250:       s = 2;
251:     }
252:     /* Prevent that short vectors are discarded in the orthogonalization */
253:     for (i=0; i<s && i<*size_D; i++) {
254:       if (d->eps->errest[d->nconv+r_s+i] > PETSC_MACHINE_EPSILON && d->eps->errest[d->nconv+r_s+i] < PETSC_MAX_REAL) {
255:         BVScaleColumn(d->eps->V,i+kv,1.0/d->eps->errest[d->nconv+r_s+i]);
256:       }
257:     }
258:   } else {
259:     *size_D = 0;
260:   }

262:   SlepcVecPoolRestoreVecs(d->auxV,n,&Bx);
263:   SlepcVecPoolRestoreVecs(d->auxV,n,&Ax);
264:   BVDestroy(&X);
265:   MatDestroy(&M);
266:   return(0);
267: }