Actual source code: davidson.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:    Skeleton of Davidson solver. Actual solvers are GD and JD.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

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

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

 24:  #include davidson.h

 26: PetscErrorCode EPSView_XD(EPS eps,PetscViewer viewer);

 28: typedef struct {
 29:   /**** Solver options ****/
 30:   PetscInt  blocksize;     /* block size */
 31:   PetscInt  initialsize;   /* initial size of V */
 32:   PetscInt  minv;          /* size of V after restarting */
 33:   PetscInt  plusk;         /* keep plusk eigenvectors from the last iteration */
 34:   PetscBool ipB;           /* true if B-ortho is used */
 35:   PetscInt  method;        /* method for improving the approximate solution */
 36:   PetscReal fix;           /* the fix parameter */
 37:   PetscBool krylovstart;   /* true if the starting subspace is a Krylov basis */
 38:   PetscBool dynamic;       /* true if dynamic stopping criterion is used */
 39:   PetscInt  cX_in_proj;    /* converged vectors in the projected problem */
 40:   PetscInt  cX_in_impr;    /* converged vectors in the projector */
 41:   Method_t  scheme;        /* method employed: GD, JD or GD2 */

 43:   /**** Solver data ****/
 44:   dvdDashboard ddb;
 45: } EPS_DAVIDSON;

 49: PetscErrorCode EPSCreate_XD(EPS eps)
 50: {
 52:   EPS_DAVIDSON   *data;

 55:   eps->ops->solve          = EPSSolve_XD;
 56:   eps->ops->setup          = EPSSetUp_XD;
 57:   eps->ops->reset          = EPSReset_XD;
 58:   eps->ops->backtransform  = EPSBackTransform_Default;
 59:   eps->ops->computevectors = EPSComputeVectors_XD;
 60:   eps->ops->view           = EPSView_XD;

 62:   PetscNewLog(eps,&data);
 63:   eps->data = (void*)data;
 64:   PetscMemzero(&data->ddb,sizeof(dvdDashboard));

 66:   /* Set default values */
 67:   EPSXDSetKrylovStart_XD(eps,PETSC_FALSE);
 68:   EPSXDSetBlockSize_XD(eps,1);
 69:   EPSXDSetRestart_XD(eps,6,0);
 70:   EPSXDSetInitialSize_XD(eps,5);
 71:   EPSJDSetFix_JD(eps,0.01);
 72:   EPSXDSetBOrth_XD(eps,PETSC_TRUE);
 73:   EPSJDSetConstCorrectionTol_JD(eps,PETSC_TRUE);
 74:   EPSXDSetWindowSizes_XD(eps,0,0);
 75:   return(0);
 76: }

 80: PetscErrorCode EPSSetUp_XD(EPS eps)
 81: {
 83:   EPS_DAVIDSON   *data = (EPS_DAVIDSON*)eps->data;
 84:   dvdDashboard   *dvd = &data->ddb;
 85:   dvdBlackboard  b;
 86:   PetscInt       min_size_V,plusk,bs,initv,i,cX_in_proj,cX_in_impr,nmat;
 87:   Mat            A,B;
 88:   KSP            ksp;
 89:   PetscBool      t,ipB,ispositive,dynamic;
 90:   HarmType_t     harm;
 91:   InitType_t     init;
 92:   PetscReal      fix;
 93:   PetscScalar    target;

 96:   /* Setup EPS options and get the problem specification */
 97:   EPSXDGetBlockSize_XD(eps,&bs);
 98:   if (bs <= 0) bs = 1;
 99:   if (eps->ncv) {
100:     if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of ncv must be at least nev");
101:   } else if (eps->mpd) eps->ncv = eps->mpd + eps->nev + bs;
102:   else if (eps->nev<500) eps->ncv = PetscMin(eps->n-bs,PetscMax(2*eps->nev,eps->nev+15))+bs;
103:   else eps->ncv = PetscMin(eps->n-bs,eps->nev+500)+bs;
104:   if (!eps->mpd) eps->mpd = eps->ncv;
105:   if (eps->mpd > eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd has to be less or equal than ncv");
106:   if (eps->mpd < 2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd has to be greater than 2");
107:   if (!eps->max_it) eps->max_it = PetscMax(100*eps->ncv,2*eps->n);
108:   if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
109:   if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Wrong value of eps->which");
110:   if (!(eps->nev + bs <= eps->ncv)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The ncv has to be greater than nev plus blocksize");
111:   if (eps->trueres) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"-eps_true_residual is temporally disable in this solver.");

113:   EPSXDGetRestart_XD(eps,&min_size_V,&plusk);
114:   if (!min_size_V) min_size_V = PetscMin(PetscMax(bs,5),eps->mpd/2);
115:   if (!(min_size_V+bs <= eps->mpd)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of minv must be less than mpd minus blocksize");
116:   EPSXDGetInitialSize_XD(eps,&initv);
117:   if (eps->mpd < initv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The initv has to be less or equal than mpd");

119:   /* Set STPrecond as the default ST */
120:   if (!((PetscObject)eps->st)->type_name) {
121:     STSetType(eps->st,STPRECOND);
122:   }
123:   STPrecondSetKSPHasMat(eps->st,PETSC_FALSE);

125:   /* Change the default sigma to inf if necessary */
126:   if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL || eps->which == EPS_LARGEST_IMAGINARY) {
127:     STSetDefaultShift(eps->st,PETSC_MAX_REAL);
128:   }

130:   /* Davidson solvers only support STPRECOND */
131:   STSetUp(eps->st);
132:   PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&t);
133:   if (!t) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"%s only works with precond spectral transformation",((PetscObject)eps)->type_name);

135:   /* Setup problem specification in dvd */
136:   STGetNumMatrices(eps->st,&nmat);
137:   STGetOperators(eps->st,0,&A);
138:   if (nmat>1) { STGetOperators(eps->st,1,&B); }
139:   EPSReset_XD(eps);
140:   PetscMemzero(dvd,sizeof(dvdDashboard));
141:   dvd->A = A; dvd->B = eps->isgeneralized? B: NULL;
142:   ispositive = eps->ispositive;
143:   dvd->sA = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN: 0) | ((ispositive && !eps->isgeneralized) ? DVD_MAT_POS_DEF: 0);
144:   /* Asume -eps_hermitian means hermitian-definite in generalized problems */
145:   if (!ispositive && !eps->isgeneralized && eps->ishermitian) ispositive = PETSC_TRUE;
146:   if (!eps->isgeneralized) dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY | DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
147:   else dvd->sB = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN: 0) | (ispositive? DVD_MAT_POS_DEF: 0);
148:   ipB = (dvd->B && data->ipB && DVD_IS(dvd->sB,DVD_MAT_HERMITIAN))?PETSC_TRUE:PETSC_FALSE;
149:   if (data->ipB && !ipB) data->ipB = PETSC_FALSE;
150:   dvd->correctXnorm = ipB;
151:   dvd->sEP = ((!eps->isgeneralized || (eps->isgeneralized && ipB))? DVD_EP_STD: 0) | (ispositive? DVD_EP_HERMITIAN: 0) | ((eps->problem_type == EPS_GHIEP && ipB) ? DVD_EP_INDEFINITE : 0);
152:   dvd->nev        = eps->nev;
153:   dvd->which      = eps->which;
154:   dvd->withTarget = PETSC_TRUE;
155:   switch (eps->which) {
156:     case EPS_TARGET_MAGNITUDE:
157:     case EPS_TARGET_IMAGINARY:
158:       dvd->target[0] = target = eps->target;
159:       dvd->target[1] = 1.0;
160:       break;
161:     case EPS_TARGET_REAL:
162:       dvd->target[0] = PetscRealPart(target = eps->target);
163:       dvd->target[1] = 1.0;
164:       break;
165:     case EPS_LARGEST_REAL:
166:     case EPS_LARGEST_MAGNITUDE:
167:     case EPS_LARGEST_IMAGINARY: /* TODO: think about this case */
168:       dvd->target[0] = 1.0;
169:       dvd->target[1] = target = 0.0;
170:       break;
171:     case EPS_SMALLEST_MAGNITUDE:
172:     case EPS_SMALLEST_REAL:
173:     case EPS_SMALLEST_IMAGINARY: /* TODO: think about this case */
174:       dvd->target[0] = target = 0.0;
175:       dvd->target[1] = 1.0;
176:       break;
177:     case EPS_WHICH_USER:
178:       STGetShift(eps->st,&target);
179:       dvd->target[0] = target;
180:       dvd->target[1] = 1.0;
181:       break;
182:     case EPS_ALL:
183:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported option: which == EPS_ALL");
184:       break;
185:     default:
186:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported value of option 'which'");
187:   }
188:   dvd->tol = (eps->tol==PETSC_DEFAULT)? SLEPC_DEFAULT_TOL: eps->tol;
189:   dvd->eps = eps;

