Actual source code: davidson.c
slepc-3.7.4 2017-05-17
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-2016, 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: static PetscBool cited = PETSC_FALSE;
27: static const char citation[] =
28: "@Article{slepc-davidson,\n"
29: " author = \"E. Romero and J. E. Roman\",\n"
30: " title = \"A parallel implementation of {Davidson} methods for large-scale eigenvalue problems in {SLEPc}\",\n"
31: " journal = \"{ACM} Trans. Math. Software\",\n"
32: " volume = \"40\",\n"
33: " number = \"2\",\n"
34: " pages = \"13:1--13:29\",\n"
35: " year = \"2014,\"\n"
36: " doi = \"http://dx.doi.org/10.1145/2543696\"\n"
37: "}\n";
41: PetscErrorCode EPSSetUp_XD(EPS eps)
42: {
44: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
45: dvdDashboard *dvd = &data->ddb;
46: dvdBlackboard b;
47: PetscInt min_size_V,plusk,bs,initv,i,cX_in_proj,cX_in_impr,nmat;
48: Mat A,B;
49: KSP ksp;
50: PetscBool t,ipB,ispositive,dynamic;
51: HarmType_t harm;
52: InitType_t init;
53: PetscReal fix;
54: PetscScalar target;
57: /* Setup EPS options and get the problem specification */
58: EPSXDGetBlockSize_XD(eps,&bs);
59: if (bs <= 0) bs = 1;
60: if (eps->ncv) {
61: if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of ncv must be at least nev");
62: } else if (eps->mpd) eps->ncv = eps->mpd + eps->nev + bs;
63: else if (eps->nev<500) eps->ncv = PetscMin(eps->n-bs,PetscMax(2*eps->nev,eps->nev+15))+bs;
64: else eps->ncv = PetscMin(eps->n-bs,eps->nev+500)+bs;
65: if (!eps->mpd) eps->mpd = eps->ncv;
66: if (eps->mpd > eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd has to be less or equal than ncv");
67: if (eps->mpd < 2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd has to be greater than 2");
68: if (!eps->max_it) eps->max_it = PetscMax(100*eps->ncv,2*eps->n);
69: if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
70: 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");
71: if (!(eps->nev + bs <= eps->ncv)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The ncv has to be greater than nev plus blocksize");
72: if (eps->trueres) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"-eps_true_residual is temporally disable in this solver.");
74: EPSXDGetRestart_XD(eps,&min_size_V,&plusk);
75: if (!min_size_V) min_size_V = PetscMin(PetscMax(bs,5),eps->mpd/2);
76: 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");
77: EPSXDGetInitialSize_XD(eps,&initv);
78: if (eps->mpd < initv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The initv has to be less or equal than mpd");
80: /* Set STPrecond as the default ST */
81: if (!((PetscObject)eps->st)->type_name) {
82: STSetType(eps->st,STPRECOND);
83: }
84: STPrecondSetKSPHasMat(eps->st,PETSC_FALSE);
86: /* Change the default sigma to inf if necessary */
87: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL || eps->which == EPS_LARGEST_IMAGINARY) {
88: STSetDefaultShift(eps->st,PETSC_MAX_REAL);
89: }
91: /* Davidson solvers only support STPRECOND */
92: STSetUp(eps->st);
93: PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&t);
94: if (!t) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"%s only works with precond spectral transformation",((PetscObject)eps)->type_name);
96: /* Setup problem specification in dvd */
97: STGetNumMatrices(eps->st,&nmat);
98: STGetOperators(eps->st,0,&A);
99: if (nmat>1) { STGetOperators(eps->st,1,&B); }
100: EPSReset_XD(eps);
101: PetscMemzero(dvd,sizeof(dvdDashboard));
102: dvd->A = A; dvd->B = eps->isgeneralized? B: NULL;
103: ispositive = eps->ispositive;
104: dvd->sA = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN: 0) | ((ispositive && !eps->isgeneralized) ? DVD_MAT_POS_DEF: 0);
105: /* Asume -eps_hermitian means hermitian-definite in generalized problems */
106: if (!ispositive && !eps->isgeneralized && eps->ishermitian) ispositive = PETSC_TRUE;
107: if (!eps->isgeneralized) dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY | DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
108: else dvd->sB = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN: 0) | (ispositive? DVD_MAT_POS_DEF: 0);
109: ipB = (dvd->B && data->ipB && DVD_IS(dvd->sB,DVD_MAT_HERMITIAN))?PETSC_TRUE:PETSC_FALSE;
110: if (data->ipB && !ipB) data->ipB = PETSC_FALSE;
111: dvd->correctXnorm = ipB;
112: 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);
113: dvd->nev = eps->nev;
114: dvd->which = eps->which;
115: dvd->withTarget = PETSC_TRUE;
116: switch (eps->which) {
117: case EPS_TARGET_MAGNITUDE:
118: case EPS_TARGET_IMAGINARY:
119: dvd->target[0] = target = eps->target;
120: dvd->target[1] = 1.0;
121: break;
122: case EPS_TARGET_REAL:
123: dvd->target[0] = PetscRealPart(target = eps->target);
124: dvd->target[1] = 1.0;
125: break;
126: case EPS_LARGEST_REAL:
127: case EPS_LARGEST_MAGNITUDE:
128: case EPS_LARGEST_IMAGINARY: /* TODO: think about this case */
129: dvd->target[0] = 1.0;
130: dvd->target[1] = target = 0.0;
131: break;
132: case EPS_SMALLEST_MAGNITUDE:
133: case EPS_SMALLEST_REAL:
134: case EPS_SMALLEST_IMAGINARY: /* TODO: think about this case */
135: dvd->target[0] = target = 0.0;
136: dvd->target[1] = 1.0;
137: break;
138: case EPS_WHICH_USER:
139: STGetShift(eps->st,&target);
140: dvd->target[0] = target;
141: dvd->target[1] = 1.0;
142: break;
143: case EPS_ALL:
144: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported option: which == EPS_ALL");
145: break;
146: default:
147: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported value of option 'which'");
148: }
149: dvd->tol = (eps->tol==PETSC_DEFAULT)? SLEPC_DEFAULT_TOL: eps->tol;
150: dvd->eps = eps;
152: /* Setup the extraction technique */
153: if (!eps->extraction) {
154: if (ipB || ispositive) eps->extraction = EPS_RITZ;
155: else {
156: switch (eps->which) {
157: case EPS_TARGET_REAL:
158: case EPS_TARGET_MAGNITUDE:
159: case EPS_TARGET_IMAGINARY:
160: case EPS_SMALLEST_MAGNITUDE:
161: case EPS_SMALLEST_REAL:
162: case EPS_SMALLEST_IMAGINARY:
163: eps->extraction = EPS_HARMONIC;
164: break;
165: case EPS_LARGEST_REAL:
166: case EPS_LARGEST_MAGNITUDE:
167: case EPS_LARGEST_IMAGINARY:
168: eps->extraction = EPS_HARMONIC_LARGEST;
169: break;
170: default:
171: eps->extraction = EPS_RITZ;
172: }
173: }
174: }
175: switch (eps->extraction) {
176: case EPS_RITZ: harm = DVD_HARM_NONE; break;
177: case EPS_HARMONIC: harm = DVD_HARM_RR; break;
178: case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
179: case EPS_HARMONIC_RIGHT: harm = DVD_HARM_REIGS; break;
180: case EPS_HARMONIC_LARGEST: harm = DVD_HARM_LEIGS; break;
181: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
182: }
184: /* Setup the type of starting subspace */
185: EPSXDGetKrylovStart_XD(eps,&t);
186: init = (!t)? DVD_INITV_CLASSIC : DVD_INITV_KRYLOV;
188: /* Setup the presence of converged vectors in the projected problem and the projector */
189: EPSXDGetWindowSizes_XD(eps,&cX_in_impr,&cX_in_proj);
190: if (cX_in_impr>0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The option pwindow is temporally disable in this solver.");
191: if (cX_in_proj>0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The option qwindow is temporally disable in this solver.");
192: if (min_size_V <= cX_in_proj) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"minv has to be greater than qwindow");
193: if (bs > 1 && cX_in_impr > 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported option: pwindow > 0 and bs > 1");
195: /* Get the fix parameter */
196: EPSXDGetFix_XD(eps,&fix);
198: /* Get whether the stopping criterion is used */
199: EPSJDGetConstCorrectionTol_JD(eps,&dynamic);
201: /* Preconfigure dvd */
202: STGetKSP(eps->st,&ksp);
203: 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->doubleexp);
205: /* Allocate memory */
206: EPSAllocateSolution(eps,0);
208: /* Setup orthogonalization */
209: EPS_SetInnerProduct(eps);
210: if (!(ipB && dvd->B)) {
211: BVSetMatrix(eps->V,NULL,PETSC_FALSE);
212: }
214: for (i=0;i<eps->ncv;i++) eps->perm[i] = i;
216: /* Configure dvd for a basic GD */
217: 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->doubleexp);
218: return(0);
219: }
223: PetscErrorCode EPSSolve_XD(EPS eps)
224: {
225: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
226: dvdDashboard *d = &data->ddb;
227: PetscInt l,k;
231: PetscCitationsRegister(citation,&cited);
232: /* Call the starting routines */
233: EPSDavidsonFLCall(d->startList,d);
235: while (eps->reason == EPS_CONVERGED_ITERATING) {
237: /* Initialize V, if it is needed */
238: BVGetActiveColumns(d->eps->V,&l,&k);
239: if (l == k) { d->initV(d); }
241: /* Find the best approximated eigenpairs in V, X */
242: d->calcPairs(d);
244: /* Test for convergence */
245: (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
246: if (eps->reason != EPS_CONVERGED_ITERATING) break;
248: /* Expand the subspace */
249: d->updateV(d);
251: /* Monitor progress */
252: eps->nconv = d->nconv;
253: eps->its++;
254: BVGetActiveColumns(d->eps->V,&l,&k);
255: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,k);
256: }
258: /* Call the ending routines */
259: EPSDavidsonFLCall(d->endList,d);
260: return(0);
261: }
265: PetscErrorCode EPSReset_XD(EPS eps)
266: {
267: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
268: dvdDashboard *dvd = &data->ddb;
272: /* Call step destructors and destroys the list */
273: EPSDavidsonFLCall(dvd->destroyList,dvd);
274: EPSDavidsonFLDestroy(&dvd->destroyList);
275: EPSDavidsonFLDestroy(&dvd->startList);
276: EPSDavidsonFLDestroy(&dvd->endList);
277: return(0);
278: }
282: PetscErrorCode EPSXDSetKrylovStart_XD(EPS eps,PetscBool krylovstart)
283: {
284: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
287: data->krylovstart = krylovstart;
288: return(0);
289: }
293: PetscErrorCode EPSXDGetKrylovStart_XD(EPS eps,PetscBool *krylovstart)
294: {
295: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
298: *krylovstart = data->krylovstart;
299: return(0);
300: }
304: PetscErrorCode EPSXDSetBlockSize_XD(EPS eps,PetscInt blocksize)
305: {
306: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
309: if (blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
310: if (blocksize <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value");
311: data->blocksize = blocksize;
312: return(0);
313: }
317: PetscErrorCode EPSXDGetBlockSize_XD(EPS eps,PetscInt *blocksize)
318: {
319: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
322: *blocksize = data->blocksize;
323: return(0);
324: }
328: PetscErrorCode EPSXDSetRestart_XD(EPS eps,PetscInt minv,PetscInt plusk)
329: {
330: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
333: if (minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 5;
334: if (minv <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value");
335: if (plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = 5;
336: if (plusk < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value");
337: data->minv = minv;
338: data->plusk = plusk;
339: return(0);
340: }
344: PetscErrorCode EPSXDGetRestart_XD(EPS eps,PetscInt *minv,PetscInt *plusk)
345: {
346: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
349: if (minv) *minv = data->minv;
350: if (plusk) *plusk = data->plusk;
351: return(0);
352: }
356: PetscErrorCode EPSXDGetInitialSize_XD(EPS eps,PetscInt *initialsize)
357: {
358: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
361: *initialsize = data->initialsize;
362: return(0);
363: }
367: PetscErrorCode EPSXDSetInitialSize_XD(EPS eps,PetscInt initialsize)
368: {
369: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
372: if (initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 5;
373: if (initialsize <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value");
374: data->initialsize = initialsize;
375: return(0);
376: }
380: PetscErrorCode EPSXDGetFix_XD(EPS eps,PetscReal *fix)
381: {
382: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
385: *fix = data->fix;
386: return(0);
387: }
391: PetscErrorCode EPSJDSetFix_JD(EPS eps,PetscReal fix)
392: {
393: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
396: if (fix == PETSC_DEFAULT || fix == PETSC_DECIDE) fix = 0.01;
397: if (fix < 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
398: data->fix = fix;
399: return(0);
400: }
404: PetscErrorCode EPSXDSetBOrth_XD(EPS eps,PetscBool borth)
405: {
406: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
409: data->ipB = borth;
410: return(0);
411: }
415: PetscErrorCode EPSXDGetBOrth_XD(EPS eps,PetscBool *borth)
416: {
417: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
420: *borth = data->ipB;
421: return(0);
422: }
426: PetscErrorCode EPSJDSetConstCorrectionTol_JD(EPS eps,PetscBool constant)
427: {
428: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
431: data->dynamic = PetscNot(constant);
432: return(0);
433: }
437: PetscErrorCode EPSJDGetConstCorrectionTol_JD(EPS eps,PetscBool *constant)
438: {
439: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
442: *constant = PetscNot(data->dynamic);
443: return(0);
444: }
448: PetscErrorCode EPSXDSetWindowSizes_XD(EPS eps,PetscInt pwindow,PetscInt qwindow)
449: {
450: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
453: if (pwindow == PETSC_DEFAULT || pwindow == PETSC_DECIDE) pwindow = 0;
454: if (pwindow < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid pwindow value");
455: if (qwindow == PETSC_DEFAULT || qwindow == PETSC_DECIDE) qwindow = 0;
456: if (qwindow < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid qwindow value");
457: data->cX_in_proj = qwindow;
458: data->cX_in_impr = pwindow;
459: return(0);
460: }
464: PetscErrorCode EPSXDGetWindowSizes_XD(EPS eps,PetscInt *pwindow,PetscInt *qwindow)
465: {
466: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
469: if (pwindow) *pwindow = data->cX_in_impr;
470: if (qwindow) *qwindow = data->cX_in_proj;
471: return(0);
472: }
476: /*
477: EPSComputeVectors_XD - Compute eigenvectors from the vectors
478: provided by the eigensolver. This version is intended for solvers
479: that provide Schur vectors from the QZ decomposition. Given the partial
480: Schur decomposition OP*V=V*T, the following steps are performed:
481: 1) compute eigenvectors of (S,T): S*Z=T*Z*D
482: 2) compute eigenvectors of OP: X=V*Z
483: */
484: PetscErrorCode EPSComputeVectors_XD(EPS eps)
485: {
487: Mat X;
488: PetscBool symm;
491: PetscObjectTypeCompare((PetscObject)eps->ds,DSHEP,&symm);
492: if (symm) return(0);
493: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
495: /* V <- V * X */
496: DSGetMat(eps->ds,DS_MAT_X,&X);
497: BVSetActiveColumns(eps->V,0,eps->nconv);
498: BVMultInPlace(eps->V,X,0,eps->nconv);
499: DSRestoreMat(eps->ds,DS_MAT_X,&X);
500: return(0);
501: }