Actual source code: dvdgd2.c

slepc-3.6.1 2015-09-03
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-2015, 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 size_X;
 30: } dvdImprovex_gd2;

 34: static PetscErrorCode dvd_improvex_gd2_d(dvdDashboard *d)
 35: {
 36:   PetscErrorCode  ierr;
 37:   dvdImprovex_gd2 *data = (dvdImprovex_gd2*)d->improveX_data;

 40:   /* Free local data and objects */
 41:   PetscFree(data);
 42:   return(0);
 43: }

 47: static PetscErrorCode dvd_improvex_gd2_gen(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
 48: {
 49:   dvdImprovex_gd2 *data = (dvdImprovex_gd2*)d->improveX_data;
 50:   PetscErrorCode  ierr;
 51:   PetscInt        i,j,n,s,ld,k,lv,kv,max_size_D;
 52:   PetscScalar     *pX,*b;
 53:   Vec             *Ax,*Bx,v,*x;
 54:   Mat             M;
 55:   BV              X;

 58:   /* Compute the number of pairs to improve */
 59:   BVGetActiveColumns(d->eps->V,&lv,&kv);
 60:   max_size_D = d->eps->ncv-kv;
 61:   n = PetscMin(PetscMin(data->size_X*2,max_size_D),(r_e-r_s)*2)/2;
 62: #if !defined(PETSC_USE_COMPLEX)
 63:   /* If the last eigenvalue is a complex conjugate pair, n is increased by one */
 64:   for (i=0; i<n; i++) {
 65:     if (d->eigi[i] != 0.0) i++;
 66:   }
 67:   if (i > n) {
 68:     n = PetscMin(PetscMin(data->size_X*2,max_size_D),(n+1)*2)/2;
 69:     if (i > n) n--;
 70:   }
 71: #endif

 73:   /* Quick exit */
 74:   if (max_size_D == 0 || r_e-r_s <= 0 || n == 0) {
 75:    *size_D = 0;
 76:     return(0);
 77:   }

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

 82:   /* Compute the eigenvectors of the selected pairs */
 83:   for (i=0;i<n;) {
 84:     k = r_s+i;
 85:     DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
 86:     DSNormalize(d->eps->ds,DS_MAT_X,r_s+i);
 87:     i = k+1; /* skip complex conjugate pairs */
 88:   }
 89:   DSGetArray(d->eps->ds,DS_MAT_X,&pX);
 90:   DSGetLeadingDimension(d->eps->ds,&ld);

 92:   SlepcVecPoolGetVecs(d->auxV,n,&Ax);
 93:   SlepcVecPoolGetVecs(d->auxV,n,&Bx);

 95:   /* Bx <- B*X(i) */
 96:   if (d->BX) {
 97:     /* Compute the norms of the eigenvectors */
 98:     if (d->correctXnorm) {
 99:       dvd_improvex_compute_X(d,r_s,r_s+n,Bx,pX,ld);
100:     } else {
101:       for (i=0;i<n;i++) d->nX[r_s+i] = 1.0;
102:     }
103:     for (i=0;i<n;i++) {
104:       BVMultVec(d->BX,1.0,0.0,Bx[i],&pX[ld*(r_s+i)]);
105:     }
106:   } else if (d->B) {
107:     SlepcVecPoolGetVecs(d->auxV,1,&x);
108:     for (i=0;i<n;i++) {
109:       /* auxV(0) <- X(i) */
110:       dvd_improvex_compute_X(d,r_s+i,r_s+i+1,x,pX,ld);
111:       /* Bx(i) <- B*auxV(0) */
112:       MatMult(d->B,x[0],Bx[i]);
113:     }
114:     SlepcVecPoolRestoreVecs(d->auxV,1,&x);
115:   } else {
116:     /* Bx <- X */
117:     dvd_improvex_compute_X(d,r_s,r_s+n,Bx,pX,ld);
118:   }

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

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

127:   for (i=0,s=0;i<n;i+=s) {
128: #if !defined(PETSC_USE_COMPLEX)
129:     if (d->eigi[r_s+i] != 0.0) {
130:        /* [Ax_i Ax_i+1 Bx_i Bx_i+1]*= [   1        0
131:                                           0        1
132:                                        -eigr_i -eigi_i
133:                                         eigi_i -eigr_i] */
134:       MatDenseGetArray(M,&b);
135:       b[0] = b[5] = 1.0/d->nX[r_s+i];
136:       b[2] = b[7] = -d->eigr[r_s+i]/d->nX[r_s+i];
137:       b[6] = -(b[3] = d->eigi[r_s+i]/d->nX[r_s+i]);
138:       b[1] = b[4] = 0.0;
139:       MatDenseRestoreArray(M,&b);
140:       BVInsertVec(X,0,Ax[i]);
141:       BVInsertVec(X,1,Ax[i+1]);
142:       BVInsertVec(X,2,Bx[i]);
143:       BVInsertVec(X,3,Bx[i+1]);
144:       BVSetActiveColumns(X,0,4);
145:       BVMultInPlace(X,M,0,2);
146:       BVCopyVec(X,0,Ax[i]);
147:       BVCopyVec(X,1,Ax[i+1]);
148:       s = 2;
149:     } else
150: #endif
151:     {
152:       /* [Ax_i Bx_i]*= [ 1/nX_i    conj(eig_i/nX_i)
153:                        -eig_i/nX_i     1/nX_i       ] */
154:       MatDenseGetArray(M,&b);
155:       b[0] = 1.0/d->nX[r_s+i];
156:       b[1] = -d->eigr[r_s+i]/d->nX[r_s+i];
157:       b[4] = PetscConj(d->eigr[r_s+i]/d->nX[r_s+i]);
158:       b[5] = 1.0/d->nX[r_s+i];
159:       MatDenseRestoreArray(M,&b);
160:       BVInsertVec(X,0,Ax[i]);
161:       BVInsertVec(X,1,Bx[i]);
162:       BVSetActiveColumns(X,0,2);
163:       BVMultInPlace(X,M,0,2);
164:       BVCopyVec(X,0,Ax[i]);
165:       BVCopyVec(X,1,Bx[i]);
166:       s = 1;
167:     }
168:     for (j=0;j<s;j++) d->nX[r_s+i+j] = 1.0;

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

173:     /* Check if the first eigenpairs are converged */
174:     if (i == 0) {
175:       d->preTestConv(d,0,s,s,&d->npreconv);
176:       if (d->npreconv > 0) break;
177:     }
178:   }

180:   /* D <- K*[Ax Bx] */
181:   if (d->npreconv == 0) {
182:     for (i=0;i<n;i++) {
183:       BVGetColumn(d->eps->V,kv+i,&v);
184:       d->improvex_precond(d,r_s+i,Ax[i],v);
185:       BVRestoreColumn(d->eps->V,kv+i,&v);
186:     }
187:     for (i=n;i<n*2;i++) {
188:       BVGetColumn(d->eps->V,kv+i,&v);
189:       d->improvex_precond(d,r_s+i-n,Bx[i-n],v);
190:       BVRestoreColumn(d->eps->V,kv+i,&v);
191:     }
192:     *size_D = 2*n;
193: #if !defined(PETSC_USE_COMPLEX)
194:     if (d->eigi[r_s] != 0.0) {
195:       s = 4;
196:     } else
197: #endif
198:     {
199:       s = 2;
200:     }
201:     /* Prevent that short vectors are discarded in the orthogonalization */
202:     for (i=0; i<s && i<*size_D; i++) {
203:       if (d->eps->errest[d->nconv+r_s+i] > PETSC_MACHINE_EPSILON && d->eps->errest[d->nconv+r_s+i] < PETSC_MAX_REAL) {
204:         BVScaleColumn(d->eps->V,i+kv,1.0/d->eps->errest[d->nconv+r_s+i]);
205:       }
206:     }
207:   } else *size_D = 0;

209:   SlepcVecPoolRestoreVecs(d->auxV,n,&Bx);
210:   SlepcVecPoolRestoreVecs(d->auxV,n,&Ax);
211:   BVDestroy(&X);
212:   MatDestroy(&M);
213:   return(0);
214: }

218: PetscErrorCode dvd_improvex_gd2(dvdDashboard *d,dvdBlackboard *b,KSP ksp,PetscInt max_bs)
219: {
220:   PetscErrorCode  ierr;
221:   dvdImprovex_gd2 *data;
222:   PC              pc;

225:   /* Setting configuration constrains */
226:   /* If the arithmetic is real and the problem is not Hermitian, then
227:      the block size is incremented in one */
228: #if !defined(PETSC_USE_COMPLEX)
229:   if (!DVD_IS(d->sEP, DVD_EP_HERMITIAN)) {
230:     max_bs++;
231:     b->max_size_P = PetscMax(b->max_size_P,2);
232:   } else
233: #endif
234:   {
235:     b->max_size_P = PetscMax(b->max_size_P,1);
236:   }
237:   b->max_size_X = PetscMax(b->max_size_X,max_bs);

239:   /* Setup the preconditioner */
240:   if (ksp) {
241:     KSPGetPC(ksp,&pc);
242:     dvd_static_precond_PC(d,b,pc);
243:   } else {
244:     dvd_static_precond_PC(d,b,0);
245:   }

247:   /* Setup the step */
248:   if (b->state >= DVD_STATE_CONF) {
249:     PetscNewLog(d->eps,&data);
250:     d->improveX_data = data;
251:     data->size_X = b->max_size_X;
252:     d->improveX = dvd_improvex_gd2_gen;

254:     EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_gd2_d);
255:   }
256:   return(0);
257: }