Actual source code: dvd_updatev.c

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

 29: static PetscErrorCode dvd_updateV_start(dvdDashboard *d);
 30: static PetscErrorCode dvd_isrestarting_fullV(dvdDashboard *d,PetscBool *r);
 31: static PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d);
 32: static PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d);
 33: static PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d);
 34: static PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d);
 35: static PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d);
 36: static PetscErrorCode dvd_updateV_testConv(dvdDashboard *d,PetscInt s,PetscInt pre,PetscInt e,PetscInt *nConv);

 38: typedef struct {
 39:   PetscInt
 40:     min_size_V,       /* restart with this number of eigenvectors */
 41:     plusk,            /* when restart, save plusk vectors from last iteration */
 42:     mpd;              /* max size of the searching subspace */
 43:   void *old_updateV_data;
 44:                       /* old updateV data */
 45:   isRestarting_type old_isRestarting;
 46:                       /* old isRestarting */
 47:   Mat oldU;           /* previous projected right igenvectors */
 48:   Mat oldV;           /* previous projected left eigenvectors */
 49:   PetscInt size_oldU; /* size of oldU */
 50:   PetscBool allResiduals;
 51:                       /* if computing all the residuals */
 52: } dvdManagV_basic;


 57: PetscErrorCode dvd_managementV_basic(dvdDashboard *d,dvdBlackboard *b,PetscInt bs,PetscInt mpd,PetscInt min_size_V,PetscInt plusk,PetscBool harm,PetscBool allResiduals)
 58: {
 59:   PetscErrorCode  ierr;
 60:   dvdManagV_basic *data;
 61: #if !defined(PETSC_USE_COMPLEX)
 62:   PetscBool       her_probl, std_probl;
 63: #endif

 66:   /* Setting configuration constrains */
 67: #if !defined(PETSC_USE_COMPLEX)
 68:   /* if the last converged eigenvalue is complex its conjugate pair is also
 69:      converged */
 70:   her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
 71:   std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
 72:   b->max_size_X = PetscMax(b->max_size_X, bs+((her_probl && std_probl)?0:1));
 73: #else
 74:   b->max_size_X = PetscMax(b->max_size_X, bs);
 75: #endif

 77:   b->max_size_V = PetscMax(b->max_size_V, mpd);
 78:   min_size_V = PetscMin(min_size_V, mpd-bs);
 79:   b->size_V = PetscMax(b->size_V, b->max_size_V + b->max_size_P + b->max_nev);
 80:   b->max_size_oldX = plusk;

 82:   /* Setup the step */
 83:   if (b->state >= DVD_STATE_CONF) {
 84:     PetscMalloc(sizeof(dvdManagV_basic),&data);
 85:     PetscLogObjectMemory((PetscObject)d->eps,sizeof(dvdManagV_basic));
 86:     data->mpd = b->max_size_V;
 87:     data->min_size_V = min_size_V;
 88:     d->bs = bs;
 89:     data->plusk = plusk;
 90:     data->allResiduals = allResiduals;

 92:     d->eigr = d->eps->eigr;
 93:     d->eigi = d->eps->eigi;
 94:     d->errest = d->eps->errest;
 95:     PetscMalloc1(d->eps->ncv,&d->real_nR);
 96:     PetscMalloc1(d->eps->ncv,&d->real_nX);
 97:     if (plusk > 0) {MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldU);}
 98:     else data->oldU = NULL;
 99:     if (harm && plusk>0) {MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldV);}
100:     else data->oldV = NULL;

102:     data->old_updateV_data = d->updateV_data;
103:     d->updateV_data = data;
104:     data->old_isRestarting = d->isRestarting;
105:     d->isRestarting = dvd_isrestarting_fullV;
106:     d->updateV = dvd_updateV_extrapol;
107:     d->preTestConv = dvd_updateV_testConv;
108:     EPSDavidsonFLAdd(&d->startList,dvd_updateV_start);
109:     EPSDavidsonFLAdd(&d->destroyList,dvd_managementV_basic_d);
110:   }
111:   return(0);
112: }

