Actual source code: dvd_updatev.c

  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-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 <slepc-private/dsimpl.h>      /*I "slepcds.h" I*/

 29: PetscErrorCode dvd_updateV_start(dvdDashboard *d);
 30: PetscBool dvd_isrestarting_fullV(dvdDashboard *d);
 31: PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d);
 32: PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d);
 33: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d);
 34: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d);
 35: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d);
 36: PetscErrorCode dvd_updateV_conv_finish(dvdDashboard *d);
 37: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d, PetscInt s, PetscInt pre,
 38:                                     PetscInt e, Vec *auxV, PetscScalar *auxS,
 39:                                     PetscInt *nConv);

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


 65: PetscErrorCode dvd_managementV_basic(dvdDashboard *d, dvdBlackboard *b,
 66:                                      PetscInt bs, PetscInt mpd,
 67:                                      PetscInt min_size_V,
 68:                                      PetscInt plusk, PetscBool harm,
 69:                                      PetscBool allResiduals)
 70: {
 71:   PetscErrorCode  ierr;
 72:   dvdManagV_basic *data;
 73: #if !defined(PETSC_USE_COMPLEX)
 74:   PetscBool       her_probl, std_probl;
 75: #endif

 78:   /* Setting configuration constrains */
 79: #if !defined(PETSC_USE_COMPLEX)
 80:   /* if the last converged eigenvalue is complex its conjugate pair is also
 81:      converged */
 82:   her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
 83:   std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
 84:   b->max_size_X = PetscMax(b->max_size_X, bs+(her_probl && std_probl)?0:1);
 85: #else
 86:   b->max_size_X = PetscMax(b->max_size_X, bs);
 87: #endif

 89:   b->max_size_V = PetscMax(b->max_size_V, mpd);
 90:   min_size_V = PetscMin(min_size_V, mpd-bs);
 91:   b->max_size_auxV = PetscMax(b->max_size_auxV, 1); /* dvd_updateV_testConv */
 92:   b->size_V = PetscMax(b->size_V, b->max_size_V + b->max_size_P + b->max_nev);
 93:   b->own_scalars+= b->size_V*2 /* eigr, eigr */ +
 94:                    b->size_V /* nR */   +
 95:                    b->size_V /* nX */   +
 96:                    b->size_V /* errest */ +
 97:                    b->max_size_V*b->max_size_V*(harm?2:1)*(plusk>0?1:0)
 98:                                                /* oldU,oldV? */;
 99:   b->max_size_oldX = plusk;

101:   /* Setup the step */
102:   if (b->state >= DVD_STATE_CONF) {
103:     PetscMalloc(sizeof(dvdManagV_basic), &data);
104:     data->mpd = b->max_size_V;
105:     data->min_size_V = min_size_V;
106:     d->bs = bs;
107:     d->max_size_X = b->max_size_X;
108:     data->plusk = plusk;
109:     data->allResiduals = allResiduals;

111:     d->size_real_eigr = b->size_V;
112:     d->real_eigr = b->free_scalars; b->free_scalars+= b->size_V;
113:     d->real_eigi = b->free_scalars; b->free_scalars+= b->size_V;
114:     d->real_nR = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(b->size_V);
115:     d->real_nX = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(b->size_V);
116:     d->real_errest = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(b->size_V);
117:     if (plusk > 0) {
118:       data->oldU = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
119:     }
120:     if (harm) {
121:       if (plusk > 0) {
122:         data->oldV = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
123:       }
124:     } else {
125:       data->oldV = PETSC_NULL;
126:     }

128:     data->old_updateV_data = d->updateV_data;
129:     d->updateV_data = data;
130:     data->old_isRestarting = d->isRestarting;
131:     d->isRestarting = dvd_isrestarting_fullV;
132:     d->updateV = dvd_updateV_extrapol;
133:     d->preTestConv = dvd_updateV_testConv;
134:     DVD_FL_ADD(d->startList, dvd_updateV_start);
135:     DVD_FL_ADD(d->endList, dvd_updateV_conv_finish);
136:     DVD_FL_ADD(d->destroyList, dvd_managementV_basic_d);
137:   }
138:   return(0);
139: }

