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);
245: if (min_size_V <= cX_in_proj) {
246: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"minv has to be greater than qwindow");
247: }
248: if (bs > 1 && cX_in_impr > 0) {
249: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported option: pwindow > 0 and bs > 1");
250: }
252: /* Setup IP */
253: if (ipB && dvd->B) {
254: IPSetMatrix(eps->ip,dvd->B);
255: } else {
256: IPSetMatrix(eps->ip,PETSC_NULL);
257: }
259: /* Get the fix parameter */
260: EPSDavidsonGetFix_Davidson(eps,&fix);
262: /* Get whether the stopping criterion is used */
263: EPSDavidsonGetConstantCorrectionTolerance_Davidson(eps,&dynamic);
265: /* Orthonormalize the DS */
266: dvd_orthV(eps->ip,PETSC_NULL,0,PETSC_NULL,0,eps->DS,0,
267: PetscAbs(eps->nds),PETSC_NULL,eps->rand);
269: /* Preconfigure dvd */
270: STGetKSP(eps->OP,&ksp);
271: dvd_schm_basic_preconf(dvd,&b,eps->mpd,min_size_V,bs,
272: initv,
273: PetscAbs(eps->nini),
274: plusk,harm,
275: ksp,init,eps->trackall,
276: ipB?DVD_ORTHOV_BOneMV:DVD_ORTHOV_I,cX_in_proj,cX_in_impr);
277:
279: /* Allocate memory */
280: nvecs = b.max_size_auxV + b.own_vecs;
281: nscalars = b.own_scalars + b.max_size_auxS;
282: PetscMalloc(nscalars*sizeof(PetscScalar),&data->wS);
283: VecDuplicateVecs(eps->t,nvecs,&data->wV);
284: data->size_wV = nvecs;
285: b.free_vecs = data->wV;
286: b.free_scalars = data->wS;
287: dvd->auxV = data->wV + b.own_vecs;
288: dvd->auxS = b.free_scalars + b.own_scalars;
289: dvd->size_auxV = b.max_size_auxV;
290: dvd->size_auxS = b.max_size_auxS;
292: eps->errest_left = PETSC_NULL;
293: PetscMalloc(eps->ncv*sizeof(PetscInt),&eps->perm);
294: for(i=0; i<eps->ncv; i++) eps->perm[i] = i;
296: /* Configure dvd for a basic GD */
297: dvd_schm_basic_conf(dvd,&b,eps->mpd,min_size_V,bs,
298: initv,
299: PetscAbs(eps->nini),plusk,
300: eps->ip,harm,dvd->withTarget,
301: target,ksp,
302: fix,init,eps->trackall,
303: ipB?DVD_ORTHOV_BOneMV:DVD_ORTHOV_I,cX_in_proj,cX_in_impr,dynamic);
304:
306: /* Associate the eigenvalues to the EPS */
307: eps->eigr = dvd->real_eigr;
308: eps->eigi = dvd->real_eigi;
309: eps->errest = dvd->real_errest;
310: eps->V = dvd->real_V;
312:
313: return(0);
314: }
318: PetscErrorCode EPSSolve_Davidson(EPS eps)
319: {
320: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
321: dvdDashboard *d = &data->ddb;
325: /* Call the starting routines */
326: DVD_FL_CALL(d->startList,d);
328: for(eps->its=0; eps->its < eps->max_it; eps->its++) {
329: /* Initialize V, if it is needed */
330: if (d->size_V == 0) { d->initV(d); }
332: /* Find the best approximated eigenpairs in V, X */
333: d->calcPairs(d);
335: /* Test for convergence */
336: if (eps->nconv >= eps->nev) break;
338: /* Expand the subspace */
339: d->updateV(d);
341: /* Monitor progress */
342: eps->nconv = d->nconv;
343: EPSMonitor(eps,eps->its+1,eps->nconv,eps->eigr,eps->eigi,eps->errest,d->size_V+d->size_cX);
344: }
346: /* Call the ending routines */
347: DVD_FL_CALL(d->endList,d);
349: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
350: else eps->reason = EPS_DIVERGED_ITS;
352: return(0);
353: }
357: PetscErrorCode EPSReset_Davidson(EPS eps)
358: {
359: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
360: dvdDashboard *dvd = &data->ddb;
364: /* Call step destructors and destroys the list */
365: DVD_FL_CALL(dvd->destroyList,dvd);
366: DVD_FL_DEL(dvd->destroyList);
367: DVD_FL_DEL(dvd->startList);
368: DVD_FL_DEL(dvd->endList);
370: if (data->size_wV > 0) {
371: VecDestroyVecs(data->size_wV,&data->wV);
372: }
373: PetscFree(data->wS);
374: PetscFree(eps->perm);
375: return(0);
376: }
380: PetscErrorCode EPSView_Davidson(EPS eps,PetscViewer viewer)
381: {
383: PetscBool isascii,opb;
384: PetscInt opi,opi0;
385: const char* name;
388: name = ((PetscObject)eps)->type_name;
389: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
390: if (!isascii) {
391: SETERRQ2(((PetscObject)eps)->comm,1,"Viewer type %s not supported for %s",((PetscObject)viewer)->type_name,name);
392: }
393:
394: EPSDavidsonGetBOrth_Davidson(eps,&opb);
395: EPSDavidsonGetBlockSize_Davidson(eps,&opi);
396: if(!opb) {
397: PetscViewerASCIIPrintf(viewer," Davidson: search subspace is I-orthogonalized\n");
398: } else {
399: PetscViewerASCIIPrintf(viewer," Davidson: search subspace is B-orthogonalized\n");
400: }
401: PetscViewerASCIIPrintf(viewer," Davidson: block size=%D\n",opi);
402: EPSDavidsonGetKrylovStart_Davidson(eps,&opb);
403: if(!opb) {
404: PetscViewerASCIIPrintf(viewer," Davidson: type of the initial subspace: non-Krylov\n");
405: } else {
406: PetscViewerASCIIPrintf(viewer," Davidson: type of the initial subspace: Krylov\n");
407: }
408: EPSDavidsonGetRestart_Davidson(eps,&opi,&opi0);
409: PetscViewerASCIIPrintf(viewer," Davidson: size of the subspace after restarting: %D\n",opi);
410: PetscViewerASCIIPrintf(viewer," Davidson: number of vectors after restarting from the previous iteration: %D\n",opi0);
411: return(0);
412: }
416: PetscErrorCode EPSDavidsonSetKrylovStart_Davidson(EPS eps,PetscBool krylovstart)
417: {
418: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
421: data->krylovstart = krylovstart;
422: return(0);
423: }
427: PetscErrorCode EPSDavidsonGetKrylovStart_Davidson(EPS eps,PetscBool *krylovstart)
428: {
429: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
432: *krylovstart = data->krylovstart;
433: return(0);
434: }
438: PetscErrorCode EPSDavidsonSetBlockSize_Davidson(EPS eps,PetscInt blocksize)
439: {
440: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
443: if(blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
444: if(blocksize <= 0)
445: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value");
446: data->blocksize = blocksize;
447: return(0);
448: }
452: PetscErrorCode EPSDavidsonGetBlockSize_Davidson(EPS eps,PetscInt *blocksize)
453: {
454: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
457: *blocksize = data->blocksize;
458: return(0);
459: }
463: PetscErrorCode EPSDavidsonSetRestart_Davidson(EPS eps,PetscInt minv,PetscInt plusk)
464: {
465: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
468: if(minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 5;
469: if(minv <= 0)
470: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value");
471: if(plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = 5;
472: if(plusk < 0)
473: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value");
474: data->minv = minv;
475: data->plusk = plusk;
476: return(0);
477: }
481: PetscErrorCode EPSDavidsonGetRestart_Davidson(EPS eps,PetscInt *minv,PetscInt *plusk)
482: {
483: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
486: *minv = data->minv;
487: *plusk = data->plusk;
488: return(0);
489: }
493: PetscErrorCode EPSDavidsonGetInitialSize_Davidson(EPS eps,PetscInt *initialsize)
494: {
495: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
498: *initialsize = data->initialsize;
499: return(0);
500: }
504: PetscErrorCode EPSDavidsonSetInitialSize_Davidson(EPS eps,PetscInt initialsize)
505: {
506: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
509: if(initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 5;
510: if(initialsize <= 0)
511: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value");
512: data->initialsize = initialsize;
513: return(0);
514: }
518: PetscErrorCode EPSDavidsonGetFix_Davidson(EPS eps,PetscReal *fix)
519: {
520: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
523: *fix = data->fix;
524: return(0);
525: }
529: PetscErrorCode EPSDavidsonSetFix_Davidson(EPS eps,PetscReal fix)
530: {
531: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
534: if(fix == PETSC_DEFAULT || fix == PETSC_DECIDE) fix = 0.01;
535: if(fix < 0.0)
536: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
537: data->fix = fix;
538: return(0);
539: }
543: PetscErrorCode EPSDavidsonSetBOrth_Davidson(EPS eps,PetscBool borth)
544: {
545: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
548: data->ipB = borth;
549: return(0);
550: }
554: PetscErrorCode EPSDavidsonGetBOrth_Davidson(EPS eps,PetscBool *borth)
555: {
556: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
559: *borth = data->ipB;
560: return(0);
561: }
565: PetscErrorCode EPSDavidsonSetConstantCorrectionTolerance_Davidson(EPS eps,PetscBool constant)
566: {
567: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
570: data->dynamic = !constant?PETSC_TRUE:PETSC_FALSE;
571: return(0);
572: }
576: PetscErrorCode EPSDavidsonGetConstantCorrectionTolerance_Davidson(EPS eps,PetscBool *constant)
577: {
578: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
581: *constant = !data->dynamic?PETSC_TRUE:PETSC_FALSE;
582: return(0);
583: }
588: PetscErrorCode EPSDavidsonSetWindowSizes_Davidson(EPS eps,PetscInt pwindow,PetscInt qwindow)
589: {
590: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
593: if(pwindow == PETSC_DEFAULT || pwindow == PETSC_DECIDE) pwindow = 0;
594: if(pwindow < 0)
595: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid pwindow value");
596: if(qwindow == PETSC_DEFAULT || qwindow == PETSC_DECIDE) qwindow = 0;
597: if(qwindow < 0)
598: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid qwindow value");
599: data->cX_in_proj = qwindow;
600: data->cX_in_impr = pwindow;
601: return(0);
602: }
606: PetscErrorCode EPSDavidsonGetWindowSizes_Davidson(EPS eps,PetscInt *pwindow,PetscInt *qwindow)
607: {
608: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
611: *pwindow = data->cX_in_impr;
612: *qwindow = data->cX_in_proj;
613: return(0);
614: }
619: /*
620: EPSComputeVectors_Davidson - Compute eigenvectors from the vectors
621: provided by the eigensolver. This version is intended for solvers
622: that provide Schur vectors from the QZ decompositon. Given the partial
623: Schur decomposition OP*V=V*T, the following steps are performed:
624: 1) compute eigenvectors of (S,T): S*Z=T*Z*D
625: 2) compute eigenvectors of OP: X=V*Z
626: If left eigenvectors are required then also do Z'*Tl=D*Z', Y=W*Z
627: */
628: PetscErrorCode EPSComputeVectors_Davidson(EPS eps)
629: {
631: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
632: dvdDashboard *d = &data->ddb;
633: PetscScalar *pX,*auxS;
634: PetscInt size_auxS,i,j;
638: if (d->cS) {
639: /* Compute the eigenvectors associated to (cS, cT) */
640: PetscMalloc(sizeof(PetscScalar)*d->nconv*d->nconv,&pX);
641: size_auxS = 6*d->nconv;
642: PetscMalloc(sizeof(PetscScalar)*size_auxS,&auxS);
643: #if defined(PETSC_USE_COMPLEX)
644: for (i=0; i<d->nconv; i++) {
645: for (j=i+1; j<d->nconv; j++) {
646: d->cS[d->ldcS*i+j] = 0.0;
647: if (d->cT) d->cT[d->ldcS*i+j] = 0.0;
648: }
649: }
650: #else
651: for (i=0; i<d->nconv; i++) {
652: if (d->cS[d->ldcS*i+i+1] != 0.0 && d->ceigi[i] != 0.0) {
653: for (j=i+2; j<d->nconv; j++) d->cS[d->ldcS*i+j] = 0.0;
654: for (j=i+2; j<d->nconv; j++) d->cS[d->ldcS*(i+1)+j] = 0.0;
655: if (d->cT) {
656: d->cT[d->ldcS*(i+1)+i] = 0.0;
657: for (j=i+1; j<d->nconv; j++) d->cT[d->ldcS*i+j] = 0.0;
658: for (j=i+2; j<d->nconv; j++) d->cT[d->ldcS*(i+1)+j] = 0.0;
659: }
660: i++;
661: } else {
662: for (j=i+1; j<d->nconv; j++) {
663: d->cS[d->ldcS*i+j] = 0.0;
664: if (d->cT) d->cT[d->ldcS*i+j] = 0.0;
665: }
666: }
667: }
668: #endif
669: dvd_compute_eigenvectors(d->nconv,d->cS,d->ldcS,d->cT,d->ldcT,
670: pX,d->nconv,PETSC_NULL,0,auxS,
671: size_auxS,PETSC_FALSE);
672:
673: /* pX[i] <- pX[i] / ||pX[i]|| */
674: SlepcDenseNorm(pX,d->nconv,d->nconv,d->nconv,d->ceigi);
675:
676: /* V <- cX * pX */
677: SlepcUpdateVectorsZ(eps->V,0.0,1.0,d->cX,d->size_cX,pX,
678: d->nconv,d->nconv,d->nconv);
679:
680: PetscFree(pX);
681: PetscFree(auxS);
682: }
684: eps->evecsavailable = PETSC_TRUE;
685: return(0);
686: }