116: static PetscErrorCode dvd_updateV_start(dvdDashboard *d)
117: {
118:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
119:   PetscInt        i;

122:   for (i=0;i<d->eps->ncv;i++) d->eigi[i] = 0.0;
123:   d->nR = d->real_nR;
124:   for (i=0;i<d->eps->ncv;i++) d->nR[i] = PETSC_MAX_REAL;
125:   d->nX = d->real_nX;
126:   for (i=0;i<d->eps->ncv;i++) d->errest[i] = PETSC_MAX_REAL;
127:   data->size_oldU = 0;
128:   d->nconv = 0;
129:   d->npreconv = 0;
130:   d->V_tra_s = d->V_tra_e = d->V_new_s = d->V_new_e = 0;
131:   d->size_D = 0;
132:   return(0);
133: }

137: static PetscErrorCode dvd_isrestarting_fullV(dvdDashboard *d,PetscBool *r)
138: {
139:   PetscErrorCode  ierr;
140:   PetscInt        l,k;
141:   PetscBool       restart;
142:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;

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

148:   /* Check old isRestarting function */
149:   if (!restart && data->old_isRestarting) {
150:     data->old_isRestarting(d,&restart);
151:   }
152:   *r = restart;
153:   return(0);
154: }

158: static PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
159: {
160:   PetscErrorCode  ierr;
161:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;

164:   /* Restore changes in dvdDashboard */
165:   d->updateV_data = data->old_updateV_data;

167:   /* Free local data */
168:   if (data->oldU) {MatDestroy(&data->oldU);}
169:   if (data->oldV) {MatDestroy(&data->oldV);}
170:   PetscFree(d->real_nR);
171:   PetscFree(d->real_nX);
172:   PetscFree(data);
173:   return(0);
174: }

178: static PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
179: {
180:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
181:   PetscInt        i;
182:   PetscBool       restart;
183:   PetscErrorCode  ierr;

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

189:   /* If the subspaces doesn't need restart, add new vector */
190:   d->isRestarting(d,&restart);
191:   if (!restart) {
192:     d->size_D = 0;
193:     dvd_updateV_update_gen(d);

195:     /* If some vector were add, exit */
196:     if (d->size_D > 0) return(0);
197:   }

199:   /* If some eigenpairs were converged, lock them  */
200:   if (d->npreconv > 0) {
201:     i = d->npreconv;
202:     dvd_updateV_conv_gen(d);

204:     /* If some eigenpair was locked, exit */
205:     if (i > d->npreconv) return(0);
206:   }

208:   /* Else, a restarting is performed */
209:   dvd_updateV_restart_gen(d);
210:   return(0);
211: }

215: static PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
216: {
217:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
218:   PetscInt        npreconv,cMT,cMTX,lV,kV,nV;
219:   PetscErrorCode  ierr;
220:   Mat             Q,Z;
221: #if !defined(PETSC_USE_COMPLEX)
222:   PetscInt        i;
223: #endif

226:   npreconv = d->npreconv;
227:   /* Constrains the converged pairs to nev */
228: #if !defined(PETSC_USE_COMPLEX)
229:   /* Tries to maintain together conjugate eigenpairs */
230:   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));
231:   npreconv = i;
232: #else
233:   npreconv = PetscMax(PetscMin(d->nev - d->nconv, npreconv), 0);
234: #endif
235:   /* Quick exit */
236:   if (npreconv == 0) return(0);

