Actual source code: davidson.c
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-2013, 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: EPSOrthType 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;
46: /**** Things to destroy ****/
47: PetscScalar *wS;
48: Vec *wV;
49: PetscInt size_wV;
50: } EPS_DAVIDSON;
54: PetscErrorCode EPSCreate_XD(EPS eps)
55: {
57: EPS_DAVIDSON *data;
60: eps->st->ops->getbilinearform = STGetBilinearForm_Default;
61: eps->ops->solve = EPSSolve_XD;
62: eps->ops->setup = EPSSetUp_XD;
63: eps->ops->reset = EPSReset_XD;
64: eps->ops->backtransform = EPSBackTransform_Default;
65: eps->ops->computevectors = EPSComputeVectors_XD;
66: eps->ops->view = EPSView_XD;
68: PetscNewLog(eps,EPS_DAVIDSON,&data);
69: eps->data = data;
70: data->wS = NULL;
71: data->wV = NULL;
72: data->size_wV = 0;
73: PetscMemzero(&data->ddb,sizeof(dvdDashboard));
75: /* Set default values */
76: EPSXDSetKrylovStart_XD(eps,PETSC_FALSE);
77: EPSXDSetBlockSize_XD(eps,1);
78: EPSXDSetRestart_XD(eps,6,0);
79: EPSXDSetInitialSize_XD(eps,5);
80: EPSJDSetFix_JD(eps,0.01);
81: EPSXDSetBOrth_XD(eps,EPS_ORTH_B);
82: EPSJDSetConstCorrectionTol_JD(eps,PETSC_TRUE);
83: EPSXDSetWindowSizes_XD(eps,0,0);
84: return(0);
85: }
89: PetscErrorCode EPSSetUp_XD(EPS eps)
90: {
92: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
93: dvdDashboard *dvd = &data->ddb;
94: dvdBlackboard b;
95: PetscInt nvecs,nscalars,min_size_V,plusk,bs,initv,i,cX_in_proj,cX_in_impr,nmat;
96: Mat A,B;
97: KSP ksp;
98: PetscBool t,ipB,ispositive,dynamic;
99: HarmType_t harm;
100: InitType_t init;
101: PetscReal fix;
102: PetscScalar target;
105: /* Setup EPS options and get the problem specification */
106: EPSXDGetBlockSize_XD(eps,&bs);
107: if (bs <= 0) bs = 1;
108: if (eps->ncv) {
109: if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of ncv must be at least nev");
110: } else if (eps->mpd) eps->ncv = eps->mpd + eps->nev + bs;
111: else if (eps->nev<500) eps->ncv = PetscMin(eps->n-bs,PetscMax(2*eps->nev,eps->nev+15))+bs;
112: else eps->ncv = PetscMin(eps->n-bs,eps->nev+500)+bs;
113: if (!eps->mpd) eps->mpd = eps->ncv;
114: if (eps->mpd > eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd has to be less or equal than ncv");
115: if (eps->mpd < 2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd has to be greater than 2");
116: if (!eps->max_it) eps->max_it = PetscMax(100*eps->ncv,2*eps->n);
117: if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
118: 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");
119: if (!(eps->nev + bs <= eps->ncv)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The ncv has to be greater than nev plus blocksize");
121: EPSXDGetRestart_XD(eps,&min_size_V,&plusk);
122: if (!min_size_V) min_size_V = PetscMin(PetscMax(bs,5),eps->mpd/2);
123: 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");
124: EPSXDGetInitialSize_XD(eps,&initv);
125: if (eps->mpd < initv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The initv has to be less or equal than mpd");
127: /* Davidson solvers do not support left eigenvectors */
128: if (eps->leftvecs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Left vectors not supported in this solver");
130: /* Set STPrecond as the default ST */
131: if (!((PetscObject)eps->st)->type_name) {
132: STSetType(eps->st,STPRECOND);
133: }
134: STPrecondSetKSPHasMat(eps->st,PETSC_FALSE);
136: /* Change the default sigma to inf if necessary */
137: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
138: eps->which == EPS_LARGEST_IMAGINARY) {
139: STSetDefaultShift(eps->st,PETSC_MAX_REAL);
140: }
142: /* Davidson solvers only support STPRECOND */
143: STSetUp(eps->st);
144: PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&t);
145: if (!t) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"%s only works with precond spectral transformation",
146: ((PetscObject)eps)->type_name);
148: /* Setup problem specification in dvd */
149: STGetNumMatrices(eps->st,&nmat);
150: STGetOperators(eps->st,0,&A);
151: if (nmat>1) { STGetOperators(eps->st,1,&B); }
152: EPSReset_XD(eps);
153: PetscMemzero(dvd,sizeof(dvdDashboard));
154: dvd->A = A; dvd->B = eps->isgeneralized? B : NULL;
155: ispositive = eps->ispositive;
156: dvd->sA = DVD_MAT_IMPLICIT |
157: (eps->ishermitian? DVD_MAT_HERMITIAN : 0) |
158: ((ispositive && !eps->isgeneralized) ? DVD_MAT_POS_DEF : 0);
159: /* Asume -eps_hermitian means hermitian-definite in generalized problems */
160: if (!ispositive && !eps->isgeneralized && eps->ishermitian) ispositive = PETSC_TRUE;
161: if (!eps->isgeneralized) dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY | DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
162: else dvd->sB = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN : 0) | (ispositive? DVD_MAT_POS_DEF : 0);
163: ipB = (dvd->B && data->ipB != EPS_ORTH_I && DVD_IS(dvd->sB,DVD_MAT_HERMITIAN))?PETSC_TRUE:PETSC_FALSE;
164: if (data->ipB != EPS_ORTH_I && !ipB) data->ipB = EPS_ORTH_I;
165: dvd->correctXnorm = ipB;
166: dvd->sEP = ((!eps->isgeneralized || (eps->isgeneralized && ipB))? DVD_EP_STD : 0) |
167: (ispositive? DVD_EP_HERMITIAN : 0) |
168: ((eps->problem_type == EPS_GHIEP && ipB) ? DVD_EP_INDEFINITE : 0);
169: dvd->nev = eps->nev;
170: dvd->which = eps->which;
171: dvd->withTarget = PETSC_TRUE;
172: switch (eps->which) {
173: case EPS_TARGET_MAGNITUDE:
174: case EPS_TARGET_IMAGINARY:
175: dvd->target[0] = target = eps->target; dvd->target[1] = 1.0;
176: break;
177: case EPS_TARGET_REAL:
178: dvd->target[0] = PetscRealPart(target = eps->target); dvd->target[1] = 1.0;
179: break;
180: case EPS_LARGEST_REAL:
181: case EPS_LARGEST_MAGNITUDE:
182: case EPS_LARGEST_IMAGINARY: /* TODO: think about this case */
183: dvd->target[0] = 1.0; dvd->target[1] = target = 0.0;
184: break;
185: case EPS_SMALLEST_MAGNITUDE:
186: case EPS_SMALLEST_REAL:
187: case EPS_SMALLEST_IMAGINARY: /* TODO: think about this case */
188: dvd->target[0] = target = 0.0; dvd->target[1] = 1.0;
189: break;
190: case EPS_WHICH_USER:
191: STGetShift(eps->st,&target);
192: dvd->target[0] = target; dvd->target[1] = 1.0;
193: break;
194: case EPS_ALL:
195: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported option: which == EPS_ALL");
196: break;
197: default:
198: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported value of option 'which'");
199: }
200: dvd->tol = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
201: dvd->eps = eps;
203: /* Setup the extraction technique */
204: if (!eps->extraction) {
205: if (ipB || ispositive) eps->extraction = EPS_RITZ;
206: else {
207: switch (eps->which) {
208: case EPS_TARGET_REAL:
209: case EPS_TARGET_MAGNITUDE:
210: case EPS_TARGET_IMAGINARY:
211: case EPS_SMALLEST_MAGNITUDE:
212: case EPS_SMALLEST_REAL:
213: case EPS_SMALLEST_IMAGINARY:
214: eps->extraction = EPS_HARMONIC;
215: break;
216: case EPS_LARGEST_REAL:
217: case EPS_LARGEST_MAGNITUDE:
218: case EPS_LARGEST_IMAGINARY:
219: eps->extraction = EPS_HARMONIC_LARGEST;
220: break;
221: default:
222: eps->extraction = EPS_RITZ;
223: }
224: }
225: }
226: switch (eps->extraction) {
227: case EPS_RITZ: harm = DVD_HARM_NONE; break;
228: case EPS_HARMONIC: harm = DVD_HARM_RR; break;
229: case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
230: case EPS_HARMONIC_RIGHT: harm = DVD_HARM_REIGS; break;
231: case EPS_HARMONIC_LARGEST: harm = DVD_HARM_LEIGS; break;
232: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
233: }
235: /* Setup the type of starting subspace */
236: EPSXDGetKrylovStart_XD(eps,&t);
237: init = (!t)? DVD_INITV_CLASSIC : DVD_INITV_KRYLOV;
239: /* Setup the presence of converged vectors in the projected problem and in the projector */
240: EPSXDGetWindowSizes_XD(eps,&cX_in_impr,&cX_in_proj);
241: if (min_size_V <= cX_in_proj) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"minv has to be greater than qwindow");
242: if (bs > 1 && cX_in_impr > 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported option: pwindow > 0 and bs > 1");
244: /* Setup IP */
245: if (ipB && dvd->B) {
246: IPSetMatrix(eps->ip,dvd->B);
247: } else {
248: IPSetMatrix(eps->ip,NULL);
249: }
251: /* Get the fix parameter */
252: EPSXDGetFix_XD(eps,&fix);
254: /* Get whether the stopping criterion is used */
255: EPSJDGetConstCorrectionTol_JD(eps,&dynamic);
257: /* Orthonormalize the deflation space */
258: dvd_orthV(eps->ip,NULL,0,NULL,0,eps->defl,0,PetscAbs(eps->nds),NULL,eps->rand);
260: /* Preconfigure dvd */
261: STGetKSP(eps->st,&ksp);
262: dvd_schm_basic_preconf(dvd,&b,eps->mpd,min_size_V,bs,
263: initv,
264: PetscAbs(eps->nini),
265: plusk,harm,
266: ksp,init,eps->trackall,
267: data->ipB,cX_in_proj,cX_in_impr,
268: data->scheme);
270: /* Allocate memory */
271: nvecs = b.max_size_auxV + b.own_vecs;
272: nscalars = b.own_scalars + b.max_size_auxS;
273: PetscMalloc(nscalars*sizeof(PetscScalar),&data->wS);
274: PetscLogObjectMemory(eps,nscalars*sizeof(PetscScalar));
275: VecDuplicateVecs(eps->t,nvecs,&data->wV);
276: PetscLogObjectParents(eps,nvecs,data->wV);
277: data->size_wV = nvecs;
278: b.free_vecs = data->wV;
279: b.free_scalars = data->wS;
280: dvd->auxV = data->wV + b.own_vecs;
281: dvd->auxS = b.free_scalars + b.own_scalars;
282: dvd->size_auxV = b.max_size_auxV;
283: dvd->size_auxS = b.max_size_auxS;
285: eps->errest_left = NULL;
286: PetscMalloc(eps->ncv*sizeof(PetscInt),&eps->perm);
287: PetscLogObjectMemory(eps,eps->ncv*sizeof(PetscInt));
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: data->ipB,cX_in_proj,cX_in_impr,dynamic,
298: data->scheme);
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;
305: return(0);
306: }
310: PetscErrorCode EPSSolve_XD(EPS eps)
311: {
312: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
313: dvdDashboard *d = &data->ddb;
317: /* Call the starting routines */
318: DVD_FL_CALL(d->startList,d);
320: for (eps->its=0;eps->its<eps->max_it;eps->its++) {
321: /* Initialize V, if it is needed */
322: if (d->size_V == 0) { d->initV(d); }
324: /* Find the best approximated eigenpairs in V, X */
325: d->calcPairs(d);
327: /* Test for convergence */
328: if (eps->nconv >= eps->nev) break;
330: /* Expand the subspace */
331: d->updateV(d);
333: /* Monitor progress */
334: eps->nconv = d->nconv;
335: EPSMonitor(eps,eps->its+1,eps->nconv,eps->eigr,eps->eigi,eps->errest,d->size_V+d->size_cX);
336: }
338: /* Call the ending routines */
339: DVD_FL_CALL(d->endList,d);
341: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
342: else eps->reason = EPS_DIVERGED_ITS;
343: return(0);
344: }
348: PetscErrorCode EPSReset_XD(EPS eps)
349: {
350: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
351: dvdDashboard *dvd = &data->ddb;
355: /* Call step destructors and destroys the list */
356: DVD_FL_CALL(dvd->destroyList,dvd);
357: DVD_FL_DEL(dvd->destroyList);
358: DVD_FL_DEL(dvd->startList);
359: DVD_FL_DEL(dvd->endList);
361: if (data->size_wV > 0) {
362: VecDestroyVecs(data->size_wV,&data->wV);
363: }
364: PetscFree(data->wS);
365: PetscFree(eps->perm);
366: return(0);
367: }
371: PetscErrorCode EPSView_XD(EPS eps,PetscViewer viewer)
372: {
374: PetscBool isascii,opb;
375: PetscInt opi,opi0;
376: Method_t meth;
377: EPSOrthType borth;
380: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
381: if (isascii) {
382: EPSXDGetMethod_XD(eps,&meth);
383: if (meth==DVD_METH_GD2) {
384: PetscViewerASCIIPrintf(viewer," Davidson: using double expansion variant (GD2)\n");
385: }
386: EPSXDGetBOrth_XD(eps,&borth);
387: switch (borth) {
388: case EPS_ORTH_I:
389: PetscViewerASCIIPrintf(viewer," Davidson: search subspace is orthogonalized\n");
390: break;
391: case EPS_ORTH_B:
392: PetscViewerASCIIPrintf(viewer," Davidson: search subspace is B-orthogonalized\n");
393: break;
394: case EPS_ORTH_BOPT:
395: PetscViewerASCIIPrintf(viewer," Davidson: search subspace is B-orthogonalized with an optimized method\n");
396: break;
397: }
398: EPSXDGetBlockSize_XD(eps,&opi);
399: PetscViewerASCIIPrintf(viewer," Davidson: block size=%D\n",opi);
400: EPSXDGetKrylovStart_XD(eps,&opb);
401: if (!opb) {
402: PetscViewerASCIIPrintf(viewer," Davidson: type of the initial subspace: non-Krylov\n");
403: } else {
404: PetscViewerASCIIPrintf(viewer," Davidson: type of the initial subspace: Krylov\n");
405: }
406: EPSXDGetRestart_XD(eps,&opi,&opi0);
407: PetscViewerASCIIPrintf(viewer," Davidson: size of the subspace after restarting: %D\n",opi);
408: PetscViewerASCIIPrintf(viewer," Davidson: number of vectors after restarting from the previous iteration: %D\n",opi0);
409: }
410: return(0);
411: }
415: PetscErrorCode EPSXDSetKrylovStart_XD(EPS eps,PetscBool krylovstart)
416: {
417: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
420: data->krylovstart = krylovstart;
421: return(0);
422: }
426: PetscErrorCode EPSXDGetKrylovStart_XD(EPS eps,PetscBool *krylovstart)
427: {
428: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
431: *krylovstart = data->krylovstart;
432: return(0);
433: }
437: PetscErrorCode EPSXDSetBlockSize_XD(EPS eps,PetscInt blocksize)
438: {
439: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
442: if (blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
443: if (blocksize <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value");
444: data->blocksize = blocksize;
445: return(0);
446: }
450: PetscErrorCode EPSXDGetBlockSize_XD(EPS eps,PetscInt *blocksize)
451: {
452: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
455: *blocksize = data->blocksize;
456: return(0);
457: }
461: PetscErrorCode EPSXDSetRestart_XD(EPS eps,PetscInt minv,PetscInt plusk)
462: {
463: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
466: if (minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 5;
467: if (minv <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value");
468: if (plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = 5;
469: if (plusk < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value");
470: data->minv = minv;
471: data->plusk = plusk;
472: return(0);
473: }
477: PetscErrorCode EPSXDGetRestart_XD(EPS eps,PetscInt *minv,PetscInt *plusk)
478: {
479: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
482: if (minv) *minv = data->minv;
483: if (plusk) *plusk = data->plusk;
484: return(0);
485: }
489: PetscErrorCode EPSXDGetInitialSize_XD(EPS eps,PetscInt *initialsize)
490: {
491: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
494: *initialsize = data->initialsize;
495: return(0);
496: }
500: PetscErrorCode EPSXDSetInitialSize_XD(EPS eps,PetscInt initialsize)
501: {
502: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
505: if (initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 5;
506: if (initialsize <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value");
507: data->initialsize = initialsize;
508: return(0);
509: }
513: PetscErrorCode EPSXDGetFix_XD(EPS eps,PetscReal *fix)
514: {
515: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
518: *fix = data->fix;
519: return(0);
520: }
524: PetscErrorCode EPSJDSetFix_JD(EPS eps,PetscReal fix)
525: {
526: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
529: if (fix == PETSC_DEFAULT || fix == PETSC_DECIDE) fix = 0.01;
530: if (fix < 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
531: data->fix = fix;
532: return(0);
533: }
537: PetscErrorCode EPSXDSetBOrth_XD(EPS eps,EPSOrthType borth)
538: {
539: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
542: data->ipB = borth;
543: return(0);
544: }
548: PetscErrorCode EPSXDGetBOrth_XD(EPS eps,EPSOrthType *borth)
549: {
550: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
553: *borth = data->ipB;
554: return(0);
555: }
559: PetscErrorCode EPSJDSetConstCorrectionTol_JD(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 EPSJDGetConstCorrectionTol_JD(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: }
581: PetscErrorCode EPSXDSetWindowSizes_XD(EPS eps,PetscInt pwindow,PetscInt qwindow)
582: {
583: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
586: if (pwindow == PETSC_DEFAULT || pwindow == PETSC_DECIDE) pwindow = 0;
587: if (pwindow < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid pwindow value");
588: if (qwindow == PETSC_DEFAULT || qwindow == PETSC_DECIDE) qwindow = 0;
589: if (qwindow < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid qwindow value");
590: data->cX_in_proj = qwindow;
591: data->cX_in_impr = pwindow;
592: return(0);
593: }
597: PetscErrorCode EPSXDGetWindowSizes_XD(EPS eps,PetscInt *pwindow,PetscInt *qwindow)
598: {
599: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
602: if (pwindow) *pwindow = data->cX_in_impr;
603: if (qwindow) *qwindow = data->cX_in_proj;
604: return(0);
605: }
609: PetscErrorCode EPSXDSetMethod(EPS eps,Method_t method)
610: {
611: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
614: data->scheme = method;
615: return(0);
616: }
620: PetscErrorCode EPSXDGetMethod_XD(EPS eps,Method_t *method)
621: {
622: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
625: *method = data->scheme;
626: return(0);
627: }
631: /*
632: EPSComputeVectors_XD - Compute eigenvectors from the vectors
633: provided by the eigensolver. This version is intended for solvers
634: that provide Schur vectors from the QZ decompositon. Given the partial
635: Schur decomposition OP*V=V*T, the following steps are performed:
636: 1) compute eigenvectors of (S,T): S*Z=T*Z*D
637: 2) compute eigenvectors of OP: X=V*Z
638: If left eigenvectors are required then also do Z'*T=D*Z', Y=W*Z
639: */
640: PetscErrorCode EPSComputeVectors_XD(EPS eps)
641: {
643: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
644: dvdDashboard *d = &data->ddb;
645: PetscScalar *pX,*cS,*cT;
646: PetscInt ld;
649: if (d->cS) {
650: /* Compute the eigenvectors associated to (cS, cT) */
651: DSSetDimensions(d->conv_ps,d->size_cS,0,0,0);
652: DSGetLeadingDimension(d->conv_ps,&ld);
653: DSGetArray(d->conv_ps,DS_MAT_A,&cS);
654: SlepcDenseCopyTriang(cS,0,ld,d->cS,0,d->ldcS,d->size_cS,d->size_cS);
655: DSRestoreArray(d->conv_ps,DS_MAT_A,&cS);
656: if (d->cT) {
657: DSGetArray(d->conv_ps,DS_MAT_B,&cT);
658: SlepcDenseCopyTriang(cT,0,ld,d->cT,0,d->ldcT,d->size_cS,d->size_cS);
659: DSRestoreArray(d->conv_ps,DS_MAT_B,&cT);
660: }
661: DSSetState(d->conv_ps,DS_STATE_RAW);
662: DSSolve(d->conv_ps,eps->eigr,eps->eigi);
663: DSVectors(d->conv_ps,DS_MAT_X,NULL,NULL);
664: DSNormalize(d->conv_ps,DS_MAT_X,-1);
666: /* V <- cX * pX */
667: DSGetArray(d->conv_ps,DS_MAT_X,&pX);
668: SlepcUpdateVectorsZ(eps->V,0.0,1.0,d->cX,d->size_cX,pX,ld,d->nconv,d->nconv);
669: DSRestoreArray(d->conv_ps,DS_MAT_X,&pX);
670: }
672: eps->evecsavailable = PETSC_TRUE;
673: return(0);
674: }