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-2010, Universidad 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

 28: PetscTruth dvd_isrestarting_fullV(dvdDashboard *d);
 29: PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d);
 30: PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d);
 31: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d);
 32: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d);
 33: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d);
 34: PetscErrorCode dvd_updateV_conv_finish(dvdDashboard *d);
 35: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d, PetscInt s, PetscInt pre,
 36:                                     PetscInt e, Vec *auxV, PetscScalar *auxS,
 37:                                     PetscInt *nConv);
 38: PetscErrorCode dvd_updateV_restartV_aux(Vec *V, PetscInt size_V,
 39:                                         PetscScalar *U, PetscInt ldU,
 40:                                         PetscScalar *pX, PetscInt ldpX,
 41:                                         PetscInt cpX, PetscScalar *oldpX,
 42:                                         PetscInt ldoldpX, PetscInt roldpX,
 43:                                         PetscInt coldpX, PetscScalar *auxS,
 44:                                         PetscInt size_auxS,
 45:                                         PetscInt *new_size_V);
 46: PetscErrorCode dvd_updateV_YtWx(PetscScalar *S, PetscInt ldS,
 47:                                 Vec *Y, PetscInt cY, Vec *y, PetscInt cy,
 48:                                 Vec *W, PetscInt cW, PetscScalar *x,
 49:                                 PetscInt ldx,
 50:                                 PetscInt rx, PetscInt cx, Vec *auxV);
 51: typedef struct {
 52:   PetscInt bs,      /* common number of approximated eigenpairs obtained */
 53:     real_max_size_V,
 54:                     /* real max size of V */
 55:     min_size_V,     /* restart with this number of eigenvectors */
 56:     plusk,          /* when restart, save plusk vectors from last iteration */
 57:     mpd;            /* max size of the searching subspace */
 58:   Vec *real_V,      /* real start vectors V */
 59:     *new_cY;        /* new left converged eigenvectors from the last iter */
 60:   void
 61:     *old_updateV_data;
 62:                     /* old updateV data */
 63:   isRestarting_type
 64:     old_isRestarting;
 65:                     /* old isRestarting */
 66:   PetscScalar
 67:     *oldU,          /* previous projected right igenvectors */
 68:     *oldV;          /* previous projected left eigenvectors */
 69:   PetscInt
 70:     ldoldU,         /* leading dimension of oldU */
 71:     size_oldU,      /* size of oldU */
 72:     size_new_cY;    /* size of new_cY */
 73:   PetscTruth
 74:     allResiduals;   /* if computing all the residuals */
 75: } dvdManagV_basic;

 79: PetscErrorCode dvd_managementV_basic(dvdDashboard *d, dvdBlackboard *b,
 80:                                      PetscInt bs, PetscInt max_size_V,
 81:                                      PetscInt mpd, PetscInt min_size_V,
 82:                                      PetscInt plusk, PetscTruth harm,
 83:                                      PetscTruth allResiduals)
 84: {
 85:   PetscErrorCode  ierr;
 86:   dvdManagV_basic *data;
 87:   PetscInt        i;


 91:   /* Setting configuration constrains */
 92:   mpd = PetscMin(mpd, max_size_V);
 93:   min_size_V = PetscMin(min_size_V, mpd-bs);
 94:   b->max_size_X = PetscMax(b->max_size_X, PetscMax(bs, min_size_V));
 95:   b->max_size_auxV = PetscMax(PetscMax(b->max_size_auxV,
 96:                                        b->max_size_X /* updateV_conv_gen */ ),
 97:                                        2 /* testConv */ );
 98:   b->max_size_auxS = PetscMax(PetscMax(PetscMax(b->max_size_auxS,
 99:                               max_size_V*b->max_size_X /* YtWx */ ),
100:                               max_size_V*2 /* SlepcDenseOrth  */ ),
101:                               max_size_V*b->max_size_X /* testConv:res_0 */ );
102:   b->max_size_V = mpd;
103:   b->size_V = max_size_V;
104:   b->own_vecs+= max_size_V*(harm==PETSC_TRUE?2:1);      /* V, W? */
105:   b->own_scalars+= max_size_V*2 /* eigr, eigr */ +
106:                    max_size_V /* nR */   +
107:                    max_size_V /* nX */   +
108:                    max_size_V /* errest */ +
109:                    2*b->max_size_V*b->max_size_V*(harm==PETSC_TRUE?2:1)
110:                                                /* MTX,MTY?,oldU,oldV? */;
111: //  b->max_size_oldX = plusk;

113:   /* Setup the step */
114:   if (b->state >= DVD_STATE_CONF) {
115:     PetscMalloc(sizeof(dvdManagV_basic), &data);
116:     data->real_max_size_V = max_size_V;
117:     data->min_size_V = min_size_V;
118:     data->real_V = b->free_vecs; b->free_vecs+= max_size_V;
119:     data->bs = bs;
120:     data->plusk = plusk;
121:     data->new_cY = PETSC_NULL;
122:     data->size_new_cY = 0;
123:     data->mpd = mpd;
124:     data->allResiduals = allResiduals;

126:     d->V = data->real_V;
127:     d->max_size_V = data->real_max_size_V;
128:     d->cX = data->real_V;
129:     d->eigr = b->free_scalars; b->free_scalars+= max_size_V;
130:     d->eigi = b->free_scalars; b->free_scalars+= max_size_V;
131: #ifdef PETSC_USE_COMPLEX
132:     for(i=0; i<max_size_V; i++) d->eigi[i] = 0.0;
133: #endif
134:     d->nR = (PetscReal*)b->free_scalars;
135:     b->free_scalars = (PetscScalar*)(d->nR + max_size_V);
136:     for(i=0; i<max_size_V; i++) d->nR[i] = PETSC_MAX;
137:     d->nX = (PetscReal*)b->free_scalars;
138:     b->free_scalars = (PetscScalar*)(d->nX + max_size_V);
139:     d->errest = (PetscReal*)b->free_scalars;
140:     b->free_scalars = (PetscScalar*)(d->errest + max_size_V);
141:     d->ceigr = d->eigr;
142:     d->ceigi = d->eigi;
143:     d->MTX = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
144:     data->oldU = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
145:     data->ldoldU = 0;
146:     data->oldV = PETSC_NULL;
147:     d->W = PETSC_NULL;
148:     d->MTY = PETSC_NULL;
149:     d->ldMTY = 0;
150:     if (harm == PETSC_TRUE) {
151:       d->W = b->free_vecs; b->free_vecs+= max_size_V;
152:       d->MTY = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
153:       data->oldV = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
154:     }

156:     data->size_oldU = 0;
157:     data->old_updateV_data = d->updateV_data;
158:     d->updateV_data = data;
159:     data->old_isRestarting = d->isRestarting;
160:     d->isRestarting = dvd_isrestarting_fullV;
161:     d->updateV = dvd_updateV_extrapol;
162:     d->preTestConv = dvd_updateV_testConv;
163:     DVD_FL_ADD(d->endList, dvd_updateV_conv_finish);
164:     DVD_FL_ADD(d->destroyList, dvd_managementV_basic_d);
165:   }

167:   return(0);
168: }