143: PetscErrorCode dvd_updateV_start(dvdDashboard *d)
144: {
145:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
146:   PetscInt        i;


150:   d->size_cX = 0;
151:   d->eigr = d->ceigr = d->real_eigr;
152:   d->eigi = d->ceigi = d->real_eigi;
153: #if defined(PETSC_USE_COMPLEX)
154:   for(i=0; i<d->size_real_V; i++) d->eigi[i] = 0.0;
155: #endif
156:   d->nR = d->real_nR;
157:   for(i=0; i<d->size_real_V; i++) d->nR[i] = PETSC_MAX_REAL;
158:   d->nX = d->real_nX;
159:   d->errest = d->real_errest;
160:   for(i=0; i<d->size_real_V; i++) d->errest[i] = PETSC_MAX_REAL;
161:   data->ldoldU = 0;
162:   data->oldV = PETSC_NULL;
163:   data->size_oldU = 0;
164:   d->nconv = 0;
165:   d->npreconv = 0;
166:   d->V_tra_s = d->V_tra_e = d->V_new_s = d->V_new_e = 0;
167:   d->size_D = 0;

169:   return(0);
170: }


175: PetscBool dvd_isrestarting_fullV(dvdDashboard *d)
176: {
177:   PetscBool       restart;
178:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;


182:   restart = (d->size_V + d->max_size_X > PetscMin(data->mpd,d->max_size_V))?
183:                 PETSC_TRUE:PETSC_FALSE;

185:   /* Check old isRestarting function */
186:   if (!restart && data->old_isRestarting)
187:     restart = data->old_isRestarting(d);

189:   PetscFunctionReturn(restart);
190: }

194: PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
195: {
196:   PetscErrorCode  ierr;
197:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;


201:   /* Restore changes in dvdDashboard */
202:   d->updateV_data = data->old_updateV_data;
203: 
204:   /* Free local data */
205:   PetscFree(data);

207:   return(0);
208: }

212: PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
213: {
214:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
215:   PetscInt        i;
216:   PetscErrorCode  ierr;


220:   d->calcpairs_selectPairs(d, data->min_size_V);

222:   /* If the subspaces doesn't need restart, add new vector */
223:   if (!d->isRestarting(d)) {
224:     d->size_D = 0;
225:     dvd_updateV_update_gen(d);

227:     /* If some vector were add, exit */
228:     if (d->size_D > 0) { return(0); }
229:   }

231:   /* If some eigenpairs were converged, lock them  */
232:   if (d->npreconv > 0) {
233:     i = d->npreconv;
234:     dvd_updateV_conv_gen(d);

236:     /* If some eigenpair was locked, exit */
237:     if (i > d->npreconv) { return(0); }
238:   }

240:   /* Else, a restarting is performed */
241:   dvd_updateV_restart_gen(d);

243:   return(0);
244: }

248: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
249: {
250:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
251:   PetscInt        npreconv,ld,cMT,cMTX;
252:   PetscErrorCode  ierr;
253:   PetscScalar     *pQ,*pZ;
254: #if !defined(PETSC_USE_COMPLEX)
255:   PetscInt        i;
256: #endif


260:   npreconv = d->npreconv;
261:   /* Constrains the converged pairs to nev */
262: #if !defined(PETSC_USE_COMPLEX)
263:   /* Tries to maintain together conjugate eigenpairs */
264:   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));
265:   npreconv = i;
266: #else
267:   npreconv = PetscMax(PetscMin(d->nev - d->nconv, npreconv), 0);
268: #endif
269:   /* Quick exit */
270:   if (npreconv == 0) { return(0); }

