Actual source code: davidson.c

  1: /*
  2:   Method: General Davidson Method (includes GD and JD)

  4:   References:
  5:     - Ernest R. Davidson. Super-matrix methods. Computer Physics Communications,
  6:       53:49–60, May 1989.

  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 10:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

 12:    This file is part of SLEPc.
 13:       
 14:    SLEPc is free software: you can redistribute it and/or modify it under  the
 15:    terms of version 3 of the GNU Lesser General Public License as published by
 16:    the Free Software Foundation.

 18:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 19:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 20:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 21:    more details.

 23:    You  should have received a copy of the GNU Lesser General  Public  License
 24:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 25:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 26: */

 28:  #include davidson.h

 30: PetscErrorCode EPSView_Davidson(EPS eps,PetscViewer viewer);

 32: typedef struct {
 33:   /**** Solver options ****/
 34:   PetscInt blocksize,     /* block size */
 35:     initialsize,          /* initial size of V */
 36:     minv,                 /* size of V after restarting */
 37:     plusk;                /* keep plusk eigenvectors from the last iteration */
 38:   PetscBool  ipB;         /* true if V'B*V=I */
 39:   PetscInt   method;      /* method for improving the approximate solution */
 40:   PetscReal  fix;         /* the fix parameter */
 41:   PetscBool  krylovstart; /* true if the starting subspace is a Krylov basis */
 42:   PetscBool  dynamic;     /* true if dynamic stopping criterion is used */
 43:   PetscInt   cX_in_proj,  /* converged vectors in the projected problem */
 44:     cX_in_impr;           /* converged vectors in the projector */

 46:   /**** Solver data ****/
 47:   dvdDashboard ddb;

 49:   /**** Things to destroy ****/
 50:   PetscScalar *wS;
 51:   Vec         *wV;
 52:   PetscInt    size_wV;
 53: } EPS_DAVIDSON;

 57: PetscErrorCode EPSCreate_Davidson(EPS eps)
 58: {
 60:   EPS_DAVIDSON   *data;


 64:   eps->OP->ops->getbilinearform  = STGetBilinearForm_Default;
 65:   eps->ops->solve                = EPSSolve_Davidson;
 66:   eps->ops->setup                = EPSSetUp_Davidson;
 67:   eps->ops->reset                = EPSReset_Davidson;
 68:   eps->ops->backtransform        = EPSBackTransform_Default;
 69:   eps->ops->computevectors       = EPSComputeVectors_Davidson;
 70:   eps->ops->view                 = EPSView_Davidson;

 72:   PetscMalloc(sizeof(EPS_DAVIDSON),&data);
 73:   eps->data = data;
 74:   data->wS = PETSC_NULL;
 75:   data->wV = PETSC_NULL;
 76:   data->size_wV = 0;
 77:   PetscMemzero(&data->ddb,sizeof(dvdDashboard));

 79:   /* Set default values */
 80:   EPSDavidsonSetKrylovStart_Davidson(eps,PETSC_FALSE);
 81:   EPSDavidsonSetBlockSize_Davidson(eps,1);
 82:   EPSDavidsonSetRestart_Davidson(eps,6,0);
 83:   EPSDavidsonSetInitialSize_Davidson(eps,5);
 84:   EPSDavidsonSetFix_Davidson(eps,0.01);
 85:   EPSDavidsonSetBOrth_Davidson(eps,PETSC_TRUE);
 86:   EPSDavidsonSetConstantCorrectionTolerance_Davidson(eps,PETSC_TRUE);
 87:   EPSDavidsonSetWindowSizes_Davidson(eps,0,0);
 88:   return(0);
 89: }

 93: PetscErrorCode EPSSetUp_Davidson(EPS eps)
 94: {
 96:   EPS_DAVIDSON   *data = (EPS_DAVIDSON*)eps->data;
 97:   dvdDashboard   *dvd = &data->ddb;
 98:   dvdBlackboard  b;
 99:   PetscInt       nvecs,nscalars,min_size_V,plusk,bs,initv,i,cX_in_proj,cX_in_impr;
100:   Mat            A,B;
101:   KSP            ksp;
102:   PetscBool      t,ipB,ispositive,dynamic;
103:   HarmType_t     harm;
104:   InitType_t     init;
105:   PetscReal      fix;
106:   PetscScalar    target;

109:   /* Setup EPS options and get the problem specification */
110:   EPSDavidsonGetBlockSize_Davidson(eps,&bs);
111:   if (bs <= 0) bs = 1;
112:   if(eps->ncv) {
113:     if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The value of ncv must be at least nev");
114:   } else if (eps->mpd) eps->ncv = eps->mpd + eps->nev + bs;
115:   else if (eps->nev<500)
116:     eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15))+bs;
117:   else
118:     eps->ncv = PetscMin(eps->n,eps->nev+500)+bs;
119:   if (!eps->mpd) eps->mpd = eps->ncv;
120:   if (eps->mpd > eps->ncv)
121:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The mpd has to be less or equal than ncv");
122:   if (eps->mpd < 2)
123:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The mpd has to be greater than 2");
124:   if (!eps->max_it) eps->max_it = PetscMax(100*eps->ncv,2*eps->n);
125:   if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
126:   if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
127:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Wrong value of eps->which");
128:   if (!(eps->nev + bs <= eps->ncv))
129:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The ncv has to be greater than nev plus blocksize!");

