Actual source code: dvdupdatev.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:   SLEPc eigensolver: "davidson"

  4:   Step: test for restarting, updateV, restartV

  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
 27: #include <slepc/private/dsimpl.h>

 29: typedef struct {
 30:   PetscInt          min_size_V;        /* restart with this number of eigenvectors */
 31:   PetscInt          plusk;             /* at restart, save plusk vectors from last iteration */
 32:   PetscInt          mpd;               /* max size of the searching subspace */
 33:   void              *old_updateV_data; /* old updateV data */
 34:   isRestarting_type old_isRestarting;  /* old isRestarting */
 35:   Mat               oldU;              /* previous projected right igenvectors */
 36:   Mat               oldV;              /* previous projected left eigenvectors */
 37:   PetscInt          size_oldU;         /* size of oldU */
 38:   PetscBool         allResiduals;      /* if computing all the residuals */
 39: } dvdManagV_basic;

 43: static PetscErrorCode dvd_updateV_start(dvdDashboard *d)
 44: {
 45:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
 46:   PetscInt        i;

 49:   for (i=0;i<d->eps->ncv;i++) d->eigi[i] = 0.0;
 50:   d->nR = d->real_nR;
 51:   for (i=0;i<d->eps->ncv;i++) d->nR[i] = PETSC_MAX_REAL;
 52:   d->nX = d->real_nX;
 53:   for (i=0;i<d->eps->ncv;i++) d->errest[i] = PETSC_MAX_REAL;
 54:   data->size_oldU = 0;
 55:   d->nconv = 0;
 56:   d->npreconv = 0;
 57:   d->V_tra_s = d->V_tra_e = d->V_new_s = d->V_new_e = 0;
 58:   d->size_D = 0;
 59:   return(0);
 60: }

 64: static PetscErrorCode dvd_isrestarting_fullV(dvdDashboard *d,PetscBool *r)
 65: {
 66:   PetscErrorCode  ierr;
 67:   PetscInt        l,k;
 68:   PetscBool       restart;
 69:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;

 72:   BVGetActiveColumns(d->eps->V,&l,&k);
 73:   restart = (k+2 > d->eps->ncv)? PETSC_TRUE: PETSC_FALSE;

 75:   /* Check old isRestarting function */
 76:   if (!restart && data->old_isRestarting) {
 77:     data->old_isRestarting(d,&restart);
 78:   }
 79:   *r = restart;
 80:   return(0);
 81: }

 85: static PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
 86: {
 87:   PetscErrorCode  ierr;
 88:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;

 91:   /* Restore changes in dvdDashboard */
 92:   d->updateV_data = data->old_updateV_data;

 94:   /* Free local data */
 95:   if (data->oldU) { MatDestroy(&data->oldU); }
 96:   if (data->oldV) { MatDestroy(&data->oldV); }
 97:   PetscFree(d->real_nR);
 98:   PetscFree(d->real_nX);
 99:   PetscFree(data);
100:   return(0);
101: }

105: static PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
106: {
107:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
108:   PetscInt        npreconv,cMT,cMTX,lV,kV,nV;
109:   PetscErrorCode  ierr;
110:   Mat             Q;
111: #if !defined(PETSC_USE_COMPLEX)
112:   PetscInt        i;
113: #endif

116:   npreconv = d->npreconv;
117:   /* Constrains the converged pairs to nev */
118: #if !defined(PETSC_USE_COMPLEX)
119:   /* Tries to maintain together conjugate eigenpairs */
120:   for (i=0; (i + (d->eigi[i]!=0.0?1:0) < npreconv) && (d->nconv + i < d->nev); i+= (d->eigi[i]!=0.0?2:1));
121:   npreconv = i;
122: #else
123:   npreconv = PetscMax(PetscMin(d->nev-d->nconv,npreconv),0);
124: #endif
125:   /* Quick exit */
126:   if (npreconv == 0) return(0);

128:   BVGetActiveColumns(d->eps->V,&lV,&kV);
129:   nV  = kV - lV; 
130:   cMT = nV - npreconv;
131:   /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
132:      If the problem is standard or hermitian, left and right vectors are the same */
133:   if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
134:     /* ps.Q <- [ps.Q(0:npreconv-1) ps.Z(npreconv:size_H-1)] */
135:     DSGetMat(d->eps->ds,DS_MAT_Q,&Q);
136:     DSCopyMat(d->eps->ds,DS_MAT_Z,0,npreconv,Q,0,npreconv,nV,cMT,PETSC_TRUE);
137:     DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q);
138:   }
139:   if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
140:     DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,nV,d->nBds,&cMTX,d->nBds);
141:   } else {
142:     DSOrthogonalize(d->eps->ds,DS_MAT_Q,nV,&cMTX);
143:   }
144:   cMT = cMTX - npreconv;