172: PetscTruth dvd_isrestarting_fullV(dvdDashboard *d)
173: {
174:   PetscTruth      restart;
175:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;


179:   /* Take into account the possibility of conjugate eigenpairs */
180: #if defined(PETSC_USE_COMPLEX)
181: #define ONE 0
182: #else
183: #define ONE 1
184: #endif

186:   restart = (d->size_V + data->bs + ONE > PetscMin(data->mpd,d->max_size_V))?
187:                 PETSC_TRUE:PETSC_FALSE;

189: #undef ONE

191:   /* Check old isRestarting function */
192:   if ((restart == PETSC_FALSE) && (data->old_isRestarting))
193:     restart = data->old_isRestarting(d);

195:   PetscFunctionReturn(restart);
196: }

200: PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
201: {
202:   PetscErrorCode  ierr;
203:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;


207:   /* Restore changes in dvdDashboard */
208:   d->updateV_data = data->old_updateV_data;
209: 
210:   /* Free local data */
211:   PetscFree(data);

213:   return(0);
214: }

218: PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
219: {
220:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
221:   PetscInt        i;
222:   PetscErrorCode  ierr;


226:   /// Temporal! Copy the converged left eigenvectors to cY
227:   if (data->size_new_cY > 0) {
228:     for (i=0; i<data->size_new_cY; i++) {
229:       VecCopy(data->new_cY[i], d->cY[d->size_cY+i]);
230:     }
231:     d->size_cY+= data->size_new_cY;
232:     data->size_new_cY = 0;
233:   }

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

237:   /* If the subspaces doesn't need restart, add new vector */
238:   if (d->isRestarting(d) == PETSC_FALSE) {
239:     i = d->size_V;
240:     dvd_updateV_update_gen(d);

242:     /* If some vector were add, exit */
243:     if (i < d->size_V) { return(0); }
244:   }

246:   /* If some eigenpairs were converged, lock them  */
247:   if (d->npreconv > 0) {
248:     i = d->npreconv;
249:     dvd_updateV_conv_gen(d);

251:     /* If some eigenpair was locked, exit */
252:     if (i > d->npreconv) { return(0); }
253:   }

255:   /* Else, a restarting is performed */
256:   dvd_updateV_restart_gen(d);

258:   return(0);
259: }

