Actual source code: davidson.c
slepc-3.22.2 2024-12-02
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: EPSCheckNotStructured(eps);
50: /* Setup EPS options and get the problem specification */
51: bs = data->blocksize;
52: if (bs <= 0) bs = 1;
53: if (eps->ncv!=PETSC_DETERMINE) {
54: PetscCheck(eps->ncv>=eps->nev,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of ncv must be at least nev");
55: } else if (eps->mpd!=PETSC_DETERMINE) eps->ncv = eps->mpd + eps->nev + bs;
56: else if (eps->n < 10) eps->ncv = eps->n+eps->nev+bs;
57: else if (eps->nev < 500) eps->ncv = PetscMax(eps->nev,PetscMin(eps->n-bs,PetscMax(2*eps->nev,eps->nev+15))+bs);
58: else eps->ncv = PetscMax(eps->nev,PetscMin(eps->n-bs,eps->nev+500)+bs);
59: if (eps->mpd==PETSC_DETERMINE) eps->mpd = PetscMin(eps->n,eps->ncv);
60: PetscCheck(eps->mpd<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd parameter has to be less than or equal to ncv");
61: PetscCheck(eps->mpd>=2,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd parameter has to be greater than 2");
62: if (eps->max_it == PETSC_DETERMINE) eps->max_it = PetscMax(100*eps->ncv,2*eps->n);
63: if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
64: PetscCheck(eps->nev+bs<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of ncv has to be greater than nev plus blocksize");
65: PetscCheck(!eps->trueres,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"-eps_true_residual is disabled in this solver.");
66: EPSCheckUnsupported(eps,EPS_FEATURE_REGION | EPS_FEATURE_TWOSIDED);
68: if (!data->minv) data->minv = (eps->n && eps->n<10)? 1: PetscMin(PetscMax(bs,6),eps->mpd/2);
69: min_size_V = data->minv;
70: PetscCheck(min_size_V+bs<=eps->mpd,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of minv must be less than mpd minus blocksize");
71: if (data->plusk == PETSC_DETERMINE) {
72: if (eps->problem_type == EPS_GHIEP || eps->nev+bs>eps->ncv) data->plusk = 0;
73: else data->plusk = 1;
74: }
75: if (!data->initialsize) data->initialsize = (eps->n && eps->n<10)? 1: 6;
76: initv = data->initialsize;
77: PetscCheck(eps->mpd>=initv,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The initv parameter has to be less than or equal to mpd");
79: /* Change the default sigma to inf if necessary */
80: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL || eps->which == EPS_LARGEST_IMAGINARY) PetscCall(STSetDefaultShift(eps->st,PETSC_MAX_REAL));
82: /* Set up preconditioner */
83: PetscCall(STSetUp(eps->st));
85: /* Setup problem specification in dvd */
86: PetscCall(STGetNumMatrices(eps->st,&nmat));
87: PetscCall(STGetMatrix(eps->st,0,&A));
88: if (nmat>1) PetscCall(STGetMatrix(eps->st,1,&B));
89: PetscCall(EPSReset_XD(eps));
90: PetscCall(PetscMemzero(dvd,sizeof(dvdDashboard)));
91: dvd->A = A; dvd->B = eps->isgeneralized? B: NULL;
92: ispositive = eps->ispositive;
93: dvd->sA = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN: 0) | ((ispositive && !eps->isgeneralized) ? DVD_MAT_POS_DEF: 0);
94: /* Assume -eps_hermitian means hermitian-definite in generalized problems */
95: if (!ispositive && !eps->isgeneralized && eps->ishermitian) ispositive = PETSC_TRUE;
96: if (!eps->isgeneralized) dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY | DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
97: else dvd->sB = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN: 0) | (ispositive? DVD_MAT_POS_DEF: 0);
98: ipB = (dvd->B && data->ipB && DVD_IS(dvd->sB,DVD_MAT_HERMITIAN))?PETSC_TRUE:PETSC_FALSE;
99: 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);
100: if (data->ipB && !ipB) data->ipB = PETSC_FALSE;
101: dvd->correctXnorm = (dvd->B && (DVD_IS(dvd->sB,DVD_MAT_HERMITIAN)||DVD_IS(dvd->sEP,DVD_EP_INDEFINITE)))?PETSC_TRUE:PETSC_FALSE;
102: dvd->nev = eps->nev;
103: dvd->which = eps->which;
104: dvd->withTarget = PETSC_TRUE;
105: switch (eps->which) {
106: case EPS_TARGET_MAGNITUDE:
107: case EPS_TARGET_IMAGINARY:
108: dvd->target[0] = target = eps->target;
109: dvd->target[1] = 1.0;
110: break;
111: case EPS_TARGET_REAL:
112: dvd->target[0] = PetscRealPart(target = eps->target);
113: dvd->target[1] = 1.0;
114: break;
115: case EPS_LARGEST_REAL:
116: case EPS_LARGEST_MAGNITUDE:
117: case EPS_LARGEST_IMAGINARY: /* TODO: think about this case */
118: dvd->target[0] = 1.0;
119: dvd->target[1] = target = 0.0;
120: break;
121: case EPS_SMALLEST_MAGNITUDE:
122: case EPS_SMALLEST_REAL:
123: case EPS_SMALLEST_IMAGINARY: /* TODO: think about this case */
124: dvd->target[0] = target = 0.0;
125: dvd->target[1] = 1.0;
126: break;
127: case EPS_WHICH_USER:
128: PetscCall(STGetShift(eps->st,&target));
129: dvd->target[0] = target;
130: dvd->target[1] = 1.0;
131: break;
132: case EPS_ALL:
133: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues");
134: default:
135: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported value of option 'which'");
136: }
137: dvd->tol = SlepcDefaultTol(eps->tol);
138: dvd->eps = eps;
140: /* Setup the extraction technique */
141: if (!eps->extraction) {
142: if (ipB || ispositive) eps->extraction = EPS_RITZ;
143: else {
144: switch (eps->which) {
145: case EPS_TARGET_REAL:
146: case EPS_TARGET_MAGNITUDE:
147: case EPS_TARGET_IMAGINARY:
148: case EPS_SMALLEST_MAGNITUDE:
149: case EPS_SMALLEST_REAL:
150: case EPS_SMALLEST_IMAGINARY:
151: eps->extraction = EPS_HARMONIC;
152: break;
153: case EPS_LARGEST_REAL:
154: case EPS_LARGEST_MAGNITUDE:
155: case EPS_LARGEST_IMAGINARY:
156: eps->extraction = EPS_HARMONIC_LARGEST;
157: break;
158: default:
159: eps->extraction = EPS_RITZ;
160: }
161: }
162: }
163: switch (eps->extraction) {
164: case EPS_RITZ: harm = DVD_HARM_NONE; break;
165: case EPS_HARMONIC: harm = DVD_HARM_RR; break;
166: case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
167: case EPS_HARMONIC_RIGHT: harm = DVD_HARM_REIGS; break;
168: case EPS_HARMONIC_LARGEST: harm = DVD_HARM_LEIGS; break;
169: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
170: }
172: /* Setup the type of starting subspace */
173: init = data->krylovstart? DVD_INITV_KRYLOV: DVD_INITV_CLASSIC;
175: /* Preconfigure dvd */
176: PetscCall(STGetKSP(eps->st,&ksp));
177: 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));
179: /* Allocate memory */
180: PetscCall(EPSAllocateSolution(eps,0));
182: /* Setup orthogonalization */
183: PetscCall(EPS_SetInnerProduct(eps));
184: if (!(ipB && dvd->B)) PetscCall(BVSetMatrix(eps->V,NULL,PETSC_FALSE));
186: /* Configure dvd for a basic GD */
187: 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));
188: PetscFunctionReturn(PETSC_SUCCESS);
189: }
191: PetscErrorCode EPSSolve_XD(EPS eps)
192: {
193: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
194: dvdDashboard *d = &data->ddb;
195: PetscInt l,k;
197: PetscFunctionBegin;
198: PetscCall(PetscCitationsRegister(citation,&cited));
199: /* Call the starting routines */
200: PetscCall(EPSDavidsonFLCall(d->startList,d));
202: while (eps->reason == EPS_CONVERGED_ITERATING) {
204: /* Initialize V, if it is needed */
205: PetscCall(BVGetActiveColumns(d->eps->V,&l,&k));
206: if (PetscUnlikely(l == k)) PetscCall(d->initV(d));
208: /* Find the best approximated eigenpairs in V, X */
209: PetscCall(d->calcPairs(d));
211: /* Test for convergence */
212: PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
213: if (eps->reason != EPS_CONVERGED_ITERATING) break;
215: /* Expand the subspace */
216: PetscCall(d->updateV(d));
218: /* Monitor progress */
219: eps->nconv = d->nconv;
220: eps->its++;
221: PetscCall(BVGetActiveColumns(d->eps->V,NULL,&k));
222: PetscCall(EPSMonitor(eps,eps->its,eps->nconv+d->npreconv,eps->eigr,eps->eigi,eps->errest,PetscMin(k,eps->nev)));
223: }
225: /* Call the ending routines */
226: PetscCall(EPSDavidsonFLCall(d->endList,d));
227: PetscFunctionReturn(PETSC_SUCCESS);
228: }
230: PetscErrorCode EPSReset_XD(EPS eps)
231: {
232: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
233: dvdDashboard *dvd = &data->ddb;
235: PetscFunctionBegin;
236: /* Call step destructors and destroys the list */
237: PetscCall(EPSDavidsonFLCall(dvd->destroyList,dvd));
238: PetscCall(EPSDavidsonFLDestroy(&dvd->destroyList));
239: PetscCall(EPSDavidsonFLDestroy(&dvd->startList));
240: PetscCall(EPSDavidsonFLDestroy(&dvd->endList));
241: PetscFunctionReturn(PETSC_SUCCESS);
242: }
244: PetscErrorCode EPSXDSetKrylovStart_XD(EPS eps,PetscBool krylovstart)
245: {
246: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
248: PetscFunctionBegin;
249: data->krylovstart = krylovstart;
250: PetscFunctionReturn(PETSC_SUCCESS);
251: }
253: PetscErrorCode EPSXDGetKrylovStart_XD(EPS eps,PetscBool *krylovstart)
254: {
255: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
257: PetscFunctionBegin;
258: *krylovstart = data->krylovstart;
259: PetscFunctionReturn(PETSC_SUCCESS);
260: }
262: PetscErrorCode EPSXDSetBlockSize_XD(EPS eps,PetscInt blocksize)
263: {
264: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
266: PetscFunctionBegin;
267: if (blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
268: PetscCheck(blocksize>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value, must be >0");
269: if (data->blocksize != blocksize) {
270: data->blocksize = blocksize;
271: eps->state = EPS_STATE_INITIAL;
272: }
273: PetscFunctionReturn(PETSC_SUCCESS);
274: }
276: PetscErrorCode EPSXDGetBlockSize_XD(EPS eps,PetscInt *blocksize)
277: {
278: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
280: PetscFunctionBegin;
281: *blocksize = data->blocksize;
282: PetscFunctionReturn(PETSC_SUCCESS);
283: }
285: PetscErrorCode EPSXDSetRestart_XD(EPS eps,PetscInt minv,PetscInt plusk)
286: {
287: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
289: PetscFunctionBegin;
290: if (minv == PETSC_DETERMINE) {
291: if (data->minv != 0) eps->state = EPS_STATE_INITIAL;
292: data->minv = 0;
293: } else if (minv != PETSC_CURRENT) {
294: PetscCheck(minv>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value, must be >0");
295: if (data->minv != minv) eps->state = EPS_STATE_INITIAL;
296: data->minv = minv;
297: }
298: if (plusk == PETSC_DETERMINE) {
299: if (data->plusk != PETSC_DETERMINE) eps->state = EPS_STATE_INITIAL;
300: data->plusk = PETSC_DETERMINE;
301: } else if (plusk != PETSC_CURRENT) {
302: PetscCheck(plusk>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value, must be >0");
303: if (data->plusk != plusk) eps->state = EPS_STATE_INITIAL;
304: data->plusk = plusk;
305: }
306: PetscFunctionReturn(PETSC_SUCCESS);
307: }
309: PetscErrorCode EPSXDGetRestart_XD(EPS eps,PetscInt *minv,PetscInt *plusk)
310: {
311: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
313: PetscFunctionBegin;
314: if (minv) *minv = data->minv;
315: if (plusk) *plusk = data->plusk;
316: PetscFunctionReturn(PETSC_SUCCESS);
317: }
319: PetscErrorCode EPSXDGetInitialSize_XD(EPS eps,PetscInt *initialsize)
320: {
321: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
323: PetscFunctionBegin;
324: *initialsize = data->initialsize;
325: PetscFunctionReturn(PETSC_SUCCESS);
326: }
328: PetscErrorCode EPSXDSetInitialSize_XD(EPS eps,PetscInt initialsize)
329: {
330: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
332: PetscFunctionBegin;
333: if (initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 0;
334: else PetscCheck(initialsize>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value, must be >0");
335: if (data->initialsize != initialsize) {
336: data->initialsize = initialsize;
337: eps->state = EPS_STATE_INITIAL;
338: }
339: PetscFunctionReturn(PETSC_SUCCESS);
340: }
342: PetscErrorCode EPSXDSetBOrth_XD(EPS eps,PetscBool borth)
343: {
344: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
346: PetscFunctionBegin;
347: data->ipB = borth;
348: PetscFunctionReturn(PETSC_SUCCESS);
349: }
351: PetscErrorCode EPSXDGetBOrth_XD(EPS eps,PetscBool *borth)
352: {
353: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
355: PetscFunctionBegin;
356: *borth = data->ipB;
357: PetscFunctionReturn(PETSC_SUCCESS);
358: }
360: /*
361: EPSComputeVectors_XD - Compute eigenvectors from the vectors
362: provided by the eigensolver. This version is intended for solvers
363: that provide Schur vectors from the QZ decomposition. Given the partial
364: Schur decomposition OP*V=V*T, the following steps are performed:
365: 1) compute eigenvectors of (S,T): S*Z=T*Z*D
366: 2) compute eigenvectors of OP: X=V*Z
367: */
368: PetscErrorCode EPSComputeVectors_XD(EPS eps)
369: {
370: Mat X;
371: PetscBool symm;
373: PetscFunctionBegin;
374: PetscCall(PetscObjectTypeCompare((PetscObject)eps->ds,DSHEP,&symm));
375: if (symm) PetscFunctionReturn(PETSC_SUCCESS);
376: PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
378: /* V <- V * X */
379: PetscCall(DSGetMat(eps->ds,DS_MAT_X,&X));
380: PetscCall(BVSetActiveColumns(eps->V,0,eps->nconv));
381: PetscCall(BVMultInPlace(eps->V,X,0,eps->nconv));
382: PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&X));
383: PetscFunctionReturn(PETSC_SUCCESS);
384: }