146:   if (d->W) {
147:     DSOrthogonalize(d->eps->ds,DS_MAT_Z,nV,&cMTX);
148:     cMT = PetscMin(cMT,cMTX - npreconv);
149:   }

151:   /* Lock the converged pairs */
152:   d->eigr+= npreconv;
153: #if !defined(PETSC_USE_COMPLEX)
154:   if (d->eigi) d->eigi+= npreconv;
155: #endif
156:   d->nconv+= npreconv;
157:   d->errest+= npreconv;
158:   /* Notify the changes in V and update the other subspaces */
159:   d->V_tra_s = npreconv;          d->V_tra_e = nV;
160:   d->V_new_s = cMT;               d->V_new_e = d->V_new_s;
161:   /* Remove oldU */
162:   data->size_oldU = 0;

164:   d->npreconv-= npreconv;
165:   return(0);
166: }

170: static PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
171: {
172:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
173:   PetscInt        lV,kV,nV,size_plusk,size_X,cMTX,cMTY;
174:   Mat             Q;
175:   PetscErrorCode  ierr;

178:   /* Select size_X desired pairs from V */
179:   BVGetActiveColumns(d->eps->V,&lV,&kV);
180:   nV = kV - lV;
181:   size_X = PetscMin(data->min_size_V,nV);

183:   /* Add plusk eigenvectors from the previous iteration */
184:   size_plusk = PetscMax(0,PetscMin(PetscMin(data->plusk,data->size_oldU),d->eps->ncv - size_X));

186:   d->size_MT = nV;
187:   /* ps.Q <- orth([pX(0:size_X-1) [oldU(0:size_plusk-1); 0] ]) */
188:   /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
189:      If the problem is standard or hermitian, left and right vectors are the same */
190:   if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
191:     DSGetMat(d->eps->ds,DS_MAT_Q,&Q);
192:     DSCopyMat(d->eps->ds,DS_MAT_Z,0,0,Q,0,0,nV,size_X,PETSC_TRUE);
193:     DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q);
194:   }
195:   if (size_plusk > 0 && DVD_IS(d->sEP,DVD_EP_INDEFINITE)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported plusk>0 in indefinite eigenvalue problems");
196:   if (size_plusk > 0) {
197:     DSCopyMat(d->eps->ds,DS_MAT_Q,0,size_X,data->oldU,0,0,nV,size_plusk,PETSC_FALSE);
198:   }
199:   if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
200:     DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,size_X,d->nBds,&cMTX,d->nBds);
201:   } else {
202:     DSOrthogonalize(d->eps->ds,DS_MAT_Q,size_X+size_plusk,&cMTX);
203:   }

205:   if (d->W && size_plusk > 0) {
206:     /* ps.Z <- orth([ps.Z(0:size_X-1) [oldV(0:size_plusk-1); 0] ]) */
207:     DSCopyMat(d->eps->ds,DS_MAT_Z,0,size_X,data->oldV,0,0,nV,size_plusk,PETSC_FALSE);
208:     DSOrthogonalize(d->eps->ds,DS_MAT_Z,size_X+size_plusk,&cMTY);
209:     cMTX = PetscMin(cMTX, cMTY);
210:   }