263: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
264: {
265:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
266:   PetscInt        i, j, npreconv, ldpX, inc_V, cMT, size_cX, size_cy;
267:   PetscErrorCode  ierr;
268:   PetscScalar     *pX;
269:   PetscReal       norm;
270:   Vec             *new_cY=0, *cX, *cy;
271:   PetscTruth      lindep;


275:   /* If the left subspace is present, constrains the converged pairs to the
276:      number of free vectors in V */
277:   if (d->cY && d->pY)
278:     npreconv = PetscMin(d->max_size_V-d->size_V, d->npreconv);
279:   else
280:     npreconv = d->npreconv;

282:   /* Constrains the converged pairs to nev */
283: #ifndef PETSC_USE_COMPLEX
284:   /* Tries to maintain together conjugate eigenpairs */
285:   for(i = 0;
286:       (i + (d->eigi[i]!=0.0?1:0) < npreconv) && (d->nconv + i < d->nev);
287:       i+= (d->eigi[i]!=0.0?2:1));
288:   npreconv = i;
289: #else
290:   npreconv = PetscMax(PetscMin(d->nev - d->nconv, npreconv), 0);
291: #endif

293:   if (d->npreconv == 0) { return(0); }

295:   /* f.cY <- [f.cY f.V*f.pY(0:npreconv-1)] */
296:   if (d->cY && d->pY) {
297:     new_cY = &d->V[d->size_V];
298:     SlepcUpdateVectorsZ(new_cY, 0.0, 1.0, d->W?d->W:d->V, d->size_V,
299:                                d->pY, d->ldpY, d->size_H, npreconv);
300: 

302:     /// Temporal! Postpone the copy of new_cY to cY
303:     data->new_cY = new_cY;
304:     data->size_new_cY = npreconv;
305:   }

307: #if !defined(PETSC_USE_COMPLEX)
308:   /* Correct the order of the conjugate eigenpairs */
309:   if (d->T) for (i=0; i<npreconv; i++) if (d->eigi[i] != 0.0) {
310:     if (d->eigi[i] < 0.0) {
311:       d->eigi[i]*= -1.0;
312:       d->eigi[i+1]*= -1.0;
313:       for (j=0; j<d->size_H; j++) d->pX[j+(i+1)*d->ldpX]*= -1.0;
314:       for (j=0; j<d->size_H; j++) d->S[j+(i+1)*d->ldS]*= -1.0;
315:       for (j=0; j<d->size_H; j++) d->T[j+(i+1)*d->ldT]*= -1.0;
316:     }
317:     i++;
318:   }
319: #endif

321:   /* BcX <- orth(BcX, B*cX),
322:      auxV = B*X, being X the last converged eigenvectors */
323:   if (d->BcX) for(i=0; i<npreconv; i++) {
324:     /* BcX <- [BcX auxV(i)] */
325:     VecCopy(d->auxV[i], d->BcX[d->size_cX+i]);
326:     IPOrthogonalize(d->ipI, 0, PETSC_NULL, d->size_cX+i, PETSC_NULL,
327:                            d->BcX, d->BcX[d->size_cX+i], PETSC_NULL,
328:                            &norm, &lindep);
329:     if(lindep == PETSC_TRUE) {
330:         SETERRQ(1, "Error during orth(BcX, B*cX(new))!");
331:     }
332:     VecScale(d->BcX[d->size_cX+i], 1.0/norm);
333:   }

335:   /* Harmonics restarts wiht right eigenvectors, and other with
336:      the left ones */
337:   pX = (d->W||!d->cY||d->BcX)?d->pX:d->pY;
338:   ldpX = (d->W||!d->cY||d->BcX)?d->ldpX:d->ldpY;

340:   /* If BcX, f.V <- orth(BcX, f.V) */
341:   if (d->BcX) cMT = 0;
342:   else        cMT = d->size_H - npreconv;

344:   /* f.MTX <- pY(npreconv:size_H-1), f.MTY <- f.pY(npreconv:size_H-1) */
345:   d->ldMTX = d->ldMTY = d->size_H;
346:   d->size_MT = d->size_H;
347:   d->MT_type = DVD_MT_ORTHO;
348:   SlepcDenseCopy(d->MTX, d->ldMTX, &pX[ldpX*npreconv], ldpX,
349:                         d->size_H, cMT);
350:   if (d->W && d->pY) {
351:     SlepcDenseCopy(d->MTY, d->ldMTY, &d->pY[d->ldpY*npreconv], d->ldpY,
352:                           d->size_H, cMT);
353:   }

355:   /* [f.cX(f.nconv) f.V] <- f.V*[f.pX(0:npreconv-1) f.pY(npreconv:f.size_V-1)] */
356:   if (&d->cX[d->nconv] == d->V) { /* cX and V are contiguous */
357:     SlepcDenseCopy(pX, ldpX, d->pX, d->ldpX, d->size_H, npreconv);
358: 
359:     SlepcUpdateVectorsZ(d->V, 0.0, 1.0, d->V, d->size_V, pX,
360:                                ldpX, d->size_H, d->size_H);
361:     d->V+= npreconv;
362:     inc_V = npreconv;
363:     d->max_size_V-= npreconv;
364:   } else {
365:     SETERRQ(1, "Untested case!");
366:     /*SlepcUpdateVectorsZ(&d->cX[d->nconv], 0.0, 1.0, d->V, d->size_V,
367:                                d->pX, d->ldpX, d->size_H, npreconv);
368:     
369:     SlepcUpdateVectorsZ(d->V, 0.0, 1.0, d->V, d->size_V,
370:                                &pX[ldpX*npreconv], ldpX,
371:                                d->size_H, cMT); 
372:     inc_V = 0;*/
373:   }
374:   d->size_cX+= npreconv;
375:   d->size_V -= npreconv;

377:   /* Udpate cS and cT, if needed */
378:   if (d->cS) {
379:     PetscInt size_cS = d->size_cX-npreconv;
380:     cX = d->cY?d->cY:d->cX; size_cX = d->cY?d->size_cY:d->size_cX;
381:     cy = d->cY?new_cY:0; size_cy = d->cY?npreconv:0;
382: 
383:     dvd_updateV_YtWx(&d->cS[d->ldcS*size_cS], d->ldcS, cX, size_cX, cy,
384:                      size_cy, d->AV, d->size_AV,
385:                      d->pX, d->ldpX,
386:                      d->size_H, npreconv, d->auxV);

388:     if (DVD_ISNOT(d->sEP, DVD_EP_STD)) {if (d->BV) {
389: 
390:       dvd_updateV_YtWx(&d->cT[d->ldcT*size_cS], d->ldcT, cX, size_cX, cy,
391:                        size_cy, d->BV, d->size_AV,
392:                        d->pX, d->ldpX,
393:                        d->size_H, npreconv, d->auxV);
394:     } else if (!d->B) {
395:       VecsMultIa(&d->cT[d->ldcT*size_cS], 0, d->ldcT, cX, 0, size_cX,
396:                         &d->cX[size_cS], 0, npreconv);
397:       VecsMultIa(&d->cT[d->ldcT*size_cS+size_cX], 0, d->ldcT, cy, 0,
398:                         size_cy, &d->cX[size_cS], 0, npreconv);
399:     } else {
400:       /* TODO: Only for nprecond==1 */
401:       MatMult(d->B, d->cX[d->size_cX-1], d->auxV[0]);
402:       VecsMultIa(&d->cT[d->ldcT*size_cS], 0, d->ldcT, cX, 0, size_cX,
403:                         &d->auxV[0], 0, npreconv);
404:       VecsMultIa(&d->cT[d->ldcT*size_cS+size_cX], 0, d->ldcT, cy, 0,
405:                         size_cy, &d->auxV[0], 0, npreconv);
406:     }}
407:   }

409:   /* f.W <- f.W * f.MTY */
410:   if (d->W) {
411:     SlepcUpdateVectorsZ(d->W, 0.0, 1.0, d->W, d->size_V+npreconv,
412:                                d->MTY, d->ldMTY, d->size_H, cMT);
413: 
414:   }

416:   /* Lock the converged pairs */
417:   d->eigr+= npreconv;
418: #ifndef PETSC_USE_COMPLEX
419:   if (d->eigi) d->eigi+= npreconv;
420: #endif
421:   d->nconv+= npreconv;
422:   d->errest+= npreconv;

424:   /* Notify the changes in V and update the other subspaces */
425:   d->V_imm_s = inc_V;             d->V_imm_e = inc_V;
426:   d->V_tra_s = 0;                 d->V_tra_e = cMT;
427:   d->V_new_s = d->V_tra_e;        d->V_new_e = d->size_V;

429:   /* Remove oldU */
430:   data->size_oldU = 0;

432:   /* f.pX <- I, f.pY <- I */
433:   d->ldpX = d->size_H;
434:   if (d->pY) d->ldpY = d->size_H;
435:   for(i=0; i<d->size_H; i++) {
436:     for(j=0; j<d->size_H; j++) {
437:       d->pX[d->ldpX*i+j] = 0.0;
438:       if (d->pY) d->pY[d->ldpY*i+j] = 0.0;
439:     }
440:     d->pX[d->ldpX*i+i] = 1.0;
441:     if (d->pY) d->pY[d->ldpY*i+i] = 1.0;
442:   }

444:   d->npreconv-= npreconv;

446:   return(0);
447: }

