Actual source code: subspace.c
slepc-3.14.0 2020-09-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: PetscErrorCode EPSSetUp_Subspace(EPS eps)
29: {
33: EPSCheckDefinite(eps);
34: EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
35: if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
36: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
37: if (eps->which!=EPS_LARGEST_MAGNITUDE && eps->which!=EPS_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only largest magnitude or target magnitude eigenvalues");
38: EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_EXTRACTION | EPS_FEATURE_TWOSIDED);
39: if (eps->converged != EPSConvergedRelative) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver only supports relative convergence test");
41: EPSAllocateSolution(eps,0);
42: EPS_SetInnerProduct(eps);
43: if (eps->ishermitian) {
44: DSSetType(eps->ds,DSHEP);
45: } else {
46: DSSetType(eps->ds,DSNHEP);
47: }
48: DSAllocate(eps->ds,eps->ncv);
49: return(0);
50: }
52: /*
53: EPSSubspaceFindGroup - Find a group of nearly equimodular eigenvalues, provided
54: in arrays wr and wi, according to the tolerance grptol. Also the 2-norms
55: of the residuals must be passed in (rsd). Arrays are processed from index
56: l to index m only. The output information is:
58: ngrp - number of entries of the group
59: ctr - (w(l)+w(l+ngrp-1))/2
60: ae - average of wr(l),...,wr(l+ngrp-1)
61: arsd - average of rsd(l),...,rsd(l+ngrp-1)
62: */
63: static PetscErrorCode EPSSubspaceFindGroup(PetscInt l,PetscInt m,PetscScalar *wr,PetscScalar *wi,PetscReal *rsd,PetscReal grptol,PetscInt *ngrp,PetscReal *ctr,PetscReal *ae,PetscReal *arsd)
64: {
65: PetscInt i;
66: PetscReal rmod,rmod1;
69: *ngrp = 0;
70: *ctr = 0;
71: rmod = SlepcAbsEigenvalue(wr[l],wi[l]);
73: for (i=l;i<m;) {
74: rmod1 = SlepcAbsEigenvalue(wr[i],wi[i]);
75: if (PetscAbsReal(rmod-rmod1) > grptol*(rmod+rmod1)) break;
76: *ctr = (rmod+rmod1)/2.0;
77: if (wi[i] != 0.0) {
78: (*ngrp)+=2;
79: i+=2;
80: } else {
81: (*ngrp)++;
82: i++;
83: }
84: }
86: *ae = 0;
87: *arsd = 0;
88: if (*ngrp) {
89: for (i=l;i<l+*ngrp;i++) {
90: (*ae) += PetscRealPart(wr[i]);
91: (*arsd) += rsd[i]*rsd[i];
92: }
93: *ae = *ae / *ngrp;
94: *arsd = PetscSqrtScalar(*arsd / *ngrp);
95: }
96: return(0);
97: }
99: /*
100: EPSSubspaceResidualNorms - Computes the column norms of residual vectors
101: OP*V(1:n,l:m) - V*T(1:m,l:m), where, on entry, OP*V has been computed and
102: stored in R. On exit, rsd(l) to rsd(m) contain the computed norms.
103: */
104: static PetscErrorCode EPSSubspaceResidualNorms(BV R,BV V,Mat T,PetscInt l,PetscInt m,PetscScalar *eigi,PetscReal *rsd)
105: {
107: PetscInt i;
110: BVMult(R,-1.0,1.0,V,T);
111: for (i=l;i<m;i++) { BVNormColumnBegin(R,i,NORM_2,rsd+i); }
112: for (i=l;i<m;i++) { BVNormColumnEnd(R,i,NORM_2,rsd+i); }
113: #if !defined(PETSC_USE_COMPLEX)
114: for (i=l;i<m-1;i++) {
115: if (eigi[i]!=0.0) {
116: rsd[i] = SlepcAbs(rsd[i],rsd[i+1]);
117: rsd[i+1] = rsd[i];
118: i++;
119: }
120: }
121: #endif
122: return(0);
123: }
125: PetscErrorCode EPSSolve_Subspace(EPS eps)
126: {
128: Mat H,Q,S,T,B;
129: BV AV,R;
130: PetscBool indef;
131: PetscInt i,k,ld,ngrp,nogrp,*itrsd,*itrsdold;
132: PetscInt nxtsrr,idsrr,idort,nxtort,nv,ncv = eps->ncv,its;
133: PetscReal arsd,oarsd,ctr,octr,ae,oae,*rsd,tcond=1.0;
134: /* Parameters */
135: PetscInt init = 5; /* Number of initial iterations */
136: PetscReal stpfac = 1.5; /* Max num of iter before next SRR step */
137: PetscReal alpha = 1.0; /* Used to predict convergence of next residual */
138: PetscReal beta = 1.1; /* Used to predict convergence of next residual */
139: PetscReal grptol = SLEPC_DEFAULT_TOL; /* Tolerance for EPSSubspaceFindGroup */
140: PetscReal cnvtol = 1e-6; /* Convergence criterion for cnv */
141: PetscInt orttol = 2; /* Number of decimal digits whose loss
142: can be tolerated in orthogonalization */
145: its = 0;
146: PetscMalloc3(ncv,&rsd,ncv,&itrsd,ncv,&itrsdold);
147: DSGetLeadingDimension(eps->ds,&ld);
148: BVDuplicate(eps->V,&AV);
149: BVDuplicate(eps->V,&R);
150: STGetOperator(eps->st,&S);
152: for (i=0;i<ncv;i++) {
153: rsd[i] = 0.0;
154: itrsd[i] = -1;
155: }
157: /* Complete the initial basis with random vectors and orthonormalize them */
158: for (k=eps->nini;k<ncv;k++) {
159: BVSetRandomColumn(eps->V,k);
160: BVOrthonormalizeColumn(eps->V,k,PETSC_TRUE,NULL,NULL);
161: }
163: while (eps->reason == EPS_CONVERGED_ITERATING) {
164: eps->its++;
165: nv = PetscMin(eps->nconv+eps->mpd,ncv);
166: DSSetDimensions(eps->ds,nv,0,eps->nconv,0);
168: /* Find group in previously computed eigenvalues */
169: EPSSubspaceFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,rsd,grptol,&nogrp,&octr,&oae,&oarsd);
171: /* AV(:,idx) = OP * V(:,idx) */
172: BVSetActiveColumns(eps->V,eps->nconv,nv);
173: BVSetActiveColumns(AV,eps->nconv,nv);
174: BVMatMult(eps->V,S,AV);
176: /* T(:,idx) = V' * AV(:,idx) */
177: BVSetActiveColumns(eps->V,0,nv);
178: DSGetMat(eps->ds,DS_MAT_A,&H);
179: BVDot(AV,eps->V,H);
180: DSRestoreMat(eps->ds,DS_MAT_A,&H);
181: DSSetState(eps->ds,DS_STATE_RAW);
183: /* Solve projected problem */
184: DSSolve(eps->ds,eps->eigr,eps->eigi);
185: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
186: DSSynchronize(eps->ds,eps->eigr,eps->eigi);
188: /* Update vectors V(:,idx) = V * U(:,idx) */
189: DSGetMat(eps->ds,DS_MAT_Q,&Q);
190: BVSetActiveColumns(AV,0,nv);
191: BVSetActiveColumns(R,0,nv);
192: BVMultInPlace(eps->V,Q,eps->nconv,nv);
193: BVMultInPlace(AV,Q,eps->nconv,nv);
194: BVCopy(AV,R);
195: MatDestroy(&Q);
197: /* Convergence check */
198: DSGetMat(eps->ds,DS_MAT_A,&T);
199: EPSSubspaceResidualNorms(R,eps->V,T,eps->nconv,nv,eps->eigi,rsd);
200: DSRestoreMat(eps->ds,DS_MAT_A,&T);
202: for (i=eps->nconv;i<nv;i++) {
203: itrsdold[i] = itrsd[i];
204: itrsd[i] = its;
205: eps->errest[i] = rsd[i];
206: }
208: for (;;) {
209: /* Find group in currently computed eigenvalues */
210: EPSSubspaceFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,eps->errest,grptol,&ngrp,&ctr,&ae,&arsd);
211: if (ngrp!=nogrp) break;
212: if (ngrp==0) break;
213: if (PetscAbsReal(ae-oae)>ctr*cnvtol*(itrsd[eps->nconv]-itrsdold[eps->nconv])) break;
214: if (arsd>ctr*eps->tol) break;
215: eps->nconv = eps->nconv + ngrp;
216: if (eps->nconv>=nv) break;
217: }
219: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
220: (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
221: if (eps->reason != EPS_CONVERGED_ITERATING) break;
223: /* Compute nxtsrr (iteration of next projection step) */
224: nxtsrr = PetscMin(eps->max_it,PetscMax((PetscInt)PetscFloorReal(stpfac*its),init));
226: if (ngrp!=nogrp || ngrp==0 || arsd>=oarsd) {
227: idsrr = nxtsrr - its;
228: } else {
229: idsrr = (PetscInt)PetscFloorReal(alpha+beta*(itrsdold[eps->nconv]-itrsd[eps->nconv])*PetscLogReal(arsd/eps->tol)/PetscLogReal(arsd/oarsd));
230: idsrr = PetscMax(1,idsrr);
231: }
232: nxtsrr = PetscMin(nxtsrr,its+idsrr);
234: /* Compute nxtort (iteration of next orthogonalization step) */
235: DSCond(eps->ds,&tcond);
236: idort = PetscMax(1,(PetscInt)PetscFloorReal(orttol/PetscMax(1,PetscLog10Real(tcond))));
237: nxtort = PetscMin(its+idort,nxtsrr);
239: /* V(:,idx) = AV(:,idx) */
240: BVSetActiveColumns(eps->V,eps->nconv,nv);
241: BVSetActiveColumns(AV,eps->nconv,nv);
242: BVCopy(AV,eps->V);
243: its++;
245: /* Orthogonalization loop */
246: do {
247: BVGetMatrix(eps->V,&B,&indef);
248: BVSetMatrix(eps->V,NULL,PETSC_FALSE);
249: while (its<nxtort) {
250: /* A(:,idx) = OP*V(:,idx) with normalization */
251: BVMatMult(eps->V,S,AV);
252: BVCopy(AV,eps->V);
253: BVNormalize(eps->V,NULL);
254: its++;
255: }
256: BVSetMatrix(eps->V,B,indef);
257: /* Orthonormalize vectors */
258: BVOrthogonalize(eps->V,NULL);
259: nxtort = PetscMin(its+idort,nxtsrr);
260: } while (its<nxtsrr);
261: }
263: PetscFree3(rsd,itrsd,itrsdold);
264: BVDestroy(&AV);
265: BVDestroy(&R);
266: STRestoreOperator(eps->st,&S);
267: /* truncate Schur decomposition and change the state to raw so that
268: DSVectors() computes eigenvectors from scratch */
269: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
270: DSSetState(eps->ds,DS_STATE_RAW);
271: return(0);
272: }
274: PetscErrorCode EPSDestroy_Subspace(EPS eps)
275: {
279: PetscFree(eps->data);
280: return(0);
281: }
283: SLEPC_EXTERN PetscErrorCode EPSCreate_Subspace(EPS eps)
284: {
286: eps->useds = PETSC_TRUE;
287: eps->categ = EPS_CATEGORY_OTHER;
289: eps->ops->solve = EPSSolve_Subspace;
290: eps->ops->setup = EPSSetUp_Subspace;
291: eps->ops->setupsort = EPSSetUpSort_Default;
292: eps->ops->destroy = EPSDestroy_Subspace;
293: eps->ops->backtransform = EPSBackTransform_Default;
294: eps->ops->computevectors = EPSComputeVectors_Schur;
295: return(0);
296: }