131:   EPSDavidsonGetRestart_Davidson(eps,&min_size_V,&plusk);
132:   if (!min_size_V) min_size_V = PetscMin(PetscMax(bs,5),eps->mpd/2);
133:   if (!(min_size_V+bs <= eps->mpd))
134:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The value of minv must be less than mpd minus blocksize");
135:   EPSDavidsonGetInitialSize_Davidson(eps,&initv);
136:   if (eps->mpd < initv)
137:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The initv has to be less or equal than mpd");

139:   /* Davidson solvers do not support left eigenvectors */
140:   if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");

142:   /* Set STPrecond as the default ST */
143:   if (!((PetscObject)eps->OP)->type_name) {
144:     STSetType(eps->OP,STPRECOND);
145:   }
146:   STPrecondSetKSPHasMat(eps->OP,PETSC_FALSE);

148:   /* Change the default sigma to inf if necessary */
149:   if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
150:       eps->which == EPS_LARGEST_IMAGINARY) {
151:     STSetDefaultShift(eps->OP,PETSC_MAX_REAL);
152:   }
153: 
154:   /* Davidson solvers only support STPRECOND */
155:   STSetUp(eps->OP);
156:   PetscTypeCompare((PetscObject)eps->OP,STPRECOND,&t);
157:   if (!t) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_SUP,"%s only works with precond spectral transformation",
158:     ((PetscObject)eps)->type_name);

