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-2012, 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: EPSOrthType ipB; /* true if B-ortho is used */
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 */
45: Method_t scheme; /* method employed: GD, JD or GD2 */
47: /**** Solver data ****/
48: dvdDashboard ddb;
50: /**** Things to destroy ****/
51: PetscScalar *wS;
52: Vec *wV;
53: PetscInt size_wV;
54: } EPS_DAVIDSON;
58: PetscErrorCode EPSCreate_Davidson(EPS eps)
59: {
61: EPS_DAVIDSON *data;
65: eps->OP->ops->getbilinearform = STGetBilinearForm_Default;
66: eps->ops->solve = EPSSolve_Davidson;
67: eps->ops->setup = EPSSetUp_Davidson;
68: eps->ops->reset = EPSReset_Davidson;
69: eps->ops->backtransform = EPSBackTransform_Default;
70: eps->ops->computevectors = EPSComputeVectors_Davidson;
71: eps->ops->view = EPSView_Davidson;
73: PetscMalloc(sizeof(EPS_DAVIDSON),&data);
74: eps->data = data;
75: data->wS = PETSC_NULL;
76: data->wV = PETSC_NULL;
77: data->size_wV = 0;
78: PetscMemzero(&data->ddb,sizeof(dvdDashboard));
80: /* Set default values */
81: EPSDavidsonSetKrylovStart_Davidson(eps,PETSC_FALSE);
82: EPSDavidsonSetBlockSize_Davidson(eps,1);
83: EPSDavidsonSetRestart_Davidson(eps,6,0);
84: EPSDavidsonSetInitialSize_Davidson(eps,5);
85: EPSDavidsonSetFix_Davidson(eps,0.01);
86: EPSDavidsonSetBOrth_Davidson(eps,EPS_ORTH_B);
87: EPSDavidsonSetConstantCorrectionTolerance_Davidson(eps,PETSC_TRUE);
88: EPSDavidsonSetWindowSizes_Davidson(eps,0,0);
89: return(0);
90: }
94: PetscErrorCode EPSSetUp_Davidson(EPS eps)
95: {
97: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
98: dvdDashboard *dvd = &data->ddb;
99: dvdBlackboard b;
100: PetscInt nvecs,nscalars,min_size_V,plusk,bs,initv,i,cX_in_proj,cX_in_impr;
101: Mat A,B;
102: KSP ksp;
103: PetscBool t,ipB,ispositive,dynamic;
104: HarmType_t harm;
105: InitType_t init;
106: PetscReal fix;
107: PetscScalar target;
110: /* Setup EPS options and get the problem specification */
111: EPSDavidsonGetBlockSize_Davidson(eps,&bs);
112: if (bs <= 0) bs = 1;
113: if(eps->ncv) {
114: if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The value of ncv must be at least nev");
115: } else if (eps->mpd) eps->ncv = eps->mpd + eps->nev + bs;
116: else if (eps->nev<500)
117: eps->ncv = PetscMin(eps->n-bs,PetscMax(2*eps->nev,eps->nev+15))+bs;
118: else
119: eps->ncv = PetscMin(eps->n-bs,eps->nev+500)+bs;
120: if (!eps->mpd) eps->mpd = eps->ncv;
121: if (eps->mpd > eps->ncv)
122: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The mpd has to be less or equal than ncv");
123: if (eps->mpd < 2)
124: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The mpd has to be greater than 2");
125: if (!eps->max_it) eps->max_it = PetscMax(100*eps->ncv,2*eps->n);
126: if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
127: if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
128: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Wrong value of eps->which");
129: if (!(eps->nev + bs <= eps->ncv))
130: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The ncv has to be greater than nev plus blocksize");
132: EPSDavidsonGetRestart_Davidson(eps,&min_size_V,&plusk);
133: if (!min_size_V) min_size_V = PetscMin(PetscMax(bs,5),eps->mpd/2);
134: if (!(min_size_V+bs <= eps->mpd))
135: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The value of minv must be less than mpd minus blocksize");
136: EPSDavidsonGetInitialSize_Davidson(eps,&initv);
137: if (eps->mpd < initv)
138: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"The initv has to be less or equal than mpd");
140: /* Davidson solvers do not support left eigenvectors */
141: if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
143: /* Set STPrecond as the default ST */
144: if (!((PetscObject)eps->OP)->type_name) {
145: STSetType(eps->OP,STPRECOND);
146: }
147: STPrecondSetKSPHasMat(eps->OP,PETSC_FALSE);
149: /* Change the default sigma to inf if necessary */
150: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
151: eps->which == EPS_LARGEST_IMAGINARY) {
152: STSetDefaultShift(eps->OP,PETSC_MAX_REAL);
153: }
154:
155: /* Davidson solvers only support STPRECOND */
156: STSetUp(eps->OP);
157: PetscObjectTypeCompare((PetscObject)eps->OP,STPRECOND,&t);
158: if (!t) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_SUP,"%s only works with precond spectral transformation",
159: ((PetscObject)eps)->type_name);
161: /* Setup problem specification in dvd */
162: STGetOperators(eps->OP,&A,&B);
163: EPSReset_Davidson(eps);
164: PetscMemzero(dvd,sizeof(dvdDashboard));
165: dvd->A = A; dvd->B = eps->isgeneralized? B : PETSC_NULL;
166: ispositive = eps->ispositive;
167: dvd->sA = DVD_MAT_IMPLICIT |
168: (eps->ishermitian? DVD_MAT_HERMITIAN : 0) |
169: ((ispositive && !eps->isgeneralized) ? DVD_MAT_POS_DEF : 0);
170: /* Asume -eps_hermitian means hermitian-definite in generalized problems */
171: if (!ispositive && !eps->isgeneralized && eps->ishermitian) ispositive = PETSC_TRUE;
172: if (!eps->isgeneralized)
173: dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY |
174: DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
175: else
176: dvd->sB = DVD_MAT_IMPLICIT |
177: (eps->ishermitian? DVD_MAT_HERMITIAN : 0) |
178: (ispositive? DVD_MAT_POS_DEF : 0);
179: ipB = (dvd->B && data->ipB != EPS_ORTH_I && DVD_IS(dvd->sB,DVD_MAT_HERMITIAN))?PETSC_TRUE:PETSC_FALSE;
180: if (data->ipB != EPS_ORTH_I && !ipB) data->ipB = EPS_ORTH_I;
181: dvd->correctXnorm = ipB;
182: dvd->sEP = ((!eps->isgeneralized || (eps->isgeneralized && ipB))? DVD_EP_STD : 0) |
183: (ispositive? DVD_EP_HERMITIAN : 0) |
184: ((eps->problem_type == EPS_GHIEP && ipB) ? DVD_EP_INDEFINITE : 0);
185: dvd->nev = eps->nev;
186: dvd->which = eps->which;
187: dvd->withTarget = PETSC_TRUE;
188: switch(eps->which) {
189: case EPS_TARGET_MAGNITUDE:
190: case EPS_TARGET_IMAGINARY:
191: dvd->target[0] = target = eps->target; dvd->target[1] = 1.0;
192: break;
194: case EPS_TARGET_REAL:
195: dvd->target[0] = PetscRealPart(target = eps->target); dvd->target[1] = 1.0;
196: break;
198: case EPS_LARGEST_REAL:
199: case EPS_LARGEST_MAGNITUDE:
200: case EPS_LARGEST_IMAGINARY: /* TODO: think about this case */
201: default:
202: dvd->target[0] = 1.0; dvd->target[1] = target = 0.0;
203: break;
204:
205: case EPS_SMALLEST_MAGNITUDE:
206: case EPS_SMALLEST_REAL:
207: case EPS_SMALLEST_IMAGINARY: /* TODO: think about this case */
208: dvd->target[0] = target = 0.0; dvd->target[1] = 1.0;
209: break;
211: case EPS_WHICH_USER:
212: STGetShift(eps->OP,&target);
213: dvd->target[0] = target; dvd->target[1] = 1.0;
214: break;
216: case EPS_ALL:
217: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported option: which == EPS_ALL");
218: break;
219: }
220: dvd->tol = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
221: dvd->eps = eps;
223: /* Setup the extraction technique */
224: if (!eps->extraction) {
225: if (ipB || ispositive) eps->extraction = EPS_RITZ;
226: else {
227: switch(eps->which) {
228: case EPS_TARGET_REAL: case EPS_TARGET_MAGNITUDE: case EPS_TARGET_IMAGINARY:
229: case EPS_SMALLEST_MAGNITUDE: case EPS_SMALLEST_REAL: case EPS_SMALLEST_IMAGINARY:
230: eps->extraction = EPS_HARMONIC;
231: break;
232: case EPS_LARGEST_REAL: case EPS_LARGEST_MAGNITUDE: case EPS_LARGEST_IMAGINARY:
233: eps->extraction = EPS_HARMONIC_LARGEST;
234: break;
235: default:
236: eps->extraction = EPS_RITZ;
237: }
238: }
239: }
240: switch(eps->extraction) {
241: case EPS_RITZ: harm = DVD_HARM_NONE; break;
242: case EPS_HARMONIC: harm = DVD_HARM_RR; break;
243: case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
244: case EPS_HARMONIC_RIGHT: harm = DVD_HARM_REIGS; break;
245: case EPS_HARMONIC_LARGEST: harm = DVD_HARM_LEIGS; break;
246: default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
247: }
249: /* Setup the type of starting subspace */
250: EPSDavidsonGetKrylovStart_Davidson(eps,&t);
251: init = (!t)? DVD_INITV_CLASSIC : DVD_INITV_KRYLOV;
253: /* Setup the presence of converged vectors in the projected problem and in the projector */
254: EPSDavidsonGetWindowSizes_Davidson(eps,&cX_in_impr,&cX_in_proj);
255: if (min_size_V <= cX_in_proj) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"minv has to be greater than qwindow");
256: if (bs > 1 && cX_in_impr > 0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported option: pwindow > 0 and bs > 1");
258: /* Setup IP */
259: if (ipB && dvd->B) {
260: IPSetMatrix(eps->ip,dvd->B);
261: } else {
262: IPSetMatrix(eps->ip,PETSC_NULL);
263: }
265: /* Get the fix parameter */
266: EPSDavidsonGetFix_Davidson(eps,&fix);
268: /* Get whether the stopping criterion is used */
269: EPSDavidsonGetConstantCorrectionTolerance_Davidson(eps,&dynamic);
271: /* Orthonormalize the deflation space */
272: dvd_orthV(eps->ip,PETSC_NULL,0,PETSC_NULL,0,eps->defl,0,
273: PetscAbs(eps->nds),PETSC_NULL,eps->rand);
275: /* Preconfigure dvd */
276: STGetKSP(eps->OP,&ksp);
277: dvd_schm_basic_preconf(dvd,&b,eps->mpd,min_size_V,bs,
278: initv,
279: PetscAbs(eps->nini),
280: plusk,harm,
281: ksp,init,eps->trackall,
282: data->ipB,cX_in_proj,cX_in_impr,
283: data->scheme);
284:
286: /* Allocate memory */
287: nvecs = b.max_size_auxV + b.own_vecs;
288: nscalars = b.own_scalars + b.max_size_auxS;
289: PetscMalloc(nscalars*sizeof(PetscScalar),&data->wS);
290: VecDuplicateVecs(eps->t,nvecs,&data->wV);
291: data->size_wV = nvecs;
292: b.free_vecs = data->wV;
293: b.free_scalars = data->wS;
294: dvd->auxV = data->wV + b.own_vecs;
295: dvd->auxS = b.free_scalars + b.own_scalars;
296: dvd->size_auxV = b.max_size_auxV;
297: dvd->size_auxS = b.max_size_auxS;
299: eps->errest_left = PETSC_NULL;
300: PetscMalloc(eps->ncv*sizeof(PetscInt),&eps->perm);
301: for(i=0; i<eps->ncv; i++) eps->perm[i] = i;
303: /* Configure dvd for a basic GD */
304: dvd_schm_basic_conf(dvd,&b,eps->mpd,min_size_V,bs,
305: initv,
306: PetscAbs(eps->nini),plusk,
307: eps->ip,harm,dvd->withTarget,
308: target,ksp,
309: fix,init,eps->trackall,
310: data->ipB,cX_in_proj,cX_in_impr,dynamic,
311: data->scheme);
312:
314: /* Associate the eigenvalues to the EPS */
315: eps->eigr = dvd->real_eigr;
316: eps->eigi = dvd->real_eigi;
317: eps->errest = dvd->real_errest;
318: eps->V = dvd->real_V;
320:
321: return(0);
322: }
326: PetscErrorCode EPSSolve_Davidson(EPS eps)
327: {
328: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
329: dvdDashboard *d = &data->ddb;
333: /* Call the starting routines */
334: DVD_FL_CALL(d->startList,d);
336: for(eps->its=0; eps->its < eps->max_it; eps->its++) {
337: /* Initialize V, if it is needed */
338: if (d->size_V == 0) { d->initV(d); }
340: /* Find the best approximated eigenpairs in V, X */
341: d->calcPairs(d);
343: /* Test for convergence */
344: if (eps->nconv >= eps->nev) break;
346: /* Expand the subspace */
347: d->updateV(d);
349: /* Monitor progress */
350: eps->nconv = d->nconv;
351: EPSMonitor(eps,eps->its+1,eps->nconv,eps->eigr,eps->eigi,eps->errest,d->size_V+d->size_cX);
352: }
354: /* Call the ending routines */
355: DVD_FL_CALL(d->endList,d);
357: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
358: else eps->reason = EPS_DIVERGED_ITS;
360: return(0);
361: }
365: PetscErrorCode EPSReset_Davidson(EPS eps)
366: {
367: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
368: dvdDashboard *dvd = &data->ddb;
372: /* Call step destructors and destroys the list */
373: DVD_FL_CALL(dvd->destroyList,dvd);
374: DVD_FL_DEL(dvd->destroyList);
375: DVD_FL_DEL(dvd->startList);
376: DVD_FL_DEL(dvd->endList);
378: if (data->size_wV > 0) {
379: VecDestroyVecs(data->size_wV,&data->wV);
380: }
381: PetscFree(data->wS);
382: PetscFree(eps->perm);
383: return(0);
384: }
388: PetscErrorCode EPSView_Davidson(EPS eps,PetscViewer viewer)
389: {
391: PetscBool isascii,opb;
392: PetscInt opi,opi0;
393: const char* name;
394: Method_t meth;
395: EPSOrthType borth;
398: name = ((PetscObject)eps)->type_name;
399: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
400: if (!isascii) SETERRQ2(((PetscObject)eps)->comm,1,"Viewer type %s not supported for %s",((PetscObject)viewer)->type_name,name);
401:
402: EPSDavidsonGetMethod_Davidson(eps,&meth);
403: if (meth==DVD_METH_GD2) {
404: PetscViewerASCIIPrintf(viewer," Davidson: using double expansion variant (GD2)\n");
405: }
406: EPSDavidsonGetBOrth_Davidson(eps,&borth);
407: switch (borth) {
408: case EPS_ORTH_I:
409: PetscViewerASCIIPrintf(viewer," Davidson: search subspace is orthogonalized\n");
410: break;
411: case EPS_ORTH_B:
412: PetscViewerASCIIPrintf(viewer," Davidson: search subspace is B-orthogonalized\n");
413: break;
414: case EPS_ORTH_BOPT:
415: PetscViewerASCIIPrintf(viewer," Davidson: search subspace is B-orthogonalized with an optimized method\n");
416: break;
417: }
418: EPSDavidsonGetBlockSize_Davidson(eps,&opi);
419: PetscViewerASCIIPrintf(viewer," Davidson: block size=%D\n",opi);
420: EPSDavidsonGetKrylovStart_Davidson(eps,&opb);
421: if (!opb) {
422: PetscViewerASCIIPrintf(viewer," Davidson: type of the initial subspace: non-Krylov\n");
423: } else {
424: PetscViewerASCIIPrintf(viewer," Davidson: type of the initial subspace: Krylov\n");
425: }
426: EPSDavidsonGetRestart_Davidson(eps,&opi,&opi0);
427: PetscViewerASCIIPrintf(viewer," Davidson: size of the subspace after restarting: %D\n",opi);
428: PetscViewerASCIIPrintf(viewer," Davidson: number of vectors after restarting from the previous iteration: %D\n",opi0);
429: return(0);
430: }
434: PetscErrorCode EPSDavidsonSetKrylovStart_Davidson(EPS eps,PetscBool krylovstart)
435: {
436: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
439: data->krylovstart = krylovstart;
440: return(0);
441: }
445: PetscErrorCode EPSDavidsonGetKrylovStart_Davidson(EPS eps,PetscBool *krylovstart)
446: {
447: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
450: *krylovstart = data->krylovstart;
451: return(0);
452: }
456: PetscErrorCode EPSDavidsonSetBlockSize_Davidson(EPS eps,PetscInt blocksize)
457: {
458: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
461: if(blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
462: if(blocksize <= 0)
463: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value");
464: data->blocksize = blocksize;
465: return(0);
466: }
470: PetscErrorCode EPSDavidsonGetBlockSize_Davidson(EPS eps,PetscInt *blocksize)
471: {
472: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
475: *blocksize = data->blocksize;
476: return(0);
477: }
481: PetscErrorCode EPSDavidsonSetRestart_Davidson(EPS eps,PetscInt minv,PetscInt plusk)
482: {
483: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
486: if(minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 5;
487: if(minv <= 0)
488: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value");
489: if(plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = 5;
490: if(plusk < 0)
491: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value");
492: data->minv = minv;
493: data->plusk = plusk;
494: return(0);
495: }
499: PetscErrorCode EPSDavidsonGetRestart_Davidson(EPS eps,PetscInt *minv,PetscInt *plusk)
500: {
501: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
504: if (minv) *minv = data->minv;
505: if (plusk) *plusk = data->plusk;
506: return(0);
507: }
511: PetscErrorCode EPSDavidsonGetInitialSize_Davidson(EPS eps,PetscInt *initialsize)
512: {
513: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
516: *initialsize = data->initialsize;
517: return(0);
518: }
522: PetscErrorCode EPSDavidsonSetInitialSize_Davidson(EPS eps,PetscInt initialsize)
523: {
524: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
527: if(initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 5;
528: if(initialsize <= 0)
529: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value");
530: data->initialsize = initialsize;
531: return(0);
532: }
536: PetscErrorCode EPSDavidsonGetFix_Davidson(EPS eps,PetscReal *fix)
537: {
538: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
541: *fix = data->fix;
542: return(0);
543: }
547: PetscErrorCode EPSDavidsonSetFix_Davidson(EPS eps,PetscReal fix)
548: {
549: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
552: if(fix == PETSC_DEFAULT || fix == PETSC_DECIDE) fix = 0.01;
553: if(fix < 0.0)
554: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
555: data->fix = fix;
556: return(0);
557: }
561: PetscErrorCode EPSDavidsonSetBOrth_Davidson(EPS eps,EPSOrthType borth)
562: {
563: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
566: data->ipB = borth;
567: return(0);
568: }
572: PetscErrorCode EPSDavidsonGetBOrth_Davidson(EPS eps,EPSOrthType *borth)
573: {
574: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
577: *borth = data->ipB;
578: return(0);
579: }
583: PetscErrorCode EPSDavidsonSetConstantCorrectionTolerance_Davidson(EPS eps,PetscBool constant)
584: {
585: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
588: data->dynamic = !constant?PETSC_TRUE:PETSC_FALSE;
589: return(0);
590: }
594: PetscErrorCode EPSDavidsonGetConstantCorrectionTolerance_Davidson(EPS eps,PetscBool *constant)
595: {
596: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
599: *constant = !data->dynamic?PETSC_TRUE:PETSC_FALSE;
600: return(0);
601: }
606: PetscErrorCode EPSDavidsonSetWindowSizes_Davidson(EPS eps,PetscInt pwindow,PetscInt qwindow)
607: {
608: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
611: if(pwindow == PETSC_DEFAULT || pwindow == PETSC_DECIDE) pwindow = 0;
612: if(pwindow < 0)
613: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid pwindow value");
614: if(qwindow == PETSC_DEFAULT || qwindow == PETSC_DECIDE) qwindow = 0;
615: if(qwindow < 0)
616: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid qwindow value");
617: data->cX_in_proj = qwindow;
618: data->cX_in_impr = pwindow;
619: return(0);
620: }
624: PetscErrorCode EPSDavidsonGetWindowSizes_Davidson(EPS eps,PetscInt *pwindow,PetscInt *qwindow)
625: {
626: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
629: if (pwindow) *pwindow = data->cX_in_impr;
630: if (qwindow) *qwindow = data->cX_in_proj;
631: return(0);
632: }
636: PetscErrorCode EPSDavidsonSetMethod_Davidson(EPS eps,Method_t method)
637: {
638: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
641: data->scheme = method;
642: return(0);
643: }
647: PetscErrorCode EPSDavidsonGetMethod_Davidson(EPS eps,Method_t *method)
648: {
649: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
652: *method = data->scheme;
653: return(0);
654: }
658: /*
659: EPSComputeVectors_Davidson - Compute eigenvectors from the vectors
660: provided by the eigensolver. This version is intended for solvers
661: that provide Schur vectors from the QZ decompositon. Given the partial
662: Schur decomposition OP*V=V*T, the following steps are performed:
663: 1) compute eigenvectors of (S,T): S*Z=T*Z*D
664: 2) compute eigenvectors of OP: X=V*Z
665: If left eigenvectors are required then also do Z'*T=D*Z', Y=W*Z
666: */
667: PetscErrorCode EPSComputeVectors_Davidson(EPS eps)
668: {
670: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
671: dvdDashboard *d = &data->ddb;
672: PetscScalar *pX,*cS,*cT;
673: PetscInt ld;
677: if (d->cS) {
678: /* Compute the eigenvectors associated to (cS, cT) */
679: DSSetDimensions(d->conv_ps,d->size_cS,PETSC_IGNORE,0,0);
680: DSGetLeadingDimension(d->conv_ps,&ld);
681: DSGetArray(d->conv_ps,DS_MAT_A,&cS);
682: SlepcDenseCopyTriang(cS,0,ld,d->cS,0,d->ldcS,d->size_cS,d->size_cS);
683: DSRestoreArray(d->conv_ps,DS_MAT_A,&cS);
684: if (d->cT) {
685: DSGetArray(d->conv_ps,DS_MAT_B,&cT);
686: SlepcDenseCopyTriang(cT,0,ld,d->cT,0,d->ldcT,d->size_cS,d->size_cS);
687: DSRestoreArray(d->conv_ps,DS_MAT_B,&cT);
688: }
689: DSSetState(d->conv_ps,DS_STATE_RAW);
690: DSSolve(d->conv_ps,eps->eigr,eps->eigi);
691: DSVectors(d->conv_ps,DS_MAT_X,PETSC_NULL,PETSC_NULL);
692: DSNormalize(d->conv_ps,DS_MAT_X,-1);
694: /* V <- cX * pX */
695: DSGetArray(d->conv_ps,DS_MAT_X,&pX);
696: SlepcUpdateVectorsZ(eps->V,0.0,1.0,d->cX,d->size_cX,pX,ld,d->nconv,d->nconv);
697: DSRestoreArray(d->conv_ps,DS_MAT_X,&pX);
698: }
700: eps->evecsavailable = PETSC_TRUE;
701: return(0);
702: }