451: PetscErrorCode dvd_updateV_conv_finish(dvdDashboard *d)
452: {
453:   PetscErrorCode  ierr;
454: #if defined(PETSC_USE_COMPLEX)
455:   PetscInt        i, j;
456:   PetscScalar     s;
457: #endif  


461:   /* Finish cS and cT */
462:   VecsMultIb(d->cS, 0, d->ldcS, d->nconv, d->nconv, d->auxS, d->V[0]);
463: 
464:   if (d->cT) {
465:     VecsMultIb(d->cT, 0, d->ldcT, d->nconv, d->nconv, d->auxS, d->V[0]);
466: 
467:   }

469:   /* Some functions need the diagonal elements in cT be real */
470: #if defined(PETSC_USE_COMPLEX)
471:   if (d->cT) for(i=0; i<d->nconv; i++) {
472:     s = PetscConj(d->cT[d->ldcT*i+i])/PetscAbsScalar(d->cT[d->ldcT*i+i]);
473:     for(j=0; j<=i; j++)
474:       d->cT[d->ldcT*i+j] = PetscRealPart(d->cT[d->ldcT*i+j]*s),
475:       d->cS[d->ldcS*i+j]*= s;
476:     VecScale(d->cX[i], s);
477:   }
478: #endif

480:   return(0);
481: }
482: 
485: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
486: {
487:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
488:   PetscInt        size_plusk, size_X, new_size_X, i, j;
489:   PetscErrorCode  ierr;


493:   /* Select size_X desired pairs from V */
494:   size_X = PetscMin(PetscMin(data->min_size_V,
495:                              d->size_V ),
496:                              d->max_size_V );

498:   /* Add plusk eigenvectors from the previous iteration */
499:   size_plusk = PetscMax(0, PetscMin(PetscMin(data->plusk,
500:                                     data->size_oldU ),
501:                                     d->max_size_V - size_X ));

503:   /* Modify the subspaces */
504:   d->ldMTX = d->size_MT = d->size_V;
505:   dvd_updateV_restartV_aux(d->V, d->size_V, d->MTX, d->ldMTX, d->pX,
506:                                   d->ldpX, size_X, data->oldU, data->ldoldU,
507:                                   data->size_oldU, size_plusk, d->auxS,
508:                                   d->size_auxS, &new_size_X);
509:   if (d->W && d->pY) {
510:     PetscInt new_size_Y;
511:     d->ldMTY = d->size_V;
512:     dvd_updateV_restartV_aux(d->W, d->size_V, d->MTY, d->ldMTY, d->pY,
513:                                     d->ldpY, size_X, data->oldV, data->ldoldU,
514:                                     data->size_oldU, new_size_X-size_X,
515:                                     d->auxS, d->size_auxS, &new_size_Y);
516: 
517:     new_size_X = PetscMin(new_size_X, new_size_Y);
518:   }

520:   /* Notify the changes in V and update the other subspaces */
521:   d->size_V = new_size_X;
522:   d->MT_type = DVD_MT_ORTHO;
523:   d->V_imm_s = 0;                 d->V_imm_e = 0;
524:   d->V_tra_s = 0;                 d->V_tra_e = new_size_X;
525:   d->V_new_s = d->V_tra_e;        d->V_new_e = d->V_tra_e;

527:   /* Remove oldU */
528:   data->size_oldU = 0;

530:   /* Remove npreconv */
531:   d->npreconv = 0;
532: 
533:   /* f.pX <- I, f.pY <- I */
534:   d->ldpX = d->size_H;
535:   if (d->pY) d->ldpY = d->size_H;
536:   for(i=0; i<d->size_H; i++) {
537:     for(j=0; j<d->size_H; j++) {
538:       d->pX[d->ldpX*i+j] = 0.0;
539:       if (d->pY) d->pY[d->ldpY*i+j] = 0.0;
540:     }
541:     d->pX[d->ldpX*i+i] = 1.0;
542:     if (d->pY) d->pY[d->ldpY*i+i] = 1.0;
543:   }

545:   return(0);
546: }