160:   /* Setup problem specification in dvd */
161:   STGetOperators(eps->OP,&A,&B);
162:   EPSReset_Davidson(eps);
163:   PetscMemzero(dvd,sizeof(dvdDashboard));
164:   dvd->A = A; dvd->B = eps->isgeneralized? B : PETSC_NULL;
165:   ispositive = eps->ispositive;
166:   dvd->sA = DVD_MAT_IMPLICIT |
167:             (eps->ishermitian? DVD_MAT_HERMITIAN : 0) |
168:             ((ispositive && !eps->isgeneralized) ? DVD_MAT_POS_DEF : 0);
169:   /* Asume -eps_hermitian means hermitian-definite in generalized problems */
170:   if (!ispositive && !eps->isgeneralized && eps->ishermitian) ispositive = PETSC_TRUE;
171:   if (!eps->isgeneralized)
172:     dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY |
173:               DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
174:   else
175:     dvd->sB = DVD_MAT_IMPLICIT |
176:               (eps->ishermitian? DVD_MAT_HERMITIAN : 0) |
177:               (ispositive? DVD_MAT_POS_DEF : 0);
178:   ipB = (dvd->B && data->ipB && DVD_IS(dvd->sB,DVD_MAT_POS_DEF))?PETSC_TRUE:PETSC_FALSE;
179:   data->ipB = ipB;
180:   dvd->correctXnorm = ipB;
181:   dvd->sEP = ((!eps->isgeneralized || (eps->isgeneralized && ipB))? DVD_EP_STD : 0) |
182:              (ispositive? DVD_EP_HERMITIAN : 0);
183:   dvd->nev = eps->nev;
184:   dvd->which = eps->which;
185:   dvd->withTarget = PETSC_TRUE;
186:   switch(eps->which) {
187:   case EPS_TARGET_MAGNITUDE:
188:   case EPS_TARGET_IMAGINARY:
189:     dvd->target[0] = target = eps->target; dvd->target[1] = 1.0;
190:     break;

192:   case EPS_TARGET_REAL:
193:     dvd->target[0] = PetscRealPart(target = eps->target); dvd->target[1] = 1.0;
194:     break;

196:   case EPS_LARGEST_REAL:
197:   case EPS_LARGEST_MAGNITUDE:
198:   case EPS_LARGEST_IMAGINARY: /* TODO: think about this case */
199:   default:
200:     dvd->target[0] = 1.0; dvd->target[1] = target = 0.0;
201:     break;
202: 
203:   case EPS_SMALLEST_MAGNITUDE:
204:   case EPS_SMALLEST_REAL:
205:   case EPS_SMALLEST_IMAGINARY: /* TODO: think about this case */
206:     dvd->target[0] = target = 0.0; dvd->target[1] = 1.0;
207:     break;

209:   case EPS_WHICH_USER:
210:     STGetShift(eps->OP,&target);
211:     dvd->target[0] = target; dvd->target[1] = 1.0;
212:     break;

214:   case EPS_ALL:
215:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported option: which == EPS_ALL");
216:     break;
217:   }
218:   dvd->tol = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
219:   dvd->eps = eps;

221:   /* Setup the extraction technique */
222:   if (!eps->extraction) {
223:     if (ipB || ispositive)
224:       eps->extraction = EPS_RITZ;
225:     else if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_IMAGINARY || eps->which == EPS_LARGEST_REAL)
226:       eps->extraction = EPS_HARMONIC_LARGEST;
227:     else
228:       eps->extraction = EPS_RITZ;
229:   }
230:   switch(eps->extraction) {
231:   case EPS_RITZ:              harm = DVD_HARM_NONE; break;
232:   case EPS_HARMONIC:          harm = DVD_HARM_RR; break;
233:   case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
234:   case EPS_HARMONIC_RIGHT:    harm = DVD_HARM_REIGS; break;
235:   case EPS_HARMONIC_LARGEST:  harm = DVD_HARM_LEIGS; break;
236:   default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
237:   }

239:   /* Setup the type of starting subspace */
240:   EPSDavidsonGetKrylovStart_Davidson(eps,&t);
241:   init = (!t)? DVD_INITV_CLASSIC : DVD_INITV_KRYLOV;

243:   /* Setup the presence of converged vectors in the projected problem and in the projector */
244:   EPSDavidsonGetWindowSizes_Davidson(eps,&cX_in_impr,&cX_in_proj);

246:   /* Setup IP */
247:   if (ipB && dvd->B) {
248:     IPSetMatrix(eps->ip,dvd->B);
249:   } else {
250:     IPSetMatrix(eps->ip,PETSC_NULL);
251:   }

253:   /* Get the fix parameter */
254:   EPSDavidsonGetFix_Davidson(eps,&fix);

256:   /* Get whether the stopping criterion is used */
257:   EPSDavidsonGetConstantCorrectionTolerance_Davidson(eps,&dynamic);

259:   /* Orthonormalize the DS */
260:   dvd_orthV(eps->ip,PETSC_NULL,0,PETSC_NULL,0,eps->DS,0,
261:                    PetscAbs(eps->nds),PETSC_NULL,eps->rand);

263:   /* Preconfigure dvd */
264:   STGetKSP(eps->OP,&ksp);
265:   dvd_schm_basic_preconf(dvd,&b,eps->mpd,min_size_V,bs,
266:                                 initv,
267:                                 PetscAbs(eps->nini),
268:                                 plusk,harm,
269:                                 ksp,init,eps->trackall,
270:                                 ipB?DVD_ORTHOV_BOneMV:DVD_ORTHOV_I,cX_in_proj,cX_in_impr);
271: 

