Actual source code: subspace.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: SLEPc eigensolver: "subspace"
13: Method: Subspace Iteration
15: Algorithm:
17: Subspace iteration with Rayleigh-Ritz projection and locking,
18: based on the SRRIT implementation.
20: References:
22: [1] "Subspace Iteration in SLEPc", SLEPc Technical Report STR-3,
23: available at https://slepc.upv.es.
24: */
26: #include <slepc/private/epsimpl.h>
28: typedef struct {
29: PetscBool estimatedrange; /* the filter range was not set by the user */
30: } EPS_SUBSPACE;
32: static PetscErrorCode EPSSetUp_Subspace_Filter(EPS eps)
33: {
34: EPS_SUBSPACE *ctx = (EPS_SUBSPACE*)eps->data;
35: PetscBool estimaterange=PETSC_TRUE;
36: PetscReal rleft,rright;
37: Mat A;
39: PetscFunctionBegin;
40: EPSCheckHermitianCondition(eps,PETSC_TRUE," with polynomial filter");
41: EPSCheckStandardCondition(eps,PETSC_TRUE," with polynomial filter");
42: PetscCheck(eps->intb<PETSC_MAX_REAL || eps->inta>PETSC_MIN_REAL,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The defined computational interval should have at least one of their sides bounded");
43: EPSCheckUnsupportedCondition(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_EXTRACTION,PETSC_TRUE," with polynomial filter");
44: PetscCall(STFilterSetInterval(eps->st,eps->inta,eps->intb));
45: if (!ctx->estimatedrange) {
46: PetscCall(STFilterGetRange(eps->st,&rleft,&rright));
47: estimaterange = (!rleft && !rright)? PETSC_TRUE: PETSC_FALSE;
48: }
49: if (estimaterange) { /* user did not set a range */
50: PetscCall(STGetMatrix(eps->st,0,&A));
51: PetscCall(MatEstimateSpectralRange_EPS(A,&rleft,&rright));
52: PetscCall(PetscInfo(eps,"Setting eigenvalue range to [%g,%g]\n",(double)rleft,(double)rright));
53: PetscCall(STFilterSetRange(eps->st,rleft,rright));
54: ctx->estimatedrange = PETSC_TRUE;
55: }
56: if (eps->ncv==PETSC_DEFAULT && eps->nev==1) eps->nev = 40; /* user did not provide nev estimation */
57: PetscCall(EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd));
58: PetscCheck(eps->ncv<=eps->nev+eps->mpd,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must not be larger than nev+mpd");
59: PetscFunctionReturn(PETSC_SUCCESS);
60: }
62: static PetscErrorCode EPSSetUp_Subspace(EPS eps)
63: {
64: PetscBool isfilt;
66: PetscFunctionBegin;
67: EPSCheckDefinite(eps);
68: if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
69: if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
70: if (eps->which==EPS_ALL) {
71: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt));
72: PetscCheck(isfilt,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing not supported in this solver");
73: PetscCall(EPSSetUp_Subspace_Filter(eps));
74: } else {
75: PetscCheck(eps->which==EPS_LARGEST_MAGNITUDE || eps->which==EPS_TARGET_MAGNITUDE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only largest magnitude or target magnitude eigenvalues");
76: PetscCall(EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd));
77: }
78: EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_EXTRACTION | EPS_FEATURE_TWOSIDED);
79: PetscCheck(eps->converged==EPSConvergedRelative,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver only supports relative convergence test");
81: PetscCall(EPSAllocateSolution(eps,0));
82: PetscCall(EPS_SetInnerProduct(eps));
83: if (eps->ishermitian) PetscCall(DSSetType(eps->ds,DSHEP));
84: else PetscCall(DSSetType(eps->ds,DSNHEP));
85: PetscCall(DSAllocate(eps->ds,eps->ncv));
86: PetscFunctionReturn(PETSC_SUCCESS);
87: }
89: static PetscErrorCode EPSSetUpSort_Subspace(EPS eps)
90: {
91: SlepcSC sc;
93: PetscFunctionBegin;
94: PetscCall(EPSSetUpSort_Default(eps));
95: if (eps->which==EPS_ALL) {
96: PetscCall(DSGetSlepcSC(eps->ds,&sc));
97: sc->rg = NULL;
98: sc->comparison = SlepcCompareLargestReal;
99: sc->comparisonctx = NULL;
100: sc->map = NULL;
101: sc->mapobj = NULL;
102: }
103: PetscFunctionReturn(PETSC_SUCCESS);
104: }
106: /*
107: EPSSubspaceFindGroup - Find a group of nearly equimodular eigenvalues, provided
108: in arrays wr and wi, according to the tolerance grptol. Also the 2-norms
109: of the residuals must be passed in (rsd). Arrays are processed from index
110: l to index m only. The output information is:
112: ngrp - number of entries of the group
113: ctr - (w(l)+w(l+ngrp-1))/2
114: ae - average of wr(l),...,wr(l+ngrp-1)
115: arsd - average of rsd(l),...,rsd(l+ngrp-1)
116: */
117: static PetscErrorCode EPSSubspaceFindGroup(PetscInt l,PetscInt m,PetscScalar *wr,PetscScalar *wi,PetscReal *rsd,PetscReal grptol,PetscInt *ngrp,PetscReal *ctr,PetscReal *ae,PetscReal *arsd)
118: {
119: PetscInt i;
120: PetscReal rmod,rmod1;
122: PetscFunctionBegin;
123: *ngrp = 0;
124: *ctr = 0;
125: rmod = SlepcAbsEigenvalue(wr[l],wi[l]);
127: for (i=l;i<m;) {
128: rmod1 = SlepcAbsEigenvalue(wr[i],wi[i]);
129: if (PetscAbsReal(rmod-rmod1) > grptol*(rmod+rmod1)) break;
130: *ctr = (rmod+rmod1)/2.0;
131: if (wi[i] == 0.0) {
132: (*ngrp)++;
133: i++;
134: } else {
135: (*ngrp)+=2;
136: i+=2;
137: }
138: }
140: *ae = 0;
141: *arsd = 0;
142: if (*ngrp) {
143: for (i=l;i<l+*ngrp;i++) {
144: (*ae) += PetscRealPart(wr[i]);
145: (*arsd) += rsd[i]*rsd[i];
146: }
147: *ae = *ae / *ngrp;
148: *arsd = PetscSqrtReal(*arsd / *ngrp);
149: }
150: PetscFunctionReturn(PETSC_SUCCESS);
151: }
153: /*
154: EPSSubspaceResidualNorms - Computes the column norms of residual vectors
155: OP*V(1:n,l:m) - V*T(1:m,l:m), where, on entry, OP*V has been computed and
156: stored in R. On exit, rsd(l) to rsd(m) contain the computed norms.
157: */
158: static PetscErrorCode EPSSubspaceResidualNorms(BV R,BV V,Mat T,PetscInt l,PetscInt m,PetscScalar *eigi,PetscReal *rsd)
159: {
160: PetscInt i;
162: PetscFunctionBegin;
163: PetscCall(BVMult(R,-1.0,1.0,V,T));
164: for (i=l;i<m;i++) PetscCall(BVNormColumnBegin(R,i,NORM_2,rsd+i));
165: for (i=l;i<m;i++) PetscCall(BVNormColumnEnd(R,i,NORM_2,rsd+i));
166: #if !defined(PETSC_USE_COMPLEX)
167: for (i=l;i<m-1;i++) {
168: if (eigi[i]!=0.0) {
169: rsd[i] = SlepcAbs(rsd[i],rsd[i+1]);
170: rsd[i+1] = rsd[i];
171: i++;
172: }
173: }
174: #endif
175: PetscFunctionReturn(PETSC_SUCCESS);
176: }
178: static PetscErrorCode EPSSolve_Subspace(EPS eps)
179: {
180: Mat H,Q,S,T,B;
181: BV AV,R;
182: PetscBool indef;
183: PetscInt i,k,ld,ngrp,nogrp,*itrsd,*itrsdold;
184: PetscInt nxtsrr,idsrr,idort,nxtort,nv,ncv = eps->ncv,its,ninside;
185: PetscReal arsd,oarsd,ctr,octr,ae,oae,*rsd,*orsd,tcond=1.0,gamma;
186: PetscScalar *oeigr,*oeigi;
187: /* Parameters */
188: PetscInt init = 5; /* Number of initial iterations */
189: PetscReal stpfac = 1.5; /* Max num of iter before next SRR step */
190: PetscReal alpha = 1.0; /* Used to predict convergence of next residual */
191: PetscReal beta = 1.1; /* Used to predict convergence of next residual */
192: PetscReal grptol = SLEPC_DEFAULT_TOL; /* Tolerance for EPSSubspaceFindGroup */
193: PetscReal cnvtol = 1e-6; /* Convergence criterion for cnv */
194: PetscInt orttol = 2; /* Number of decimal digits whose loss
195: can be tolerated in orthogonalization */
197: PetscFunctionBegin;
198: its = 0;
199: PetscCall(PetscMalloc6(ncv,&rsd,ncv,&orsd,ncv,&oeigr,ncv,&oeigi,ncv,&itrsd,ncv,&itrsdold));
200: PetscCall(DSGetLeadingDimension(eps->ds,&ld));
201: PetscCall(BVDuplicate(eps->V,&AV));
202: PetscCall(BVDuplicate(eps->V,&R));
203: PetscCall(STGetOperator(eps->st,&S));
205: for (i=0;i<ncv;i++) {
206: rsd[i] = 0.0;
207: itrsd[i] = -1;
208: }
210: /* Complete the initial basis with random vectors and orthonormalize them */
211: for (k=eps->nini;k<ncv;k++) {
212: PetscCall(BVSetRandomColumn(eps->V,k));
213: PetscCall(BVOrthonormalizeColumn(eps->V,k,PETSC_TRUE,NULL,NULL));
214: }
216: while (eps->reason == EPS_CONVERGED_ITERATING) {
217: eps->its++;
218: nv = PetscMin(eps->nconv+eps->mpd,ncv);
219: PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,0));
221: for (i=eps->nconv;i<nv;i++) {
222: oeigr[i] = eps->eigr[i];
223: oeigi[i] = eps->eigi[i];
224: orsd[i] = rsd[i];
225: }
227: /* AV(:,idx) = OP * V(:,idx) */
228: PetscCall(BVSetActiveColumns(eps->V,eps->nconv,nv));
229: PetscCall(BVSetActiveColumns(AV,eps->nconv,nv));
230: PetscCall(BVMatMult(eps->V,S,AV));
232: /* T(:,idx) = V' * AV(:,idx) */
233: PetscCall(BVSetActiveColumns(eps->V,0,nv));
234: PetscCall(DSGetMat(eps->ds,DS_MAT_A,&H));
235: PetscCall(BVDot(AV,eps->V,H));
236: PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&H));
237: PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
239: /* Solve projected problem */
240: PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi));
241: PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL));
242: PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));
244: /* Update vectors V(:,idx) = V * U(:,idx) */
245: PetscCall(DSGetMat(eps->ds,DS_MAT_Q,&Q));
246: PetscCall(BVSetActiveColumns(AV,0,nv));
247: PetscCall(BVSetActiveColumns(R,0,nv));
248: PetscCall(BVMultInPlace(eps->V,Q,eps->nconv,nv));
249: PetscCall(BVMultInPlace(AV,Q,eps->nconv,nv));
250: PetscCall(DSRestoreMat(eps->ds,DS_MAT_Q,&Q));
251: PetscCall(BVCopy(AV,R));
253: /* Convergence check */
254: PetscCall(DSGetMat(eps->ds,DS_MAT_A,&T));
255: PetscCall(EPSSubspaceResidualNorms(R,eps->V,T,eps->nconv,nv,eps->eigi,rsd));
256: PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&T));
258: if (eps->which==EPS_ALL && eps->its>1) { /* adjust eigenvalue count */
259: ninside = 0;
260: PetscCall(STFilterGetThreshold(eps->st,&gamma));
261: for (i=eps->nconv;i<nv;i++) {
262: if (PetscRealPart(eps->eigr[i]) < gamma) break;
263: ninside++;
264: }
265: eps->nev = eps->nconv+ninside;
266: }
267: for (i=eps->nconv;i<nv;i++) {
268: itrsdold[i] = itrsd[i];
269: itrsd[i] = its;
270: eps->errest[i] = rsd[i];
271: }
273: for (;;) {
274: /* Find clusters of computed eigenvalues */
275: PetscCall(EPSSubspaceFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,eps->errest,grptol,&ngrp,&ctr,&ae,&arsd));
276: PetscCall(EPSSubspaceFindGroup(eps->nconv,nv,oeigr,oeigi,orsd,grptol,&nogrp,&octr,&oae,&oarsd));
277: if (ngrp!=nogrp) break;
278: if (ngrp==0) break;
279: if (PetscAbsReal(ae-oae)>ctr*cnvtol*(itrsd[eps->nconv]-itrsdold[eps->nconv])) break;
280: if (arsd>ctr*eps->tol) break;
281: eps->nconv = eps->nconv + ngrp;
282: if (eps->nconv>=nv) break;
283: }
285: PetscCall(EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv));
286: PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
287: if (eps->reason != EPS_CONVERGED_ITERATING) break;
289: /* Compute nxtsrr (iteration of next projection step) */
290: nxtsrr = PetscMin(eps->max_it,PetscMax((PetscInt)PetscFloorReal(stpfac*its),init));
292: if (ngrp!=nogrp || ngrp==0 || arsd>=oarsd) {
293: idsrr = nxtsrr - its;
294: } else {
295: idsrr = (PetscInt)PetscFloorReal(alpha+beta*(itrsdold[eps->nconv]-itrsd[eps->nconv])*PetscLogReal(arsd/eps->tol)/PetscLogReal(arsd/oarsd));
296: idsrr = PetscMax(1,idsrr);
297: }
298: nxtsrr = PetscMin(nxtsrr,its+idsrr);
300: /* Compute nxtort (iteration of next orthogonalization step) */
301: PetscCall(DSCond(eps->ds,&tcond));
302: idort = PetscMax(1,(PetscInt)PetscFloorReal(orttol/PetscMax(1,PetscLog10Real(tcond))));
303: nxtort = PetscMin(its+idort,nxtsrr);
304: PetscCall(PetscInfo(eps,"Updated iteration counts: nxtort=%" PetscInt_FMT ", nxtsrr=%" PetscInt_FMT "\n",nxtort,nxtsrr));
306: /* V(:,idx) = AV(:,idx) */
307: PetscCall(BVSetActiveColumns(eps->V,eps->nconv,nv));
308: PetscCall(BVSetActiveColumns(AV,eps->nconv,nv));
309: PetscCall(BVCopy(AV,eps->V));
310: its++;
312: /* Orthogonalization loop */
313: do {
314: PetscCall(BVGetMatrix(eps->V,&B,&indef));
315: PetscCall(BVSetMatrix(eps->V,NULL,PETSC_FALSE));
316: while (its<nxtort) {
317: /* A(:,idx) = OP*V(:,idx) with normalization */
318: PetscCall(BVMatMult(eps->V,S,AV));
319: PetscCall(BVCopy(AV,eps->V));
320: PetscCall(BVNormalize(eps->V,NULL));
321: its++;
322: }
323: PetscCall(BVSetMatrix(eps->V,B,indef));
324: /* Orthonormalize vectors */
325: PetscCall(BVOrthogonalize(eps->V,NULL));
326: nxtort = PetscMin(its+idort,nxtsrr);
327: } while (its<nxtsrr);
328: }
330: PetscCall(PetscFree6(rsd,orsd,oeigr,oeigi,itrsd,itrsdold));
331: PetscCall(BVDestroy(&AV));
332: PetscCall(BVDestroy(&R));
333: PetscCall(STRestoreOperator(eps->st,&S));
334: PetscCall(DSTruncate(eps->ds,eps->nconv,PETSC_TRUE));
335: PetscFunctionReturn(PETSC_SUCCESS);
336: }
338: static PetscErrorCode EPSDestroy_Subspace(EPS eps)
339: {
340: PetscFunctionBegin;
341: PetscCall(PetscFree(eps->data));
342: PetscFunctionReturn(PETSC_SUCCESS);
343: }
345: SLEPC_EXTERN PetscErrorCode EPSCreate_Subspace(EPS eps)
346: {
347: EPS_SUBSPACE *ctx;
349: PetscFunctionBegin;
350: PetscCall(PetscNew(&ctx));
351: eps->data = (void*)ctx;
353: eps->useds = PETSC_TRUE;
354: eps->categ = EPS_CATEGORY_OTHER;
356: eps->ops->solve = EPSSolve_Subspace;
357: eps->ops->setup = EPSSetUp_Subspace;
358: eps->ops->setupsort = EPSSetUpSort_Subspace;
359: eps->ops->destroy = EPSDestroy_Subspace;
360: eps->ops->backtransform = EPSBackTransform_Default;
361: eps->ops->computevectors = EPSComputeVectors_Schur;
362: PetscFunctionReturn(PETSC_SUCCESS);
363: }