191:   /* Setup the extraction technique */
192:   if (!eps->extraction) {
193:     if (ipB || ispositive) eps->extraction = EPS_RITZ;
194:     else {
195:       switch (eps->which) {
196:         case EPS_TARGET_REAL:
197:         case EPS_TARGET_MAGNITUDE:
198:         case EPS_TARGET_IMAGINARY:
199:         case EPS_SMALLEST_MAGNITUDE:
200:         case EPS_SMALLEST_REAL:
201:         case EPS_SMALLEST_IMAGINARY:
202:           eps->extraction = EPS_HARMONIC;
203:           break;
204:         case EPS_LARGEST_REAL:
205:         case EPS_LARGEST_MAGNITUDE:
206:         case EPS_LARGEST_IMAGINARY:
207:           eps->extraction = EPS_HARMONIC_LARGEST;
208:           break;
209:         default:
210:           eps->extraction = EPS_RITZ;
211:       }
212:     }
213:   }
214:   switch (eps->extraction) {
215:     case EPS_RITZ:              harm = DVD_HARM_NONE; break;
216:     case EPS_HARMONIC:          harm = DVD_HARM_RR; break;
217:     case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
218:     case EPS_HARMONIC_RIGHT:    harm = DVD_HARM_REIGS; break;
219:     case EPS_HARMONIC_LARGEST:  harm = DVD_HARM_LEIGS; break;
220:     default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
221:   }