212:   /* Notify the changes in V and update the other subspaces */
213:   d->V_tra_s = 0;                     d->V_tra_e = cMTX;
214:   d->V_new_s = d->V_tra_e;            d->V_new_e = d->V_new_s;

216:   /* Remove oldU */
217:   data->size_oldU = 0;

219:   /* Remove npreconv */
220:   d->npreconv = 0;
221:   return(0);
222: }

226: static PetscErrorCode dvd_updateV_testConv(dvdDashboard *d,PetscInt s,PetscInt pre,PetscInt e,PetscInt *nConv)
227: {
228:   PetscInt        i,j,b;
229:   PetscReal       norm;
230:   PetscErrorCode  ierr;
231:   PetscBool       conv, c;
232:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;

235:   if (nConv) *nConv = s;
236:   for (i=s,conv=PETSC_TRUE;(conv || data->allResiduals) && (i < e);i+=b) {
237: #if !defined(PETSC_USE_COMPLEX)
238:     b = d->eigi[i]!=0.0?2:1;
239: #else
240:     b = 1;
241: #endif
242:     if (i+b-1 >= pre) {
243:       d->calcpairs_residual(d,i,i+b);
244:     }
245:     /* Test the Schur vector */
246:     for (j=0,c=PETSC_TRUE;j<b && c;j++) {
247:       norm = d->nR[i+j]/d->nX[i+j];
248:       c = d->testConv(d,d->eigr[i+j],d->eigi[i+j],norm,&d->errest[i+j]);
249:     }
250:     if (conv && c) { if (nConv) *nConv = i+b; }
251:     else conv = PETSC_FALSE;
252:   }
253:   pre = PetscMax(pre,i);

255: #if !defined(PETSC_USE_COMPLEX)
256:   /* Enforce converged conjugate complex eigenpairs */
257:   if (nConv) {
258:     for (j=0;j<*nConv;j++) if (d->eigi[j] != 0.0) j++;
259:     if (j>*nConv) (*nConv)--;
260:   }
261: #endif
262:   for (i=pre;i<e;i++) d->errest[i] = d->nR[i] = PETSC_MAX_REAL;
263:   return(0);
264: }

268: static PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
269: {
270:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
271:   PetscInt        size_D,s,lV,kV,nV;
272:   PetscErrorCode  ierr;

275:   /* Select the desired pairs */
276:   BVGetActiveColumns(d->eps->V,&lV,&kV);
277:   nV = kV - lV;
278:   size_D = PetscMin(PetscMin(PetscMin(d->bs,nV),d->eps->ncv-nV),nV);
279:   if (size_D == 0) {
280:     d->initV(d);
281:     d->calcPairs(d);
282:   }

284:   /* Fill V with D */
285:   d->improveX(d,0,size_D,&size_D);

287:   /* If D is empty, exit */
288:   d->size_D = size_D;
289:   if (size_D == 0) return(0);

291:   /* Get the residual of all pairs */
292: #if !defined(PETSC_USE_COMPLEX)
293:   s = (d->eigi[0]!=0.0)? 2: 1;
294: #else
295:   s = 1;
296: #endif
297:   BVGetActiveColumns(d->eps->V,&lV,&kV);
298:   nV = kV - lV;
299:   dvd_updateV_testConv(d,s,s,data->allResiduals?nV:size_D,NULL);

301:   /* Notify the changes in V */
302:   d->V_tra_s = 0;                 d->V_tra_e = 0;
303:   d->V_new_s = nV;                d->V_new_e = nV+size_D;

305:   /* Save the projected eigenvectors */
306:   if (data->plusk > 0) {
307:     MatZeroEntries(data->oldU);
308:     data->size_oldU = nV;
309:     DSCopyMat(d->eps->ds,DS_MAT_Q,0,0,data->oldU,0,0,nV,nV,PETSC_TRUE);
310:     if (d->W) {
311:       MatZeroEntries(data->oldV);
312:       DSCopyMat(d->eps->ds,DS_MAT_Z,0,0,data->oldV,0,0,nV,nV,PETSC_TRUE);
313:     }
314:   }
315:   return(0);
316: }