551: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
552: {
553:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
554:   PetscInt        size_D;
555:   PetscErrorCode  ierr;


559:   /* Select the desired pairs */
560:   size_D = PetscMin(PetscMin(PetscMin(data->bs,
561:                                       d->size_V ),
562:                                       d->max_size_V-d->size_V ),
563:                                       d->size_H );
564:   if (size_D == 0) {
565:     PetscPrintf(PETSC_COMM_WORLD, "MON: D:%d H:%d\n", size_D, d->size_H);
566:     d->initV(d);
567:     d->calcPairs(d);
568:     //SETERRQ(1, "D == 0!\n");
569:     //PetscFunctionReturn(1);
570:   }

572: //  PetscPrintf(PETSC_COMM_WORLD, "EIGS: ");
573: //  for(i=0; i<d->size_H; i++) PetscPrintf(PETSC_COMM_WORLD, "%d:%g ", i, d->eigr[i]);
574: //  PetscPrintf(PETSC_COMM_WORLD, "\n");

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

580:   /* If D is empty, exit */
581:   if (size_D == 0) { return(0); }

583:   /* Get the converged pairs */
584:   dvd_updateV_testConv(d, 0, size_D,
585:     data->allResiduals==PETSC_TRUE?d->size_V:size_D, d->auxV, d->auxS,
586:     &d->npreconv);

588:   /* Notify the changes in V */
589:   d->size_V+= size_D;
590:   d->size_D = size_D;
591:   d->V_imm_s = 0;                 d->V_imm_e = d->size_V-size_D;
592:   d->V_tra_s = 0;                 d->V_tra_e = 0;
593:   d->V_new_s = d->size_V-size_D;  d->V_new_e = d->size_V;

595:   /* Save the projected eigenvectors */
596:   if (data->plusk > 0) {
597:     data->ldoldU = data->size_oldU = d->size_H;
598:     SlepcDenseCopy(data->oldU, data->ldoldU, d->pX, d->ldpX, d->size_H,
599:                           d->size_H);
600:     if (d->pY && d->W) {
601:       SlepcDenseCopy(data->oldV, data->ldoldU, d->pY, d->ldpY, d->size_H,
602:                             d->size_H);
603:     }
604:   }

606:   return(0);
607: }