272:   npreconv+= d->cX_in_H;
273:   DSGetLeadingDimension(d->ps,&ld);
274:   d->size_MT = d->size_H;
275:   cMT = d->size_H - npreconv;
276:   /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
277:      If the problem is standard or hermitian, left and right vectors are the same */
278:   if (!(d->W||!d->cY||d->BcX||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
279:     /* ps.Q <- [ps.Q(0:npreconv-1) ps.Z(npreconv:size_H-1)] */
280:     DSGetArray(d->ps,DS_MAT_Q,&pQ);
281:     DSGetArray(d->ps,DS_MAT_Z,&pZ);
282:     SlepcDenseCopy(&pQ[ld*npreconv],ld,&pZ[ld*npreconv],ld,d->size_H,cMT);
283:     DSRestoreArray(d->ps,DS_MAT_Q,&pQ);
284:     DSRestoreArray(d->ps,DS_MAT_Z,&pZ);
285:   }
286:   if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
287:     DSPseudoOrthogonalize(d->ps,DS_MAT_Q,d->size_H,d->nBV-d->cX_in_H,&cMTX,d->nBpX);
288:   } else {
289:     DSOrthogonalize(d->ps,DS_MAT_Q,d->size_H,&cMTX);
290:   }
291:   cMT = cMTX - npreconv;

293:   if (d->W) {
294:     DSOrthogonalize(d->ps,DS_MAT_Z,d->size_H,&cMTX);
295:     cMT = PetscMin(cMT,cMTX - npreconv);
296:   }

298:   /* Lock the converged pairs */
299:   d->eigr+= npreconv-d->cX_in_H;
300: #if !defined(PETSC_USE_COMPLEX)
301:   if (d->eigi) d->eigi+= npreconv-d->cX_in_H;
302: #endif
303:   d->nconv+= npreconv-d->cX_in_H;
304:   d->errest+= npreconv-d->cX_in_H;
305:   /* Notify the changes in V and update the other subspaces */
306:   d->V_tra_s = npreconv;          d->V_tra_e = d->size_H;
307:   d->V_new_s = cMT;               d->V_new_e = d->V_new_s;
308:   /* Remove oldU */
309:   data->size_oldU = 0;

311:   d->npreconv-= npreconv-d->cX_in_H;
312:   return(0);
313: }