238:   BVGetActiveColumns(d->eps->V,&lV,&kV);
239:   nV = kV - lV; 
240:   cMT = nV - npreconv;
241:   /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
242:      If the problem is standard or hermitian, left and right vectors are the same */
243:   if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
244:     /* ps.Q <- [ps.Q(0:npreconv-1) ps.Z(npreconv:size_H-1)] */
245:     DSGetMat(d->eps->ds,DS_MAT_Q,&Q);
246:     DSGetMat(d->eps->ds,DS_MAT_Z,&Z);
247:     SlepcMatDenseCopy(Z,0,npreconv,Q,0,npreconv,nV,cMT);
248:     DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q);
249:     DSRestoreMat(d->eps->ds,DS_MAT_Z,&Z);
250:   }
251:   if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
252:     DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,nV,d->nBds,&cMTX,d->nBds);
253:   } else {
254:     DSOrthogonalize(d->eps->ds,DS_MAT_Q,nV,&cMTX);
255:   }
256:   cMT = cMTX - npreconv;

258:   if (d->W) {
259:     DSOrthogonalize(d->eps->ds,DS_MAT_Z,nV,&cMTX);
260:     cMT = PetscMin(cMT,cMTX - npreconv);
261:   }

263:   /* Lock the converged pairs */
264:   d->eigr+= npreconv;
265: #if !defined(PETSC_USE_COMPLEX)
266:   if (d->eigi) d->eigi+= npreconv;
267: #endif
268:   d->nconv+= npreconv;
269:   d->errest+= npreconv;
270:   /* Notify the changes in V and update the other subspaces */
271:   d->V_tra_s = npreconv;          d->V_tra_e = nV;
272:   d->V_new_s = cMT;               d->V_new_e = d->V_new_s;
273:   /* Remove oldU */
274:   data->size_oldU = 0;

276:   d->npreconv-= npreconv;
277:   return(0);
278: }

282: static PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
283: {
284:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
285:   PetscInt        lV,kV,nV,size_plusk,size_X,cMTX,cMTY;
286:   Mat             Q,Z;
287:   PetscErrorCode  ierr;

290:   /* Select size_X desired pairs from V */
291:   BVGetActiveColumns(d->eps->V,&lV,&kV);
292:   nV = kV - lV;
293:   size_X = PetscMin(data->min_size_V,nV);

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

298:   d->size_MT = nV;
299:   /* ps.Q <- orth([pX(0:size_X-1) [oldU(0:size_plusk-1); 0] ]) */
300:   /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
301:      If the problem is standard or hermitian, left and right vectors are the same */
302:   DSGetMat(d->eps->ds,DS_MAT_Q,&Q);
303:   if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
304:     DSGetMat(d->eps->ds,DS_MAT_Z,&Z);
305:     SlepcMatDenseCopy(Z,0,0,Q,0,0,nV,size_X);
306:     DSRestoreMat(d->eps->ds,DS_MAT_Z,&Z);
307:   }
308:   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");
309:   if (size_plusk > 0) {
310:     SlepcMatDenseCopy(data->oldU,0,0,Q,0,size_X,nV,size_plusk);
311:   }
312:   DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q);
313:   if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
314:     DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,size_X,d->nBds,&cMTX,d->nBds);
315:   } else {
316:     DSOrthogonalize(d->eps->ds,DS_MAT_Q,size_X+size_plusk,&cMTX);
317:   }

319:   if (d->W && size_plusk > 0) {
320:     /* ps.Z <- orth([ps.Z(0:size_X-1) [oldV(0:size_plusk-1); 0] ]) */
321:     DSGetMat(d->eps->ds,DS_MAT_Z,&Z);
322:     SlepcMatDenseCopy(data->oldV,0,0,Z,0,size_X,nV,size_plusk);
323:     DSRestoreMat(d->eps->ds,DS_MAT_Z,&Z);
324:     DSOrthogonalize(d->eps->ds,DS_MAT_Z,size_X+size_plusk,&cMTY);
325:     cMTX = PetscMin(cMTX, cMTY);
326:   }

328:   /* Notify the changes in V and update the other subspaces */
329:   d->V_tra_s = 0;                     d->V_tra_e = cMTX;
330:   d->V_new_s = d->V_tra_e;            d->V_new_e = d->V_new_s;

332:   /* Remove oldU */
333:   data->size_oldU = 0;

