Actual source code: subspace.c
slepc-3.7.2 2016-07-19
1: /*
3: SLEPc eigensolver: "subspace"
5: Method: Subspace Iteration
7: Algorithm:
9: Subspace iteration with Rayleigh-Ritz projection and locking,
10: based on the SRRIT implementation.
12: References:
14: [1] "Subspace Iteration in SLEPc", SLEPc Technical Report STR-3,
15: available at http://slepc.upv.es.
17: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18: SLEPc - Scalable Library for Eigenvalue Problem Computations
19: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
21: This file is part of SLEPc.
23: SLEPc is free software: you can redistribute it and/or modify it under the
24: terms of version 3 of the GNU Lesser General Public License as published by
25: the Free Software Foundation.
27: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
28: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
29: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
30: more details.
32: You should have received a copy of the GNU Lesser General Public License
33: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
34: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35: */
37: #include <slepc/private/epsimpl.h>
41: PetscErrorCode EPSSetUp_Subspace(EPS eps)
42: {
46: EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
47: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
48: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
49: if (eps->which!=EPS_LARGEST_MAGNITUDE && eps->which!=EPS_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
50: if (!eps->extraction) {
51: EPSSetExtraction(eps,EPS_RITZ);
52: } else if (eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
53: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
55: EPSAllocateSolution(eps,0);
56: EPS_SetInnerProduct(eps);
57: if (eps->ishermitian) {
58: DSSetType(eps->ds,DSHEP);
59: } else {
60: DSSetType(eps->ds,DSNHEP);
61: }
62: DSAllocate(eps->ds,eps->ncv);
63: EPSSetWorkVecs(eps,1);
65: if (eps->isgeneralized && eps->ishermitian && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method does not work for indefinite problems");
66: return(0);
67: }
71: /*
72: EPSSubspaceFindGroup - Find a group of nearly equimodular eigenvalues, provided
73: in arrays wr and wi, according to the tolerance grptol. Also the 2-norms
74: of the residuals must be passed in (rsd). Arrays are processed from index
75: l to index m only. The output information is:
77: ngrp - number of entries of the group
78: ctr - (w(l)+w(l+ngrp-1))/2
79: ae - average of wr(l),...,wr(l+ngrp-1)
80: arsd - average of rsd(l),...,rsd(l+ngrp-1)
81: */
82: static PetscErrorCode EPSSubspaceFindGroup(PetscInt l,PetscInt m,PetscScalar *wr,PetscScalar *wi,PetscReal *rsd,PetscReal grptol,PetscInt *ngrp,PetscReal *ctr,PetscReal *ae,PetscReal *arsd)
83: {
84: PetscInt i;
85: PetscReal rmod,rmod1;
88: *ngrp = 0;
89: *ctr = 0;
90: rmod = SlepcAbsEigenvalue(wr[l],wi[l]);
92: for (i=l;i<m;) {
93: rmod1 = SlepcAbsEigenvalue(wr[i],wi[i]);
94: if (PetscAbsReal(rmod-rmod1) > grptol*(rmod+rmod1)) break;
95: *ctr = (rmod+rmod1)/2.0;
96: if (wi[i] != 0.0) {
97: (*ngrp)+=2;
98: i+=2;
99: } else {
100: (*ngrp)++;
101: i++;
102: }
103: }
105: *ae = 0;
106: *arsd = 0;
107: if (*ngrp) {
108: for (i=l;i<l+*ngrp;i++) {
109: (*ae) += PetscRealPart(wr[i]);
110: (*arsd) += rsd[i]*rsd[i];
111: }
112: *ae = *ae / *ngrp;
113: *arsd = PetscSqrtScalar(*arsd / *ngrp);
114: }
115: return(0);
116: }
120: /*
121: EPSSubspaceResidualNorms - Computes the column norms of residual vectors
122: OP*V(1:n,l:m) - V*T(1:m,l:m), where, on entry, OP*V has been computed and
123: stored in AV. ldt is the leading dimension of T. On exit, rsd(l) to
124: rsd(m) contain the computed norms.
125: */
126: static PetscErrorCode EPSSubspaceResidualNorms(BV V,BV AV,PetscScalar *T,PetscInt l,PetscInt m,PetscInt ldt,Vec w,PetscReal *rsd)
127: {
129: PetscInt i,k;
130: PetscScalar t;
133: for (i=l;i<m;i++) {
134: if (i==m-1 || T[i+1+ldt*i]==0.0) k=i+1;
135: else k=i+2;
136: BVSetActiveColumns(V,0,k);
137: BVCopyVec(AV,i,w);
138: BVMultVec(V,-1.0,1.0,w,T+ldt*i);
139: VecDot(w,w,&t);
140: rsd[i] = PetscRealPart(t);
141: }
142: for (i=l;i<m;i++) {
143: if (i == m-1) {
144: rsd[i] = PetscSqrtReal(rsd[i]);
145: } else if (T[i+1+(ldt*i)]==0.0) {
146: rsd[i] = PetscSqrtReal(rsd[i]);
147: } else {
148: rsd[i] = PetscSqrtReal((rsd[i]+rsd[i+1])/2.0);
149: rsd[i+1] = rsd[i];
150: i++;
151: }
152: }
153: return(0);
154: }
158: PetscErrorCode EPSSolve_Subspace(EPS eps)
159: {
161: Vec v,av,w=eps->work[0];
162: Mat H,Q;
163: BV AV;
164: PetscInt i,k,ld,ngrp,nogrp,*itrsd,*itrsdold;
165: PetscInt nxtsrr,idsrr,idort,nxtort,nv,ncv = eps->ncv,its;
166: PetscScalar *T;
167: PetscReal arsd,oarsd,ctr,octr,ae,oae,*rsd,norm,tcond=1.0;
168: PetscBool breakdown;
169: /* Parameters */
170: PetscInt init = 5; /* Number of initial iterations */
171: PetscReal stpfac = 1.5; /* Max num of iter before next SRR step */
172: PetscReal alpha = 1.0; /* Used to predict convergence of next residual */
173: PetscReal beta = 1.1; /* Used to predict convergence of next residual */
174: PetscReal grptol = 1e-8; /* Tolerance for EPSSubspaceFindGroup */
175: PetscReal cnvtol = 1e-6; /* Convergence criterion for cnv */
176: PetscInt orttol = 2; /* Number of decimal digits whose loss
177: can be tolerated in orthogonalization */
180: its = 0;
181: PetscMalloc3(ncv,&rsd,ncv,&itrsd,ncv,&itrsdold);
182: DSGetLeadingDimension(eps->ds,&ld);
183: BVDuplicate(eps->V,&AV);
185: for (i=0;i<ncv;i++) {
186: rsd[i] = 0.0;
187: itrsd[i] = -1;
188: }
190: /* Complete the initial basis with random vectors and orthonormalize them */
191: k = eps->nini;
192: while (k<ncv) {
193: BVSetRandomColumn(eps->V,k);
194: BVOrthogonalizeColumn(eps->V,k,NULL,&norm,&breakdown);
195: if (norm>0.0 && !breakdown) {
196: BVScaleColumn(eps->V,k,1.0/norm);
197: k++;
198: }
199: }
201: while (eps->reason == EPS_CONVERGED_ITERATING) {
202: eps->its++;
203: nv = PetscMin(eps->nconv+eps->mpd,ncv);
204: DSSetDimensions(eps->ds,nv,0,eps->nconv,0);
206: /* Find group in previously computed eigenvalues */
207: EPSSubspaceFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,rsd,grptol,&nogrp,&octr,&oae,&oarsd);
209: /* AV(:,idx) = OP * V(:,idx) */
210: for (i=eps->nconv;i<nv;i++) {
211: BVGetColumn(eps->V,i,&v);
212: BVGetColumn(AV,i,&av);
213: STApply(eps->st,v,av);
214: BVRestoreColumn(eps->V,i,&v);
215: BVRestoreColumn(AV,i,&av);
216: }
218: /* T(:,idx) = V' * AV(:,idx) */
219: BVSetActiveColumns(eps->V,0,nv);
220: BVSetActiveColumns(AV,eps->nconv,nv);
221: DSGetMat(eps->ds,DS_MAT_A,&H);
222: BVDot(AV,eps->V,H);
223: DSRestoreMat(eps->ds,DS_MAT_A,&H);
224: DSSetState(eps->ds,DS_STATE_RAW);
226: /* Solve projected problem */
227: DSSolve(eps->ds,eps->eigr,eps->eigi);
228: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
230: /* Update vectors V(:,idx) = V * U(:,idx) */
231: DSGetMat(eps->ds,DS_MAT_Q,&Q);
232: BVSetActiveColumns(AV,0,nv);
233: BVMultInPlace(eps->V,Q,eps->nconv,nv);
234: BVMultInPlace(AV,Q,eps->nconv,nv);
235: MatDestroy(&Q);
237: /* Convergence check */
238: DSGetArray(eps->ds,DS_MAT_A,&T);
239: EPSSubspaceResidualNorms(eps->V,AV,T,eps->nconv,nv,ld,w,rsd);
240: DSRestoreArray(eps->ds,DS_MAT_A,&T);
242: for (i=eps->nconv;i<nv;i++) {
243: itrsdold[i] = itrsd[i];
244: itrsd[i] = its;
245: eps->errest[i] = rsd[i];
246: }
248: for (;;) {
249: /* Find group in currently computed eigenvalues */
250: EPSSubspaceFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,eps->errest,grptol,&ngrp,&ctr,&ae,&arsd);
251: if (ngrp!=nogrp) break;
252: if (ngrp==0) break;
253: if (PetscAbsReal(ae-oae)>ctr*cnvtol*(itrsd[eps->nconv]-itrsdold[eps->nconv])) break;
254: if (arsd>ctr*eps->tol) break;
255: eps->nconv = eps->nconv + ngrp;
256: if (eps->nconv>=nv) break;
257: }
259: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
260: (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
261: if (eps->reason != EPS_CONVERGED_ITERATING) break;
263: /* Compute nxtsrr (iteration of next projection step) */
264: nxtsrr = PetscMin(eps->max_it,PetscMax((PetscInt)PetscFloorReal(stpfac*its),init));
266: if (ngrp!=nogrp || ngrp==0 || arsd>=oarsd) {
267: idsrr = nxtsrr - its;
268: } else {
269: idsrr = (PetscInt)PetscFloorReal(alpha+beta*(itrsdold[eps->nconv]-itrsd[eps->nconv])*PetscLogReal(arsd/eps->tol)/PetscLogReal(arsd/oarsd));
270: idsrr = PetscMax(1,idsrr);
271: }
272: nxtsrr = PetscMin(nxtsrr,its+idsrr);
274: /* Compute nxtort (iteration of next orthogonalization step) */
275: DSCond(eps->ds,&tcond);
276: idort = PetscMax(1,(PetscInt)PetscFloorReal(orttol/PetscMax(1,PetscLog10Real(tcond))));
277: nxtort = PetscMin(its+idort,nxtsrr);
279: /* V(:,idx) = AV(:,idx) */
280: BVSetActiveColumns(eps->V,eps->nconv,nv);
281: BVSetActiveColumns(AV,eps->nconv,nv);
282: BVCopy(AV,eps->V);
283: its++;
285: /* Orthogonalization loop */
286: do {
287: while (its<nxtort) {
289: /* A(:,idx) = OP*V(:,idx) with normalization */
290: for (i=eps->nconv;i<nv;i++) {
291: BVGetColumn(eps->V,i,&v);
292: STApply(eps->st,v,w);
293: VecCopy(w,v);
294: BVRestoreColumn(eps->V,i,&v);
295: BVNormColumn(eps->V,i,NORM_INFINITY,&norm);
296: BVScaleColumn(eps->V,i,1/norm);
297: }
298: its++;
299: }
300: /* Orthonormalize vectors */
301: for (i=eps->nconv;i<nv;i++) {
302: BVOrthogonalizeColumn(eps->V,i,NULL,&norm,&breakdown);
303: if (breakdown) {
304: BVSetRandomColumn(eps->V,i);
305: BVOrthogonalizeColumn(eps->V,i,NULL,&norm,&breakdown);
306: }
307: BVScaleColumn(eps->V,i,1/norm);
308: }
309: nxtort = PetscMin(its+idort,nxtsrr);
310: } while (its<nxtsrr);
311: }
313: PetscFree3(rsd,itrsd,itrsdold);
314: BVDestroy(&AV);
315: /* truncate Schur decomposition and change the state to raw so that
316: DSVectors() computes eigenvectors from scratch */
317: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
318: DSSetState(eps->ds,DS_STATE_RAW);
319: return(0);
320: }
324: PetscErrorCode EPSDestroy_Subspace(EPS eps)
325: {
329: PetscFree(eps->data);
330: return(0);
331: }
335: PETSC_EXTERN PetscErrorCode EPSCreate_Subspace(EPS eps)
336: {
338: eps->ops->setup = EPSSetUp_Subspace;
339: eps->ops->solve = EPSSolve_Subspace;
340: eps->ops->destroy = EPSDestroy_Subspace;
341: eps->ops->backtransform = EPSBackTransform_Default;
342: eps->ops->computevectors = EPSComputeVectors_Schur;
343: return(0);
344: }