612: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d, PetscInt s, PetscInt pre,
613:                                     PetscInt e, Vec *auxV, PetscScalar *auxS,
614:                                                       PetscInt *nConv)
615: {
616:   PetscInt        i;
617: #ifndef PETSC_USE_COMPLEX
618:   PetscInt        j;
619: #endif
620:   PetscReal       norm;
621:   PetscErrorCode  ierr;
622:   PetscTruth      conv, c;
623:   dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;

626: 
627:   *nConv = s;
628:   for(i=s, conv=PETSC_TRUE;
629:       (conv == PETSC_TRUE || data->allResiduals == PETSC_TRUE) && (i < e);
630:       i++) {
631:     if (i >= pre) {
632:       d->calcpairs_X(d, i, i+1, &auxV[1]);
633:       if ((d->B) && DVD_IS(d->sEP, DVD_EP_STD)) {
634:         MatMult(d->B, auxV[1], auxV[0]);
635:       }
636:       d->calcpairs_residual(d, i, i+1, auxV, auxS, auxV[1]);
637: 
638:     }
639:     norm = d->nR[i]/d->nX[i];
640:     c = d->testConv(d, d->eigr[i], d->eigi[i], norm, &d->errest[i]);
641:     if (conv == PETSC_TRUE && c == PETSC_TRUE) *nConv = i+1;
642:     else conv = PETSC_FALSE;
643:   }
644:   pre = PetscMax(pre, i);

646: #ifndef PETSC_USE_COMPLEX
647:   /* Enforce converged conjugate conjugate complex eigenpairs */
648:   for(j=0; j<*nConv; j++) if(d->eigi[j] != 0.0) j++;
649:   if(j > *nConv) (*nConv)--;
650: #endif
651:   for(i=pre; i<e; i++) d->errest[i] = d->nR[i] = PETSC_MAX;
652: 
653:   return(0);
654: }