273:   /* Allocate memory */
274:   nvecs = b.max_size_auxV + b.own_vecs;
275:   nscalars = b.own_scalars + b.max_size_auxS;
276:   PetscMalloc(nscalars*sizeof(PetscScalar),&data->wS);
277:   VecDuplicateVecs(eps->t,nvecs,&data->wV);
278:   data->size_wV = nvecs;
279:   b.free_vecs = data->wV;
280:   b.free_scalars = data->wS;
281:   dvd->auxV = data->wV + b.own_vecs;
282:   dvd->auxS = b.free_scalars + b.own_scalars;
283:   dvd->size_auxV = b.max_size_auxV;
284:   dvd->size_auxS = b.max_size_auxS;

286:   eps->errest_left = PETSC_NULL;
287:   PetscMalloc(eps->ncv*sizeof(PetscInt),&eps->perm);
288:   for(i=0; i<eps->ncv; i++) eps->perm[i] = i;

290:   /* Configure dvd for a basic GD */
291:   dvd_schm_basic_conf(dvd,&b,eps->mpd,min_size_V,bs,
292:                              initv,
293:                              PetscAbs(eps->nini),plusk,
294:                              eps->ip,harm,dvd->withTarget,
295:                              target,ksp,
296:                              fix,init,eps->trackall,
297:                              ipB?DVD_ORTHOV_BOneMV:DVD_ORTHOV_I,cX_in_proj,cX_in_impr,dynamic);
298: 

300:   /* Associate the eigenvalues to the EPS */
301:   eps->eigr = dvd->real_eigr;
302:   eps->eigi = dvd->real_eigi;
303:   eps->errest = dvd->real_errest;
304:   eps->V = dvd->real_V;

306: 
307:   return(0);
308: }

312: PetscErrorCode EPSSolve_Davidson(EPS eps)
313: {
314:   EPS_DAVIDSON   *data = (EPS_DAVIDSON*)eps->data;
315:   dvdDashboard   *d = &data->ddb;

319:   /* Call the starting routines */
320:   DVD_FL_CALL(d->startList,d);

322:   for(eps->its=0; eps->its < eps->max_it; eps->its++) {
323:     /* Initialize V, if it is needed */
324:     if (d->size_V == 0) { d->initV(d); }

326:     /* Find the best approximated eigenpairs in V, X */
327:     d->calcPairs(d);

329:     /* Test for convergence */
330:     if (eps->nconv >= eps->nev) break;

332:     /* Expand the subspace */
333:     d->updateV(d);

335:     /* Monitor progress */
336:     eps->nconv = d->nconv;
337:     EPSMonitor(eps,eps->its+1,eps->nconv,eps->eigr,eps->eigi,eps->errest,d->size_V+d->size_cX);
338:   }

340:   /* Call the ending routines */
341:   DVD_FL_CALL(d->endList,d);

343:   if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
344:   else eps->reason = EPS_DIVERGED_ITS;

346:   return(0);
347: }

351: PetscErrorCode EPSReset_Davidson(EPS eps)
352: {
353:   EPS_DAVIDSON   *data = (EPS_DAVIDSON*)eps->data;
354:   dvdDashboard   *dvd = &data->ddb;

358:   /* Call step destructors and destroys the list */
359:   DVD_FL_CALL(dvd->destroyList,dvd);
360:   DVD_FL_DEL(dvd->destroyList);
361:   DVD_FL_DEL(dvd->startList);
362:   DVD_FL_DEL(dvd->endList);

364:   if (data->size_wV > 0) {
365:     VecDestroyVecs(data->size_wV,&data->wV);
366:   }
367:   PetscFree(data->wS);
368:   PetscFree(eps->perm);
369:   return(0);
370: }