223:   /* Setup the type of starting subspace */
224:   EPSXDGetKrylovStart_XD(eps,&t);
225:   init = (!t)? DVD_INITV_CLASSIC : DVD_INITV_KRYLOV;

227:   /* Setup the presence of converged vectors in the projected problem and the projector */
228:   EPSXDGetWindowSizes_XD(eps,&cX_in_impr,&cX_in_proj);
229:   if (cX_in_impr>0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The option pwindow is temporally disable in this solver.");
230:   if (cX_in_proj>0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The option qwindow is temporally disable in this solver.");
231:   if (min_size_V <= cX_in_proj) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"minv has to be greater than qwindow");
232:   if (bs > 1 && cX_in_impr > 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported option: pwindow > 0 and bs > 1");

234:   /* Get the fix parameter */
235:   EPSXDGetFix_XD(eps,&fix);

237:   /* Get whether the stopping criterion is used */
238:   EPSJDGetConstCorrectionTol_JD(eps,&dynamic);

240:   /* Preconfigure dvd */
241:   STGetKSP(eps->st,&ksp);
242:   dvd_schm_basic_preconf(dvd,&b,eps->mpd,min_size_V,bs,initv,PetscAbs(eps->nini),plusk,harm,ksp,init,eps->trackall,data->ipB,cX_in_proj,cX_in_impr,data->scheme);

244:   /* Allocate memory */
245:   EPSAllocateSolution(eps,0);

247:   /* Setup orthogonalization */
248:   EPS_SetInnerProduct(eps);
249:   if (!(ipB && dvd->B)) {
250:     BVSetMatrix(eps->V,NULL,PETSC_FALSE);
251:   }

253:   for (i=0;i<eps->ncv;i++) eps->perm[i] = i;

255:   /* Configure dvd for a basic GD */
256:   dvd_schm_basic_conf(dvd,&b,eps->mpd,min_size_V,bs,initv,PetscAbs(eps->nini),plusk,harm,dvd->withTarget,target,ksp,fix,init,eps->trackall,data->ipB,cX_in_proj,cX_in_impr,dynamic,data->scheme);
257:   return(0);
258: }

262: PetscErrorCode EPSSolve_XD(EPS eps)
263: {
264:   EPS_DAVIDSON   *data = (EPS_DAVIDSON*)eps->data;
265:   dvdDashboard   *d = &data->ddb;
266:   PetscInt       l,k;

270:   /* Call the starting routines */
271:   EPSDavidsonFLCall(d->startList,d);

273:   for (eps->its=0;eps->its<eps->max_it;eps->its++) {
274:     /* Initialize V, if it is needed */
275:     BVGetActiveColumns(d->eps->V,&l,&k);
276:     if (l == k) { d->initV(d); }

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

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

284:     /* Expand the subspace */
285:     d->updateV(d);

287:     /* Monitor progress */
288:     eps->nconv = d->nconv;
289:     BVGetActiveColumns(d->eps->V,&l,&k);
290:     EPSMonitor(eps,eps->its+1,eps->nconv,eps->eigr,eps->eigi,eps->errest,k);
291:   }

293:   /* Call the ending routines */
294:   EPSDavidsonFLCall(d->endList,d);

296:   if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
297:   else eps->reason = EPS_DIVERGED_ITS;
298:   return(0);
299: }

303: PetscErrorCode EPSReset_XD(EPS eps)
304: {
305:   EPS_DAVIDSON   *data = (EPS_DAVIDSON*)eps->data;
306:   dvdDashboard   *dvd = &data->ddb;

310:   /* Call step destructors and destroys the list */
311:   EPSDavidsonFLCall(dvd->destroyList,dvd);
312:   EPSDavidsonFLDestroy(&dvd->destroyList);
313:   EPSDavidsonFLDestroy(&dvd->startList);
314:   EPSDavidsonFLDestroy(&dvd->endList);
315:   return(0);
316: }

320: PetscErrorCode EPSView_XD(EPS eps,PetscViewer viewer)
321: {
323:   PetscBool      isascii,opb;
324:   PetscInt       opi,opi0;
325:   Method_t       meth;
326:   PetscBool      borth;

329:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
330:   if (isascii) {
331:     EPSXDGetMethod_XD(eps,&meth);
332:     if (meth==DVD_METH_GD2) {
333:       PetscViewerASCIIPrintf(viewer,"  Davidson: using double expansion variant (GD2)\n");
334:     }
335:     EPSXDGetBOrth_XD(eps,&borth);
336:     if (borth) {
337:       PetscViewerASCIIPrintf(viewer,"  Davidson: search subspace is B-orthogonalized\n");
338:     } else {
339:       PetscViewerASCIIPrintf(viewer,"  Davidson: search subspace is orthogonalized\n");
340:     }
341:     EPSXDGetBlockSize_XD(eps,&opi);
342:     PetscViewerASCIIPrintf(viewer,"  Davidson: block size=%D\n",opi);
343:     EPSXDGetKrylovStart_XD(eps,&opb);
344:     if (!opb) {
345:       PetscViewerASCIIPrintf(viewer,"  Davidson: type of the initial subspace: non-Krylov\n");
346:     } else {
347:       PetscViewerASCIIPrintf(viewer,"  Davidson: type of the initial subspace: Krylov\n");
348:     }
349:     EPSXDGetRestart_XD(eps,&opi,&opi0);
350:     PetscViewerASCIIPrintf(viewer,"  Davidson: size of the subspace after restarting: %D\n",opi);
351:     PetscViewerASCIIPrintf(viewer,"  Davidson: number of vectors after restarting from the previous iteration: %D\n",opi0);
352:   }
353:   return(0);
354: }

358: PetscErrorCode EPSXDSetKrylovStart_XD(EPS eps,PetscBool krylovstart)
359: {
360:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

363:   data->krylovstart = krylovstart;
364:   return(0);
365: }

369: PetscErrorCode EPSXDGetKrylovStart_XD(EPS eps,PetscBool *krylovstart)
370: {
371:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

374:   *krylovstart = data->krylovstart;
375:   return(0);
376: }

380: PetscErrorCode EPSXDSetBlockSize_XD(EPS eps,PetscInt blocksize)
381: {
382:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

385:   if (blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
386:   if (blocksize <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value");
387:   data->blocksize = blocksize;
388:   return(0);
389: }

393: PetscErrorCode EPSXDGetBlockSize_XD(EPS eps,PetscInt *blocksize)
394: {
395:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

398:   *blocksize = data->blocksize;
399:   return(0);
400: }

404: PetscErrorCode EPSXDSetRestart_XD(EPS eps,PetscInt minv,PetscInt plusk)
405: {
406:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

409:   if (minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 5;
410:   if (minv <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value");
411:   if (plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = 5;
412:   if (plusk < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value");
413:   data->minv = minv;
414:   data->plusk = plusk;
415:   return(0);
416: }

420: PetscErrorCode EPSXDGetRestart_XD(EPS eps,PetscInt *minv,PetscInt *plusk)
421: {
422:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

425:   if (minv) *minv = data->minv;
426:   if (plusk) *plusk = data->plusk;
427:   return(0);
428: }

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

437:   *initialsize = data->initialsize;
438:   return(0);
439: }

443: PetscErrorCode EPSXDSetInitialSize_XD(EPS eps,PetscInt initialsize)
444: {
445:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

448:   if (initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 5;
449:   if (initialsize <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value");
450:   data->initialsize = initialsize;
451:   return(0);
452: }

456: PetscErrorCode EPSXDGetFix_XD(EPS eps,PetscReal *fix)
457: {
458:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

461:   *fix = data->fix;
462:   return(0);
463: }

467: PetscErrorCode EPSJDSetFix_JD(EPS eps,PetscReal fix)
468: {
469:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

472:   if (fix == PETSC_DEFAULT || fix == PETSC_DECIDE) fix = 0.01;
473:   if (fix < 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
474:   data->fix = fix;
475:   return(0);
476: }

480: PetscErrorCode EPSXDSetBOrth_XD(EPS eps,PetscBool borth)
481: {
482:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

485:   data->ipB = borth;
486:   return(0);
487: }

491: PetscErrorCode EPSXDGetBOrth_XD(EPS eps,PetscBool *borth)
492: {
493:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

496:   *borth = data->ipB;
497:   return(0);
498: }

502: PetscErrorCode EPSJDSetConstCorrectionTol_JD(EPS eps,PetscBool constant)
503: {
504:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

507:   data->dynamic = (!constant)? PETSC_TRUE: PETSC_FALSE;
508:   return(0);
509: }

513: PetscErrorCode EPSJDGetConstCorrectionTol_JD(EPS eps,PetscBool *constant)
514: {
515:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

518:   *constant = (!data->dynamic)? PETSC_TRUE: PETSC_FALSE;
519:   return(0);
520: }

524: PetscErrorCode EPSXDSetWindowSizes_XD(EPS eps,PetscInt pwindow,PetscInt qwindow)
525: {
526:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

529:   if (pwindow == PETSC_DEFAULT || pwindow == PETSC_DECIDE) pwindow = 0;
530:   if (pwindow < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid pwindow value");
531:   if (qwindow == PETSC_DEFAULT || qwindow == PETSC_DECIDE) qwindow = 0;
532:   if (qwindow < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid qwindow value");
533:   data->cX_in_proj = qwindow;
534:   data->cX_in_impr = pwindow;
535:   return(0);
536: }

540: PetscErrorCode EPSXDGetWindowSizes_XD(EPS eps,PetscInt *pwindow,PetscInt *qwindow)
541: {
542:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

545:   if (pwindow) *pwindow = data->cX_in_impr;
546:   if (qwindow) *qwindow = data->cX_in_proj;
547:   return(0);
548: }

552: PetscErrorCode EPSXDSetMethod(EPS eps,Method_t method)
553: {
554:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

557:   data->scheme = method;
558:   return(0);
559: }

563: PetscErrorCode EPSXDGetMethod_XD(EPS eps,Method_t *method)
564: {
565:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

568:   *method = data->scheme;
569:   return(0);
570: }

574: /*
575:   EPSComputeVectors_XD - Compute eigenvectors from the vectors
576:   provided by the eigensolver. This version is intended for solvers
577:   that provide Schur vectors from the QZ decomposition. Given the partial
578:   Schur decomposition OP*V=V*T, the following steps are performed:
579:       1) compute eigenvectors of (S,T): S*Z=T*Z*D
580:       2) compute eigenvectors of OP: X=V*Z
581:  */
582: PetscErrorCode EPSComputeVectors_XD(EPS eps)
583: {
585:   Mat            X;
586:   PetscBool      symm;

589:   PetscObjectTypeCompareAny((PetscObject)eps->ds,&symm,DSHEP,"");
590:   if (symm) return(0);
591:   DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
592:   DSNormalize(eps->ds,DS_MAT_X,-1);

594:   /* V <- V * X */
595:   DSGetMat(eps->ds,DS_MAT_X,&X);
596:   BVSetActiveColumns(eps->V,0,eps->nconv);
597:   BVMultInPlace(eps->V,X,0,eps->nconv);
598:   DSRestoreMat(eps->ds,DS_MAT_X,&X);
599:   return(0);
600: }