657: /*
658:   U <- [pX(0:size_X-1) gs(pX(0:size_X-1), oldpX(0:size_plusk-1))]
659:   V <- V * U,
660:   where
661:   new_size_V, return the new size of V
662:   auxS, auxiliar vector of size 2*ldpX, at least
663:   size_auxS, the size of auxS
664: */
667: PetscErrorCode dvd_updateV_restartV_aux(Vec *V, PetscInt size_V,
668:                                         PetscScalar *U, PetscInt ldU,
669:                                         PetscScalar *pX, PetscInt ldpX,
670:                                         PetscInt cpX, PetscScalar *oldpX,
671:                                         PetscInt ldoldpX, PetscInt roldpX,
672:                                         PetscInt coldpX, PetscScalar *auxS,
673:                                         PetscInt size_auxS,
674:                                         PetscInt *new_size_V)
675: {
676:   PetscInt        i, j;
677:   PetscErrorCode  ierr;


681:   /* Add the size_X best eigenpairs */
682:   SlepcDenseCopy(U, ldU, pX, ldpX, size_V, cpX);

684:   /* Add plusk eigenvectors from the previous iteration */
685:   SlepcDenseCopy(&U[cpX*ldU], ldU, oldpX, ldoldpX, roldpX, coldpX);
686: 
687:   for(i=cpX; i<cpX+coldpX; i++)
688:     for(j=roldpX; j<size_V; j++)
689:         U[i*ldU+j] = 0.0;

691:   /* U <- orth(U) */
692:   /// Temporal! Correct sentence: U <- orth(U(0:size_X-1), U(size_X:size_X+size_plusk))
693:   if (coldpX > 0) {
694:     SlepcDenseOrth(U, ldU, size_V, cpX+coldpX, auxS, size_auxS,
695:                           new_size_V);
696:   } else
697:     *new_size_V = cpX;

699:   /* V <- V * U */
700:   SlepcUpdateVectorsZ(V, 0.0, 1.0, V, size_V, U, ldU, size_V,
701:                              *new_size_V);

703:   return(0);
704: }