374: PetscErrorCode EPSView_Davidson(EPS eps,PetscViewer viewer)
375: {
377:   PetscBool      isascii,opb;
378:   PetscInt       opi,opi0;
379:   const char*    name;

382:   name = ((PetscObject)eps)->type_name;
383:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
384:   if (!isascii) {
385:     SETERRQ2(((PetscObject)eps)->comm,1,"Viewer type %s not supported for %s",((PetscObject)viewer)->type_name,name);
386:   }
387: 
388:   EPSDavidsonGetBOrth_Davidson(eps,&opb);
389:   EPSDavidsonGetBlockSize_Davidson(eps,&opi);
390:   if(!opb) {
391:     PetscViewerASCIIPrintf(viewer,"  Davidson: search subspace is I-orthogonalized\n");
392:   } else {
393:     PetscViewerASCIIPrintf(viewer,"  Davidson: search subspace is B-orthogonalized\n");
394:   }
395:   PetscViewerASCIIPrintf(viewer,"  Davidson: block size=%D\n",opi);
396:   EPSDavidsonGetKrylovStart_Davidson(eps,&opb);
397:   if(!opb) {
398:     PetscViewerASCIIPrintf(viewer,"  Davidson: type of the initial subspace: non-Krylov\n");
399:   } else {
400:     PetscViewerASCIIPrintf(viewer,"  Davidson: type of the initial subspace: Krylov\n");
401:   }
402:   EPSDavidsonGetRestart_Davidson(eps,&opi,&opi0);
403:   PetscViewerASCIIPrintf(viewer,"  Davidson: size of the subspace after restarting: %D\n",opi);
404:   PetscViewerASCIIPrintf(viewer,"  Davidson: number of vectors after restarting from the previous iteration: %D\n",opi0);
405:   return(0);
406: }

410: PetscErrorCode EPSDavidsonSetKrylovStart_Davidson(EPS eps,PetscBool krylovstart)
411: {
412:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

415:   data->krylovstart = krylovstart;
416:   return(0);
417: }

421: PetscErrorCode EPSDavidsonGetKrylovStart_Davidson(EPS eps,PetscBool *krylovstart)
422: {
423:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

426:   *krylovstart = data->krylovstart;
427:   return(0);
428: }