317: PetscErrorCode dvd_updateV_conv_finish(dvdDashboard *d)
318: {
319:   PetscErrorCode  ierr;
320:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
321: #if defined(PETSC_USE_COMPLEX)
322:   PetscInt        i, j;
323:   PetscScalar     s;
324: #endif  


328:   /* Some functions need the diagonal elements in cT be real */
329: #if defined(PETSC_USE_COMPLEX)
330:   if (d->cT) for(i=0; i<d->nconv; i++) {
331:     s = PetscConj(d->cT[d->ldcT*i+i])/PetscAbsScalar(d->cT[d->ldcT*i+i]);
332:     for(j=0; j<=i; j++)
333:       d->cT[d->ldcT*i+j] = PetscRealPart(d->cT[d->ldcT*i+j]*s),
334:       d->cS[d->ldcS*i+j]*= s;
335:     VecScale(d->cX[i], s);
336:   }
337: #endif
338:   d->calcpairs_selectPairs(d, data->min_size_V);

340:   return(0);
341: }
342: 
345: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
346: {
347:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
348:   PetscInt        size_plusk,size_X,i,j,ld,cMTX,cMTY;
349:   PetscScalar     *pQ,*pZ;
350:   PetscErrorCode  ierr;


354:   /* Select size_X desired pairs from V */
355:   size_X = PetscMin(PetscMin(data->min_size_V,
356:                              d->size_V ),
357:                              d->max_size_V );

359:   /* Add plusk eigenvectors from the previous iteration */
360:   size_plusk = PetscMax(0, PetscMin(PetscMin(data->plusk,
361:                                     data->size_oldU ),
362:                                     d->max_size_V - size_X ));

364:   DSGetLeadingDimension(d->ps,&ld);
365:   d->size_MT = d->size_H;
366:   /* ps.Q <- orth([pX(0:size_X-1) [oldU(0:size_plusk-1); 0] ]) */
367:   /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
368:      If the problem is standard or hermitian, left and right vectors are the same */
369:   DSGetArray(d->ps,DS_MAT_Q,&pQ);
370:   if (!(d->W||!d->cY||d->BcX||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
371:     DSGetArray(d->ps,DS_MAT_Z,&pZ);
372:     SlepcDenseCopy(pQ,ld,pZ,ld,d->size_H,size_X);
373:     DSRestoreArray(d->ps,DS_MAT_Z,&pZ);
374:   }
375:   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");
376:   if (size_plusk > 0) {
377:     SlepcDenseCopy(&pQ[ld*size_X],ld,data->oldU,data->ldoldU,data->size_oldU,size_plusk);
378:     for(i=size_X; i<size_X+size_plusk; i++) {
379:       for(j=data->size_oldU; j<d->size_H; j++) {
380:         pQ[j*ld+i] = 0.0;
381:       }
382:     }
383:   }
384:   DSRestoreArray(d->ps,DS_MAT_Q,&pQ);
385:   if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
386:     DSPseudoOrthogonalize(d->ps,DS_MAT_Q,size_X,d->nBV-d->cX_in_H,&cMTX,d->nBpX);
387:   } else {
388:     DSOrthogonalize(d->ps,DS_MAT_Q,size_X+size_plusk,&cMTX);
389:   }

391:   if (d->W) {
392:     /* ps.Z <- orth([ps.Z(0:size_X-1) [oldV(0:size_plusk-1); 0] ]) */
393:     if (size_plusk > 0) {
394:       DSGetArray(d->ps,DS_MAT_Z,&pZ);
395:       SlepcDenseCopy(&pZ[ld*size_X],ld,data->oldV,data->ldoldU,data->size_oldU,size_plusk);
396:       for(i=size_X; i<size_X+size_plusk; i++) {
397:         for(j=data->size_oldU; j<d->size_H; j++) {
398:           pZ[j*ld+i] = 0.0;
399:         }
400:       }
401:     }
402:     DSRestoreArray(d->ps,DS_MAT_Z,&pZ);
403:     DSOrthogonalize(d->ps,DS_MAT_Z,size_X+size_plusk,&cMTY);
404:     cMTX = PetscMin(cMTX, cMTY);
405:   }

407:   /* Notify the changes in V and update the other subspaces */
408:   d->V_tra_s = d->cX_in_H;                  d->V_tra_e = cMTX;
409:   d->V_new_s = d->V_tra_e-d->cX_in_H; d->V_new_e = d->V_new_s;

411:   /* Remove oldU */
412:   data->size_oldU = 0;

414:   /* Remove npreconv */
415:   d->npreconv = 0;
416: 
417:   return(0);
418: }