320: static PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
321: {
322:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
323:   PetscInt        i;
324:   PetscBool       restart;
325:   PetscErrorCode  ierr;

328:   /* TODO: restrict select pairs to each case */
329:   d->calcpairs_selectPairs(d, data->min_size_V);

331:   /* If the subspaces doesn't need restart, add new vector */
332:   d->isRestarting(d,&restart);
333:   if (!restart) {
334:     d->size_D = 0;
335:     dvd_updateV_update_gen(d);

337:     /* If some vector were add, exit */
338:     if (d->size_D > 0) return(0);
339:   }

341:   /* If some eigenpairs were converged, lock them  */
342:   if (d->npreconv > 0) {
343:     i = d->npreconv;
344:     dvd_updateV_conv_gen(d);

346:     /* If some eigenpair was locked, exit */
347:     if (i > d->npreconv) return(0);
348:   }

350:   /* Else, a restarting is performed */
351:   dvd_updateV_restart_gen(d);
352:   return(0);
353: }

357: PetscErrorCode dvd_managementV_basic(dvdDashboard *d,dvdBlackboard *b,PetscInt bs,PetscInt mpd,PetscInt min_size_V,PetscInt plusk,PetscBool harm,PetscBool allResiduals)
358: {
359:   PetscErrorCode  ierr;
360:   dvdManagV_basic *data;
361: #if !defined(PETSC_USE_COMPLEX)
362:   PetscBool       her_probl,std_probl;
363: #endif

366:   /* Setting configuration constrains */
367: #if !defined(PETSC_USE_COMPLEX)
368:   /* if the last converged eigenvalue is complex its conjugate pair is also
369:      converged */
370:   her_probl = DVD_IS(d->sEP,DVD_EP_HERMITIAN)? PETSC_TRUE: PETSC_FALSE;
371:   std_probl = DVD_IS(d->sEP,DVD_EP_STD)? PETSC_TRUE: PETSC_FALSE;
372:   b->max_size_X = PetscMax(b->max_size_X,bs+((her_probl && std_probl)?0:1));
373: #else
374:   b->max_size_X = PetscMax(b->max_size_X,bs);
375: #endif

377:   b->max_size_V = PetscMax(b->max_size_V,mpd);
378:   min_size_V = PetscMin(min_size_V,mpd-bs);
379:   b->size_V = PetscMax(b->size_V,b->max_size_V+b->max_size_P+b->max_nev);
380:   b->max_size_oldX = plusk;

382:   /* Setup the step */
383:   if (b->state >= DVD_STATE_CONF) {
384:     PetscNewLog(d->eps,&data);
385:     data->mpd = b->max_size_V;
386:     data->min_size_V = min_size_V;
387:     d->bs = bs;
388:     data->plusk = plusk;
389:     data->allResiduals = allResiduals;

391:     d->eigr = d->eps->eigr;
392:     d->eigi = d->eps->eigi;
393:     d->errest = d->eps->errest;
394:     PetscMalloc1(d->eps->ncv,&d->real_nR);
395:     PetscMalloc1(d->eps->ncv,&d->real_nX);
396:     if (plusk > 0) { MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldU); }
397:     else data->oldU = NULL;
398:     if (harm && plusk>0) { MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldV); }
399:     else data->oldV = NULL;

401:     data->old_updateV_data = d->updateV_data;
402:     d->updateV_data = data;
403:     data->old_isRestarting = d->isRestarting;
404:     d->isRestarting = dvd_isrestarting_fullV;
405:     d->updateV = dvd_updateV_extrapol;
406:     d->preTestConv = dvd_updateV_testConv;
407:     EPSDavidsonFLAdd(&d->startList,dvd_updateV_start);
408:     EPSDavidsonFLAdd(&d->destroyList,dvd_managementV_basic_d);
409:   }
410:   return(0);
411: }