335:   /* Remove npreconv */
336:   d->npreconv = 0;
337:   return(0);
338: }

342: static PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
343: {
344:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
345:   PetscInt        size_D,s,lV,kV,nV;
346:   Mat             Q,Z;
347:   PetscErrorCode  ierr;

350:   /* Select the desired pairs */
351:   BVGetActiveColumns(d->eps->V,&lV,&kV);
352:   nV = kV - lV;
353:   size_D = PetscMin(PetscMin(PetscMin(d->bs,nV),d->eps->ncv-nV),nV);
354:   if (size_D == 0) {
355:     d->initV(d);
356:     d->calcPairs(d);
357:   }

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

362:   /* If D is empty, exit */
363:   d->size_D = size_D;
364:   if (size_D == 0) return(0);

366:   /* Get the residual of all pairs */
367: #if !defined(PETSC_USE_COMPLEX)
368:   s = d->eigi[0]!=0.0?2:1;
369: #else
370:   s = 1;
371: #endif
372:   BVGetActiveColumns(d->eps->V,&lV,&kV);
373:   nV = kV - lV;
374:   dvd_updateV_testConv(d,s,s,data->allResiduals?nV:size_D,NULL);

376:   /* Notify the changes in V */
377:   d->V_tra_s = 0;                 d->V_tra_e = 0;
378:   d->V_new_s = nV;                d->V_new_e = nV+size_D;

380:   /* Save the projected eigenvectors */
381:   if (data->plusk > 0) {
382:     MatZeroEntries(data->oldU);
383:     data->size_oldU = nV;
384:     DSGetMat(d->eps->ds,DS_MAT_Q,&Q);
385:     SlepcMatDenseCopy(Q,0,0,data->oldU,0,0,nV,nV);
386:     DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q);
387:     if (d->W) {
388:       MatZeroEntries(data->oldV);
389:       DSGetMat(d->eps->ds,DS_MAT_Z,&Z);
390:       SlepcMatDenseCopy(Z,0,0,data->oldV,0,0,nV,nV);
391:       DSRestoreMat(d->eps->ds,DS_MAT_Z,&Z);
392:     }
393:   }
394:   return(0);
395: }

399: static PetscErrorCode dvd_updateV_testConv(dvdDashboard *d,PetscInt s,PetscInt pre,PetscInt e,PetscInt *nConv)
400: {
401:   PetscInt        i,j,b;
402:   PetscReal       norm;
403:   PetscErrorCode  ierr;
404:   PetscBool       conv, c;
405:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;

408:   if (nConv) *nConv = s;
409:   for (i=s, conv=PETSC_TRUE;
410:       (conv || data->allResiduals) && (i < e);
411:       i+=b) {
412: #if !defined(PETSC_USE_COMPLEX)
413:     b = d->eigi[i]!=0.0?2:1;
414: #else
415:     b = 1;
416: #endif
417:     if (i+b-1 >= pre) {
418:       d->calcpairs_residual(d,i,i+b);
419:     }
420:     /* Test the Schur vector */
421:     for (j=0,c=PETSC_TRUE; j<b && c; j++) {
422:       norm = d->nR[i+j]/d->nX[i+j];
423:       c = d->testConv(d,d->eigr[i+j],d->eigi[i+j],norm,&d->errest[i+j]);
424:     }
425:     if (conv && c) { if (nConv) *nConv = i+b; }
426:     else conv = PETSC_FALSE;
427:   }
428:   pre = PetscMax(pre, i);

430: #if !defined(PETSC_USE_COMPLEX)
431:   /* Enforce converged conjugate complex eigenpairs */
432:   if (nConv) {
433:     for (j=0;j<*nConv;j++) if (d->eigi[j] != 0.0) j++;
434:     if (j>*nConv) (*nConv)--;
435:   }
436: #endif
437:   for (i=pre;i<e;i++) d->errest[i] = d->nR[i] = PETSC_MAX_REAL;
438:   return(0);
439: }