423: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
424: {
425:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
426:   PetscInt        size_D,ld,s;
427:   PetscScalar     *pQ,*pZ;
428:   PetscErrorCode  ierr;


432:   /* Select the desired pairs */
433:   size_D = PetscMin(PetscMin(PetscMin(d->bs,
434:                                       d->size_V ),
435:                                       d->max_size_V-d->size_V ),
436:                                       d->size_H );
437:   if (size_D == 0) {
438:     PetscInfo2(d->eps, "MON: D:%D H:%D\n", size_D, d->size_H);
439:     d->initV(d);
440:     d->calcPairs(d);
441:   }

443:   /* Fill V with D */
444:   d->improveX(d, d->V+d->size_V, d->max_size_V-d->size_V, 0, size_D, &size_D);

446:   /* If D is empty, exit */
447:   d->size_D = size_D;
448:   if (size_D == 0) { return(0); }

450:   /* Get the residual of all pairs */
451: #if !defined(PETSC_USE_COMPLEX)
452:   s = d->eigi[0]!=0.0?2:1;
453: #else
454:   s = 1;
455: #endif
456:   dvd_updateV_testConv(d,s,s,data->allResiduals?d->size_V:size_D,d->auxV,d->auxS,PETSC_NULL);

458:   /* Notify the changes in V */
459:   d->V_tra_s = 0;                 d->V_tra_e = 0;
460:   d->V_new_s = d->size_V;         d->V_new_e = d->size_V+size_D;

462:   /* Save the projected eigenvectors */
463:   if (data->plusk > 0) {
464:     data->ldoldU = data->size_oldU = d->size_H;
465:     DSGetLeadingDimension(d->ps,&ld);
466:     DSGetArray(d->ps,DS_MAT_Q,&pQ);
467:     SlepcDenseCopy(data->oldU,data->ldoldU,pQ,ld,d->size_H,d->size_H);
468:     DSRestoreArray(d->ps,DS_MAT_Q,&pQ);
469:     if (d->cY) {
470:       DSGetArray(d->ps,DS_MAT_Z,&pZ);
471:       SlepcDenseCopy(data->oldV,data->ldoldU,pZ,ld,d->size_H,d->size_H);
472:       DSRestoreArray(d->ps,DS_MAT_Z,&pZ);
473:     }
474:   }

476:   return(0);
477: }


482: /* auxV: (by calcpairs_residual_eig) */
483: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d, PetscInt s, PetscInt pre,
484:                                     PetscInt e, Vec *auxV, PetscScalar *auxS,
485:                                     PetscInt *nConv)
486: {
487:   PetscInt        i,j,b;
488:   PetscReal       norm;
489:   PetscErrorCode  ierr;
490:   PetscBool       conv, c;
491:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;

494: 
495:   if (nConv) *nConv = s;
496:   for(i=s, conv=PETSC_TRUE;
497:       (conv || data->allResiduals) && (i < e);
498:       i+=b) {
499: #if !defined(PETSC_USE_COMPLEX)
500:     b = d->eigi[i]!=0.0?2:1;
501: #else
502:     b = 1;
503: #endif
504:     if (i+b-1 >= pre) {
505:       d->calcpairs_residual(d, i, i+b, auxV);
506: 
507:     }
508:     /* Test the Schur vector */
509:     for (j=0,c=PETSC_TRUE; j<b && c; j++) {
510:       norm = d->nR[i+j]/d->nX[i+j];
511:       c = d->testConv(d, d->eigr[i+j], d->eigi[i+j], norm, &d->errest[i+j]);
512:     }
513:     /* Test the eigenvector */
514:     if (d->eps->trueres && conv && c) {
515:       d->calcpairs_residual_eig(d,i,i+b,auxV);
516:       for (j=0,c=PETSC_TRUE; j<b && c; j++) {
517:         norm = d->nR[i+j]/d->nX[i+j];
518:         c = d->testConv(d, d->eigr[i+j], d->eigi[i+j], norm, &d->errest[i+j]);
519:       }
520:     }
521:     if (conv && c) { if (nConv) *nConv = i+b; }
522:     else conv = PETSC_FALSE;
523:   }
524:   pre = PetscMax(pre, i);

526: #if !defined(PETSC_USE_COMPLEX)
527:   /* Enforce converged conjugate complex eigenpairs */
528:   if (nConv) {
529:     for(j=0; j<*nConv; j++) if(d->eigi[j] != 0.0) j++;
530:     if(j > *nConv) (*nConv)--;
531:   }
532: #endif
533:   for(i=pre; i<e; i++) d->errest[i] = d->nR[i] = PETSC_MAX_REAL;
534: 
535:   return(0);
536: }