Actual source code: davidson.c
slepc-3.5.2 2014-10-10
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-2014, 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: initialsize, /* initial size of V */
32: minv, /* size of V after restarting */
33: 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: 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 ||
127: eps->which == EPS_LARGEST_IMAGINARY) {
128: STSetDefaultShift(eps->st,PETSC_MAX_REAL);
129: }
131: /* Davidson solvers only support STPRECOND */
132: STSetUp(eps->st);
133: PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&t);
134: if (!t) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"%s only works with precond spectral transformation",
135: ((PetscObject)eps)->type_name);
137: /* Setup problem specification in dvd */
138: STGetNumMatrices(eps->st,&nmat);
139: STGetOperators(eps->st,0,&A);
140: if (nmat>1) { STGetOperators(eps->st,1,&B); }
141: EPSReset_XD(eps);
142: PetscMemzero(dvd,sizeof(dvdDashboard));
143: dvd->A = A; dvd->B = eps->isgeneralized? B : NULL;
144: ispositive = eps->ispositive;
145: dvd->sA = DVD_MAT_IMPLICIT |
146: (eps->ishermitian? DVD_MAT_HERMITIAN : 0) |
147: ((ispositive && !eps->isgeneralized) ? DVD_MAT_POS_DEF : 0);
148: /* Asume -eps_hermitian means hermitian-definite in generalized problems */
149: if (!ispositive && !eps->isgeneralized && eps->ishermitian) ispositive = PETSC_TRUE;
150: if (!eps->isgeneralized) dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY | DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
151: else dvd->sB = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN : 0) | (ispositive? DVD_MAT_POS_DEF : 0);
152: ipB = (dvd->B && data->ipB && DVD_IS(dvd->sB,DVD_MAT_HERMITIAN))?PETSC_TRUE:PETSC_FALSE;
153: if (data->ipB && !ipB) data->ipB = PETSC_FALSE;
154: dvd->correctXnorm = ipB;
155: dvd->sEP = ((!eps->isgeneralized || (eps->isgeneralized && ipB))? DVD_EP_STD : 0) |
156: (ispositive? DVD_EP_HERMITIAN : 0) |
157: ((eps->problem_type == EPS_GHIEP && ipB) ? DVD_EP_INDEFINITE : 0);
158: dvd->nev = eps->nev;
159: dvd->which = eps->which;
160: dvd->withTarget = PETSC_TRUE;
161: switch (eps->which) {
162: case EPS_TARGET_MAGNITUDE:
163: case EPS_TARGET_IMAGINARY:
164: dvd->target[0] = target = eps->target; dvd->target[1] = 1.0;
165: break;
166: case EPS_TARGET_REAL:
167: dvd->target[0] = PetscRealPart(target = eps->target); dvd->target[1] = 1.0;
168: break;
169: case EPS_LARGEST_REAL:
170: case EPS_LARGEST_MAGNITUDE:
171: case EPS_LARGEST_IMAGINARY: /* TODO: think about this case */
172: dvd->target[0] = 1.0; dvd->target[1] = target = 0.0;
173: break;
174: case EPS_SMALLEST_MAGNITUDE:
175: case EPS_SMALLEST_REAL:
176: case EPS_SMALLEST_IMAGINARY: /* TODO: think about this case */
177: dvd->target[0] = target = 0.0; dvd->target[1] = 1.0;
178: break;
179: case EPS_WHICH_USER:
180: STGetShift(eps->st,&target);
181: dvd->target[0] = target; dvd->target[1] = 1.0;
182: break;
183: case EPS_ALL:
184: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported option: which == EPS_ALL");
185: break;
186: default:
187: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported value of option 'which'");
188: }
189: dvd->tol = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
190: dvd->eps = eps;
192: /* Setup the extraction technique */
193: if (!eps->extraction) {
194: if (ipB || ispositive) eps->extraction = EPS_RITZ;
195: else {
196: switch (eps->which) {
197: case EPS_TARGET_REAL:
198: case EPS_TARGET_MAGNITUDE:
199: case EPS_TARGET_IMAGINARY:
200: case EPS_SMALLEST_MAGNITUDE:
201: case EPS_SMALLEST_REAL:
202: case EPS_SMALLEST_IMAGINARY:
203: eps->extraction = EPS_HARMONIC;
204: break;
205: case EPS_LARGEST_REAL:
206: case EPS_LARGEST_MAGNITUDE:
207: case EPS_LARGEST_IMAGINARY:
208: eps->extraction = EPS_HARMONIC_LARGEST;
209: break;
210: default:
211: eps->extraction = EPS_RITZ;
212: }
213: }
214: }
215: switch (eps->extraction) {
216: case EPS_RITZ: harm = DVD_HARM_NONE; break;
217: case EPS_HARMONIC: harm = DVD_HARM_RR; break;
218: case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
219: case EPS_HARMONIC_RIGHT: harm = DVD_HARM_REIGS; break;
220: case EPS_HARMONIC_LARGEST: harm = DVD_HARM_LEIGS; break;
221: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
222: }
224: /* Setup the type of starting subspace */
225: EPSXDGetKrylovStart_XD(eps,&t);
226: init = (!t)? DVD_INITV_CLASSIC : DVD_INITV_KRYLOV;
228: /* Setup the presence of converged vectors in the projected problem and in the projector */
229: EPSXDGetWindowSizes_XD(eps,&cX_in_impr,&cX_in_proj);
230: if (cX_in_impr>0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The option pwindow is temporally disable in this solver.");
231: if (cX_in_proj>0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The option qwindow is temporally disable in this solver.");
232: if (min_size_V <= cX_in_proj) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"minv has to be greater than qwindow");
233: if (bs > 1 && cX_in_impr > 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported option: pwindow > 0 and bs > 1");
235: /* Get the fix parameter */
236: EPSXDGetFix_XD(eps,&fix);
238: /* Get whether the stopping criterion is used */
239: EPSJDGetConstCorrectionTol_JD(eps,&dynamic);
241: /* Preconfigure dvd */
242: STGetKSP(eps->st,&ksp);
243: dvd_schm_basic_preconf(dvd,&b,eps->mpd,min_size_V,bs,
244: initv,
245: PetscAbs(eps->nini),
246: plusk,harm,
247: ksp,init,eps->trackall,
248: data->ipB,cX_in_proj,cX_in_impr,
249: data->scheme);
251: /* Allocate memory */
252: EPSAllocateSolution(eps,0);
254: /* Setup orthogonalization */
255: EPS_SetInnerProduct(eps);
256: if (!(ipB && dvd->B)) {
257: BVSetMatrix(eps->V,NULL,PETSC_FALSE);
258: }
260: for (i=0;i<eps->ncv;i++) eps->perm[i] = i;
262: /* Configure dvd for a basic GD */
263: dvd_schm_basic_conf(dvd,&b,eps->mpd,min_size_V,bs,
264: initv,
265: PetscAbs(eps->nini),plusk,
266: harm,dvd->withTarget,
267: target,ksp,
268: fix,init,eps->trackall,
269: data->ipB,cX_in_proj,cX_in_impr,dynamic,
270: data->scheme);
271: return(0);
272: }
276: PetscErrorCode EPSSolve_XD(EPS eps)
277: {
278: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
279: dvdDashboard *d = &data->ddb;
280: PetscInt l,k;
284: /* Call the starting routines */
285: EPSDavidsonFLCall(d->startList,d);
287: for (eps->its=0;eps->its<eps->max_it;eps->its++) {
288: /* Initialize V, if it is needed */
289: BVGetActiveColumns(d->eps->V,&l,&k);
290: if (l == k) { d->initV(d); }
292: /* Find the best approximated eigenpairs in V, X */
293: d->calcPairs(d);
295: /* Test for convergence */
296: if (eps->nconv >= eps->nev) break;
298: /* Expand the subspace */
299: d->updateV(d);
301: /* Monitor progress */
302: eps->nconv = d->nconv;
303: BVGetActiveColumns(d->eps->V,&l,&k);
304: EPSMonitor(eps,eps->its+1,eps->nconv,eps->eigr,eps->eigi,eps->errest,k);
305: }
307: /* Call the ending routines */
308: EPSDavidsonFLCall(d->endList,d);
310: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
311: else eps->reason = EPS_DIVERGED_ITS;
312: return(0);
313: }
317: PetscErrorCode EPSReset_XD(EPS eps)
318: {
319: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
320: dvdDashboard *dvd = &data->ddb;
324: /* Call step destructors and destroys the list */
325: EPSDavidsonFLCall(dvd->destroyList,dvd);
326: EPSDavidsonFLDestroy(&dvd->destroyList);
327: EPSDavidsonFLDestroy(&dvd->startList);
328: EPSDavidsonFLDestroy(&dvd->endList);
329: return(0);
330: }
334: PetscErrorCode EPSView_XD(EPS eps,PetscViewer viewer)
335: {
337: PetscBool isascii,opb;
338: PetscInt opi,opi0;
339: Method_t meth;
340: PetscBool borth;
343: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
344: if (isascii) {
345: EPSXDGetMethod_XD(eps,&meth);
346: if (meth==DVD_METH_GD2) {
347: PetscViewerASCIIPrintf(viewer," Davidson: using double expansion variant (GD2)\n");
348: }
349: EPSXDGetBOrth_XD(eps,&borth);
350: if (borth) {
351: PetscViewerASCIIPrintf(viewer," Davidson: search subspace is B-orthogonalized\n");
352: } else {
353: PetscViewerASCIIPrintf(viewer," Davidson: search subspace is orthogonalized\n");
354: }
355: EPSXDGetBlockSize_XD(eps,&opi);
356: PetscViewerASCIIPrintf(viewer," Davidson: block size=%D\n",opi);
357: EPSXDGetKrylovStart_XD(eps,&opb);
358: if (!opb) {
359: PetscViewerASCIIPrintf(viewer," Davidson: type of the initial subspace: non-Krylov\n");
360: } else {
361: PetscViewerASCIIPrintf(viewer," Davidson: type of the initial subspace: Krylov\n");
362: }
363: EPSXDGetRestart_XD(eps,&opi,&opi0);
364: PetscViewerASCIIPrintf(viewer," Davidson: size of the subspace after restarting: %D\n",opi);
365: PetscViewerASCIIPrintf(viewer," Davidson: number of vectors after restarting from the previous iteration: %D\n",opi0);
366: }
367: return(0);
368: }
372: PetscErrorCode EPSXDSetKrylovStart_XD(EPS eps,PetscBool krylovstart)
373: {
374: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
377: data->krylovstart = krylovstart;
378: return(0);
379: }
383: PetscErrorCode EPSXDGetKrylovStart_XD(EPS eps,PetscBool *krylovstart)
384: {
385: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
388: *krylovstart = data->krylovstart;
389: return(0);
390: }
394: PetscErrorCode EPSXDSetBlockSize_XD(EPS eps,PetscInt blocksize)
395: {
396: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
399: if (blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
400: if (blocksize <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value");
401: data->blocksize = blocksize;
402: return(0);
403: }
407: PetscErrorCode EPSXDGetBlockSize_XD(EPS eps,PetscInt *blocksize)
408: {
409: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
412: *blocksize = data->blocksize;
413: return(0);
414: }
418: PetscErrorCode EPSXDSetRestart_XD(EPS eps,PetscInt minv,PetscInt plusk)
419: {
420: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
423: if (minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 5;
424: if (minv <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value");
425: if (plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = 5;
426: if (plusk < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value");
427: data->minv = minv;
428: data->plusk = plusk;
429: return(0);
430: }
434: PetscErrorCode EPSXDGetRestart_XD(EPS eps,PetscInt *minv,PetscInt *plusk)
435: {
436: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
439: if (minv) *minv = data->minv;
440: if (plusk) *plusk = data->plusk;
441: return(0);
442: }
446: PetscErrorCode EPSXDGetInitialSize_XD(EPS eps,PetscInt *initialsize)
447: {
448: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
451: *initialsize = data->initialsize;
452: return(0);
453: }
457: PetscErrorCode EPSXDSetInitialSize_XD(EPS eps,PetscInt initialsize)
458: {
459: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
462: if (initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 5;
463: if (initialsize <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value");
464: data->initialsize = initialsize;
465: return(0);
466: }
470: PetscErrorCode EPSXDGetFix_XD(EPS eps,PetscReal *fix)
471: {
472: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
475: *fix = data->fix;
476: return(0);
477: }
481: PetscErrorCode EPSJDSetFix_JD(EPS eps,PetscReal fix)
482: {
483: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
486: if (fix == PETSC_DEFAULT || fix == PETSC_DECIDE) fix = 0.01;
487: if (fix < 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
488: data->fix = fix;
489: return(0);
490: }
494: PetscErrorCode EPSXDSetBOrth_XD(EPS eps,PetscBool borth)
495: {
496: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
499: data->ipB = borth;
500: return(0);
501: }
505: PetscErrorCode EPSXDGetBOrth_XD(EPS eps,PetscBool *borth)
506: {
507: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
510: *borth = data->ipB;
511: return(0);
512: }
516: PetscErrorCode EPSJDSetConstCorrectionTol_JD(EPS eps,PetscBool constant)
517: {
518: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
521: data->dynamic = !constant?PETSC_TRUE:PETSC_FALSE;
522: return(0);
523: }
527: PetscErrorCode EPSJDGetConstCorrectionTol_JD(EPS eps,PetscBool *constant)
528: {
529: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
532: *constant = !data->dynamic?PETSC_TRUE:PETSC_FALSE;
533: return(0);
534: }
538: PetscErrorCode EPSXDSetWindowSizes_XD(EPS eps,PetscInt pwindow,PetscInt qwindow)
539: {
540: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
543: if (pwindow == PETSC_DEFAULT || pwindow == PETSC_DECIDE) pwindow = 0;
544: if (pwindow < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid pwindow value");
545: if (qwindow == PETSC_DEFAULT || qwindow == PETSC_DECIDE) qwindow = 0;
546: if (qwindow < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid qwindow value");
547: data->cX_in_proj = qwindow;
548: data->cX_in_impr = pwindow;
549: return(0);
550: }
554: PetscErrorCode EPSXDGetWindowSizes_XD(EPS eps,PetscInt *pwindow,PetscInt *qwindow)
555: {
556: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
559: if (pwindow) *pwindow = data->cX_in_impr;
560: if (qwindow) *qwindow = data->cX_in_proj;
561: return(0);
562: }
566: PetscErrorCode EPSXDSetMethod(EPS eps,Method_t method)
567: {
568: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
571: data->scheme = method;
572: return(0);
573: }
577: PetscErrorCode EPSXDGetMethod_XD(EPS eps,Method_t *method)
578: {
579: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
582: *method = data->scheme;
583: return(0);
584: }
588: /*
589: EPSComputeVectors_XD - Compute eigenvectors from the vectors
590: provided by the eigensolver. This version is intended for solvers
591: that provide Schur vectors from the QZ decomposition. Given the partial
592: Schur decomposition OP*V=V*T, the following steps are performed:
593: 1) compute eigenvectors of (S,T): S*Z=T*Z*D
594: 2) compute eigenvectors of OP: X=V*Z
595: */
596: PetscErrorCode EPSComputeVectors_XD(EPS eps)
597: {
599: Mat X;
600: PetscBool symm;
603: PetscObjectTypeCompareAny((PetscObject)eps->ds,&symm,DSHEP,"");
604: if (symm) return(0);
605: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
606: DSNormalize(eps->ds,DS_MAT_X,-1);
608: /* V <- V * X */
609: DSGetMat(eps->ds,DS_MAT_X,&X);
610: BVSetActiveColumns(eps->V,0,eps->nconv);
611: BVMultInPlace(eps->V,X,0,eps->nconv);
612: DSRestoreMat(eps->ds,DS_MAT_X,&X);
613: return(0);
614: }