706: /*
707:   Compute S = [ Y' * W * x
708:                 y' * W * x ]
709:   where
710:   ldS, the leading dimension of S,
711:   cY, number of vectors of Y,
712:   ldy, the leading dimension of y,
713:   ry,cy, rows and columns of y,
714:   cW, number of vectors of W,
715:   ldH, the leading dimension of H,
716:   rH,cH, rows and columns of H,
717:   ldx, the leading dimension of y,
718:   rx,cx, rows and columns of x,
719:   r, a reduction,
720:   sr, a permanent space,
721:   auxV, array of auxiliar vectors of size cx (at the end, auxV <- W*x),
722:   auxS, auxiliar scalar vector of size rH*cx.
723: */
726: PetscErrorCode dvd_updateV_YtWx(PetscScalar *S, PetscInt ldS,
727:                                 Vec *Y, PetscInt cY, Vec *y, PetscInt cy,
728:                                 Vec *W, PetscInt cW, PetscScalar *x,
729:                                 PetscInt ldx,
730:                                 PetscInt rx, PetscInt cx, Vec *auxV)
731: {
732:   PetscErrorCode  ierr;


736:   /* auxV <- W * x */
737:   SlepcUpdateVectorsZ(auxV, 0.0, 1.0, W, cW, x, ldx, rx, cx);
738: 

740:   /* S(0:cY-1, 0:cx-1) <- Y' * auxV */
741:   VecsMultIa(S, 0, ldS, Y, 0, cY, auxV, 0, cx);

743:   /* S(cY:cY+cy-1, 0:cx-1) <- y' * auxV */
744:   VecsMultIa(&S[cY], 0, ldS, y, 0, cy, auxV, 0, cx);

746:   return(0);
747: }