432: PetscErrorCode EPSDavidsonSetBlockSize_Davidson(EPS eps,PetscInt blocksize)
433: {
434:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

437:   if(blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
438:   if(blocksize <= 0)
439:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value");
440:   data->blocksize = blocksize;
441:   return(0);
442: }

446: PetscErrorCode EPSDavidsonGetBlockSize_Davidson(EPS eps,PetscInt *blocksize)
447: {
448:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

451:   *blocksize = data->blocksize;
452:   return(0);
453: }

457: PetscErrorCode EPSDavidsonSetRestart_Davidson(EPS eps,PetscInt minv,PetscInt plusk)
458: {
459:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

462:   if(minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 5;
463:   if(minv <= 0)
464:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value");
465:   if(plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = 5;
466:   if(plusk < 0)
467:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value");
468:   data->minv = minv;
469:   data->plusk = plusk;
470:   return(0);
471: }

475: PetscErrorCode EPSDavidsonGetRestart_Davidson(EPS eps,PetscInt *minv,PetscInt *plusk)
476: {
477:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

480:   *minv = data->minv;
481:   *plusk = data->plusk;
482:   return(0);
483: }

487: PetscErrorCode EPSDavidsonGetInitialSize_Davidson(EPS eps,PetscInt *initialsize)
488: {
489:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

492:   *initialsize = data->initialsize;
493:   return(0);
494: }

498: PetscErrorCode EPSDavidsonSetInitialSize_Davidson(EPS eps,PetscInt initialsize)
499: {
500:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

503:   if(initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 5;
504:   if(initialsize <= 0)
505:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value");
506:   data->initialsize = initialsize;
507:   return(0);
508: }

512: PetscErrorCode EPSDavidsonGetFix_Davidson(EPS eps,PetscReal *fix)
513: {
514:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

517:   *fix = data->fix;
518:   return(0);
519: }

523: PetscErrorCode EPSDavidsonSetFix_Davidson(EPS eps,PetscReal fix)
524: {
525:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

528:   if(fix == PETSC_DEFAULT || fix == PETSC_DECIDE) fix = 0.01;
529:   if(fix < 0.0)
530:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
531:   data->fix = fix;
532:   return(0);
533: }

537: PetscErrorCode EPSDavidsonSetBOrth_Davidson(EPS eps,PetscBool borth)
538: {
539:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

542:   data->ipB = borth;
543:   return(0);
544: }

548: PetscErrorCode EPSDavidsonGetBOrth_Davidson(EPS eps,PetscBool *borth)
549: {
550:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

553:   *borth = data->ipB;
554:   return(0);
555: }

559: PetscErrorCode EPSDavidsonSetConstantCorrectionTolerance_Davidson(EPS eps,PetscBool constant)
560: {
561:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

564:   data->dynamic = !constant?PETSC_TRUE:PETSC_FALSE;
565:   return(0);
566: }

570: PetscErrorCode EPSDavidsonGetConstantCorrectionTolerance_Davidson(EPS eps,PetscBool *constant)
571: {
572:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

575:   *constant = !data->dynamic?PETSC_TRUE:PETSC_FALSE;
576:   return(0);
577: }


582: PetscErrorCode EPSDavidsonSetWindowSizes_Davidson(EPS eps,PetscInt pwindow,PetscInt qwindow)
583: {
584:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

587:   if(pwindow == PETSC_DEFAULT || pwindow == PETSC_DECIDE) pwindow = 0;
588:   if(pwindow < 0)
589:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid pwindow value");
590:   if(qwindow == PETSC_DEFAULT || qwindow == PETSC_DECIDE) qwindow = 0;
591:   if(qwindow < 0)
592:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid qwindow value");
593:   data->cX_in_proj = qwindow;
594:   data->cX_in_impr = pwindow;
595:   return(0);
596: }

600: PetscErrorCode EPSDavidsonGetWindowSizes_Davidson(EPS eps,PetscInt *pwindow,PetscInt *qwindow)
601: {
602:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

605:   *pwindow = data->cX_in_impr;
606:   *qwindow = data->cX_in_proj;
607:   return(0);
608: }


613: /*
614:   EPSComputeVectors_Davidson - Compute eigenvectors from the vectors
615:   provided by the eigensolver. This version is intended for solvers 
616:   that provide Schur vectors from the QZ decompositon. Given the partial
617:   Schur decomposition OP*V=V*T, the following steps are performed:
618:       1) compute eigenvectors of (S,T): S*Z=T*Z*D
619:       2) compute eigenvectors of OP: X=V*Z
620:   If left eigenvectors are required then also do Z'*Tl=D*Z', Y=W*Z
621:  */
622: PetscErrorCode EPSComputeVectors_Davidson(EPS eps)
623: {
625:   EPS_DAVIDSON   *data = (EPS_DAVIDSON*)eps->data;
626:   dvdDashboard   *d = &data->ddb;
627:   PetscScalar    *pX,*auxS;
628:   PetscInt       size_auxS;


632:   if (d->cS) {
633:     /* Compute the eigenvectors associated to (cS, cT) */
634:     PetscMalloc(sizeof(PetscScalar)*d->nconv*d->nconv,&pX);
635:     size_auxS = 11*d->nconv + 4*d->nconv*d->nconv;
636:     PetscMalloc(sizeof(PetscScalar)*size_auxS,&auxS);
637:     dvd_compute_eigenvectors(d->nconv,d->cS,d->ldcS,d->cT,d->ldcT,
638:                                     pX,d->nconv,PETSC_NULL,0,auxS,
639:                                     size_auxS,PETSC_FALSE);
640: 
641:     /* pX[i] <- pX[i] / ||pX[i]|| */
642:     SlepcDenseNorm(pX,d->nconv,d->nconv,d->nconv,d->ceigi);
643: 
644:     /* V <- cX * pX */
645:     SlepcUpdateVectorsZ(eps->V,0.0,1.0,d->cX,d->size_cX,pX,
646:                                d->nconv,d->nconv,d->nconv);
647: 
648:     PetscFree(pX);
649:     PetscFree(auxS);
650:   }

652:   eps->evecsavailable = PETSC_TRUE;
653:   return(0);
654: }