Actual source code: davidson.c
slepc-3.21.1 2024-04-26
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: Skeleton of Davidson solver. Actual solvers are GD and JD.
13: References:
15: [1] E. Romero and J.E. Roman, "A parallel implementation of Davidson
16: methods for large-scale eigenvalue problems in SLEPc", ACM Trans.
17: Math. Software 40(2):13, 2014.
18: */
20: #include "davidson.h"
22: static PetscBool cited = PETSC_FALSE;
23: static const char citation[] =
24: "@Article{slepc-davidson,\n"
25: " author = \"E. Romero and J. E. Roman\",\n"
26: " title = \"A parallel implementation of {Davidson} methods for large-scale eigenvalue problems in {SLEPc}\",\n"
27: " journal = \"{ACM} Trans. Math. Software\",\n"
28: " volume = \"40\",\n"
29: " number = \"2\",\n"
30: " pages = \"13:1--13:29\",\n"
31: " year = \"2014,\"\n"
32: " doi = \"https://doi.org/10.1145/2543696\"\n"
33: "}\n";
35: PetscErrorCode EPSSetUp_XD(EPS eps)
36: {
37: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
38: dvdDashboard *dvd = &data->ddb;
39: dvdBlackboard b;
40: PetscInt min_size_V,bs,initv,nmat;
41: Mat A,B;
42: KSP ksp;
43: PetscBool ipB,ispositive;
44: HarmType_t harm;
45: InitType_t init;
46: PetscScalar target;
48: PetscFunctionBegin;
49: /* Setup EPS options and get the problem specification */
50: bs = data->blocksize;
51: if (bs <= 0) bs = 1;
52: if (eps->ncv!=PETSC_DEFAULT && eps->ncv!=PETSC_DECIDE) {
53: PetscCheck(eps->ncv>=eps->nev,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of ncv must be at least nev");
54: } else if (eps->mpd!=PETSC_DEFAULT && eps->mpd!=PETSC_DECIDE) eps->ncv = eps->mpd + eps->nev + bs;
55: else if (eps->n < 10) eps->ncv = eps->n+eps->nev+bs;
56: else if (eps->nev < 500) eps->ncv = PetscMax(eps->nev,PetscMin(eps->n-bs,PetscMax(2*eps->nev,eps->nev+15))+bs);
57: else eps->ncv = PetscMax(eps->nev,PetscMin(eps->n-bs,eps->nev+500)+bs);
58: if (eps->mpd==PETSC_DEFAULT || eps->mpd==PETSC_DECIDE) eps->mpd = PetscMin(eps->n,eps->ncv);
59: PetscCheck(eps->mpd<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd parameter has to be less than or equal to ncv");
60: PetscCheck(eps->mpd>=2,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd parameter has to be greater than 2");
61: if (eps->max_it == PETSC_DEFAULT || eps->max_it == PETSC_DECIDE) eps->max_it = PetscMax(100*eps->ncv,2*eps->n);
62: if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
63: PetscCheck(eps->nev+bs<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of ncv has to be greater than nev plus blocksize");
64: PetscCheck(!eps->trueres,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"-eps_true_residual is disabled in this solver.");
65: EPSCheckUnsupported(eps,EPS_FEATURE_REGION | EPS_FEATURE_TWOSIDED);
67: if (!data->minv) data->minv = (eps->n && eps->n<10)? 1: PetscMin(PetscMax(bs,6),eps->mpd/2);
68: min_size_V = data->minv;
69: PetscCheck(min_size_V+bs<=eps->mpd,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of minv must be less than mpd minus blocksize");
70: if (data->plusk == PETSC_DEFAULT) {
71: if (eps->problem_type == EPS_GHIEP || eps->nev+bs>eps->ncv) data->plusk = 0;
72: else data->plusk = 1;
73: }
74: if (!data->initialsize) data->initialsize = (eps->n && eps->n<10)? 1: 6;
75: initv = data->initialsize;
76: PetscCheck(eps->mpd>=initv,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The initv parameter has to be less than or equal to mpd");
78: /* Change the default sigma to inf if necessary */
79: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL || eps->which == EPS_LARGEST_IMAGINARY) PetscCall(STSetDefaultShift(eps->st,PETSC_MAX_REAL));
81: /* Set up preconditioner */
82: PetscCall(STSetUp(eps->st));
84: /* Setup problem specification in dvd */
85: PetscCall(STGetNumMatrices(eps->st,&nmat));
86: PetscCall(STGetMatrix(eps->st,0,&A));
87: if (nmat>1) PetscCall(STGetMatrix(eps->st,1,&B));
88: PetscCall(EPSReset_XD(eps));
89: PetscCall(PetscMemzero(dvd,sizeof(dvdDashboard)));
90: dvd->A = A; dvd->B = eps->isgeneralized? B: NULL;
91: ispositive = eps->ispositive;
92: dvd->sA = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN: 0) | ((ispositive && !eps->isgeneralized) ? DVD_MAT_POS_DEF: 0);
93: /* Assume -eps_hermitian means hermitian-definite in generalized problems */
94: if (!ispositive && !eps->isgeneralized && eps->ishermitian) ispositive = PETSC_TRUE;
95: if (!eps->isgeneralized) dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY | DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
96: else dvd->sB = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN: 0) | (ispositive? DVD_MAT_POS_DEF: 0);
97: ipB = (dvd->B && data->ipB && DVD_IS(dvd->sB,DVD_MAT_HERMITIAN))?PETSC_TRUE:PETSC_FALSE;
98: 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);
99: if (data->ipB && !ipB) data->ipB = PETSC_FALSE;
100: dvd->correctXnorm = (dvd->B && (DVD_IS(dvd->sB,DVD_MAT_HERMITIAN)||DVD_IS(dvd->sEP,DVD_EP_INDEFINITE)))?PETSC_TRUE:PETSC_FALSE;
101: dvd->nev = eps->nev;
102: dvd->which = eps->which;
103: dvd->withTarget = PETSC_TRUE;
104: switch (eps->which) {
105: case EPS_TARGET_MAGNITUDE:
106: case EPS_TARGET_IMAGINARY:
107: dvd->target[0] = target = eps->target;
108: dvd->target[1] = 1.0;
109: break;
110: case EPS_TARGET_REAL:
111: dvd->target[0] = PetscRealPart(target = eps->target);
112: dvd->target[1] = 1.0;
113: break;
114: case EPS_LARGEST_REAL:
115: case EPS_LARGEST_MAGNITUDE:
116: case EPS_LARGEST_IMAGINARY: /* TODO: think about this case */
117: dvd->target[0] = 1.0;
118: dvd->target[1] = target = 0.0;
119: break;
120: case EPS_SMALLEST_MAGNITUDE:
121: case EPS_SMALLEST_REAL:
122: case EPS_SMALLEST_IMAGINARY: /* TODO: think about this case */
123: dvd->target[0] = target = 0.0;
124: dvd->target[1] = 1.0;
125: break;
126: case EPS_WHICH_USER:
127: PetscCall(STGetShift(eps->st,&target));
128: dvd->target[0] = target;
129: dvd->target[1] = 1.0;
130: break;
131: case EPS_ALL:
132: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues");
133: default:
134: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported value of option 'which'");
135: }
136: dvd->tol = SlepcDefaultTol(eps->tol);
137: dvd->eps = eps;
139: /* Setup the extraction technique */
140: if (!eps->extraction) {
141: if (ipB || ispositive) eps->extraction = EPS_RITZ;
142: else {
143: switch (eps->which) {
144: case EPS_TARGET_REAL:
145: case EPS_TARGET_MAGNITUDE:
146: case EPS_TARGET_IMAGINARY:
147: case EPS_SMALLEST_MAGNITUDE:
148: case EPS_SMALLEST_REAL:
149: case EPS_SMALLEST_IMAGINARY:
150: eps->extraction = EPS_HARMONIC;
151: break;
152: case EPS_LARGEST_REAL:
153: case EPS_LARGEST_MAGNITUDE:
154: case EPS_LARGEST_IMAGINARY:
155: eps->extraction = EPS_HARMONIC_LARGEST;
156: break;
157: default:
158: eps->extraction = EPS_RITZ;
159: }
160: }
161: }
162: switch (eps->extraction) {
163: case EPS_RITZ: harm = DVD_HARM_NONE; break;
164: case EPS_HARMONIC: harm = DVD_HARM_RR; break;
165: case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
166: case EPS_HARMONIC_RIGHT: harm = DVD_HARM_REIGS; break;
167: case EPS_HARMONIC_LARGEST: harm = DVD_HARM_LEIGS; break;
168: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
169: }
171: /* Setup the type of starting subspace */
172: init = data->krylovstart? DVD_INITV_KRYLOV: DVD_INITV_CLASSIC;
174: /* Preconfigure dvd */
175: PetscCall(STGetKSP(eps->st,&ksp));
176: PetscCall(dvd_schm_basic_preconf(dvd,&b,eps->mpd,min_size_V,bs,initv,PetscAbs(eps->nini),data->plusk,harm,ksp,init,eps->trackall,data->ipB,data->doubleexp));
178: /* Allocate memory */
179: PetscCall(EPSAllocateSolution(eps,0));
181: /* Setup orthogonalization */
182: PetscCall(EPS_SetInnerProduct(eps));
183: if (!(ipB && dvd->B)) PetscCall(BVSetMatrix(eps->V,NULL,PETSC_FALSE));
185: /* Configure dvd for a basic GD */
186: PetscCall(dvd_schm_basic_conf(dvd,&b,eps->mpd,min_size_V,bs,initv,PetscAbs(eps->nini),data->plusk,harm,dvd->withTarget,target,ksp,data->fix,init,eps->trackall,data->ipB,data->dynamic,data->doubleexp));
187: PetscFunctionReturn(PETSC_SUCCESS);
188: }
190: PetscErrorCode EPSSolve_XD(EPS eps)
191: {
192: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
193: dvdDashboard *d = &data->ddb;
194: PetscInt l,k;
196: PetscFunctionBegin;
197: PetscCall(PetscCitationsRegister(citation,&cited));
198: /* Call the starting routines */
199: PetscCall(EPSDavidsonFLCall(d->startList,d));
201: while (eps->reason == EPS_CONVERGED_ITERATING) {
203: /* Initialize V, if it is needed */
204: PetscCall(BVGetActiveColumns(d->eps->V,&l,&k));
205: if (PetscUnlikely(l == k)) PetscCall(d->initV(d));
207: /* Find the best approximated eigenpairs in V, X */
208: PetscCall(d->calcPairs(d));
210: /* Test for convergence */
211: PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
212: if (eps->reason != EPS_CONVERGED_ITERATING) break;
214: /* Expand the subspace */
215: PetscCall(d->updateV(d));
217: /* Monitor progress */
218: eps->nconv = d->nconv;
219: eps->its++;
220: PetscCall(BVGetActiveColumns(d->eps->V,NULL,&k));
221: PetscCall(EPSMonitor(eps,eps->its,eps->nconv+d->npreconv,eps->eigr,eps->eigi,eps->errest,PetscMin(k,eps->nev)));
222: }
224: /* Call the ending routines */
225: PetscCall(EPSDavidsonFLCall(d->endList,d));
226: PetscFunctionReturn(PETSC_SUCCESS);
227: }
229: PetscErrorCode EPSReset_XD(EPS eps)
230: {
231: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
232: dvdDashboard *dvd = &data->ddb;
234: PetscFunctionBegin;
235: /* Call step destructors and destroys the list */
236: PetscCall(EPSDavidsonFLCall(dvd->destroyList,dvd));
237: PetscCall(EPSDavidsonFLDestroy(&dvd->destroyList));
238: PetscCall(EPSDavidsonFLDestroy(&dvd->startList));
239: PetscCall(EPSDavidsonFLDestroy(&dvd->endList));
240: PetscFunctionReturn(PETSC_SUCCESS);
241: }
243: PetscErrorCode EPSXDSetKrylovStart_XD(EPS eps,PetscBool krylovstart)
244: {
245: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
247: PetscFunctionBegin;
248: data->krylovstart = krylovstart;
249: PetscFunctionReturn(PETSC_SUCCESS);
250: }
252: PetscErrorCode EPSXDGetKrylovStart_XD(EPS eps,PetscBool *krylovstart)
253: {
254: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
256: PetscFunctionBegin;
257: *krylovstart = data->krylovstart;
258: PetscFunctionReturn(PETSC_SUCCESS);
259: }
261: PetscErrorCode EPSXDSetBlockSize_XD(EPS eps,PetscInt blocksize)
262: {
263: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
265: PetscFunctionBegin;
266: if (blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
267: PetscCheck(blocksize>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value, must be >0");
268: if (data->blocksize != blocksize) {
269: data->blocksize = blocksize;
270: eps->state = EPS_STATE_INITIAL;
271: }
272: PetscFunctionReturn(PETSC_SUCCESS);
273: }
275: PetscErrorCode EPSXDGetBlockSize_XD(EPS eps,PetscInt *blocksize)
276: {
277: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
279: PetscFunctionBegin;
280: *blocksize = data->blocksize;
281: PetscFunctionReturn(PETSC_SUCCESS);
282: }
284: PetscErrorCode EPSXDSetRestart_XD(EPS eps,PetscInt minv,PetscInt plusk)
285: {
286: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
288: PetscFunctionBegin;
289: if (minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 0;
290: else PetscCheck(minv>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value, must be >0");
291: if (plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = PETSC_DEFAULT;
292: else PetscCheck(plusk>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value, must be >0");
293: if (data->minv != minv || data->plusk != plusk) {
294: data->minv = minv;
295: data->plusk = plusk;
296: eps->state = EPS_STATE_INITIAL;
297: }
298: PetscFunctionReturn(PETSC_SUCCESS);
299: }
301: PetscErrorCode EPSXDGetRestart_XD(EPS eps,PetscInt *minv,PetscInt *plusk)
302: {
303: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
305: PetscFunctionBegin;
306: if (minv) *minv = data->minv;
307: if (plusk) *plusk = data->plusk;
308: PetscFunctionReturn(PETSC_SUCCESS);
309: }
311: PetscErrorCode EPSXDGetInitialSize_XD(EPS eps,PetscInt *initialsize)
312: {
313: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
315: PetscFunctionBegin;
316: *initialsize = data->initialsize;
317: PetscFunctionReturn(PETSC_SUCCESS);
318: }
320: PetscErrorCode EPSXDSetInitialSize_XD(EPS eps,PetscInt initialsize)
321: {
322: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
324: PetscFunctionBegin;
325: if (initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 0;
326: else PetscCheck(initialsize>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value, must be >0");
327: if (data->initialsize != initialsize) {
328: data->initialsize = initialsize;
329: eps->state = EPS_STATE_INITIAL;
330: }
331: PetscFunctionReturn(PETSC_SUCCESS);
332: }
334: PetscErrorCode EPSXDSetBOrth_XD(EPS eps,PetscBool borth)
335: {
336: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
338: PetscFunctionBegin;
339: data->ipB = borth;
340: PetscFunctionReturn(PETSC_SUCCESS);
341: }
343: PetscErrorCode EPSXDGetBOrth_XD(EPS eps,PetscBool *borth)
344: {
345: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
347: PetscFunctionBegin;
348: *borth = data->ipB;
349: PetscFunctionReturn(PETSC_SUCCESS);
350: }
352: /*
353: EPSComputeVectors_XD - Compute eigenvectors from the vectors
354: provided by the eigensolver. This version is intended for solvers
355: that provide Schur vectors from the QZ decomposition. Given the partial
356: Schur decomposition OP*V=V*T, the following steps are performed:
357: 1) compute eigenvectors of (S,T): S*Z=T*Z*D
358: 2) compute eigenvectors of OP: X=V*Z
359: */
360: PetscErrorCode EPSComputeVectors_XD(EPS eps)
361: {
362: Mat X;
363: PetscBool symm;
365: PetscFunctionBegin;
366: PetscCall(PetscObjectTypeCompare((PetscObject)eps->ds,DSHEP,&symm));
367: if (symm) PetscFunctionReturn(PETSC_SUCCESS);
368: PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
370: /* V <- V * X */
371: PetscCall(DSGetMat(eps->ds,DS_MAT_X,&X));
372: PetscCall(BVSetActiveColumns(eps->V,0,eps->nconv));
373: PetscCall(BVMultInPlace(eps->V,X,0,eps->nconv));
374: PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&X));
375: PetscFunctionReturn(PETSC_SUCCESS);
376: }