Actual source code: subspace.c
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://www.grycap.upv.es/slepc.
17: Last update: Feb 2009
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: SLEPc - Scalable Library for Eigenvalue Problem Computations
21: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
23: This file is part of SLEPc.
24:
25: SLEPc is free software: you can redistribute it and/or modify it under the
26: terms of version 3 of the GNU Lesser General Public License as published by
27: the Free Software Foundation.
29: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
30: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
31: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
32: more details.
34: You should have received a copy of the GNU Lesser General Public License
35: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
36: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37: */
39: #include <private/epsimpl.h> /*I "slepceps.h" I*/
40: #include <slepcblaslapack.h>
42: PetscErrorCode EPSSolve_Subspace(EPS);
44: typedef struct {
45: Vec *AV;
46: } EPS_SUBSPACE;
50: PetscErrorCode EPSSetUp_Subspace(EPS eps)
51: {
53: EPS_SUBSPACE *ctx = (EPS_SUBSPACE *)eps->data;
56: if (eps->ncv) { /* ncv set */
57: if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
58: }
59: else if (eps->mpd) { /* mpd set */
60: eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
61: }
62: else { /* neither set: defaults depend on nev being small or large */
63: if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
64: else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
65: }
66: if (!eps->mpd) eps->mpd = eps->ncv;
67: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
68: if (!eps->which) { EPSDefaultSetWhich(eps); }
69: if (eps->which!=EPS_LARGEST_MAGNITUDE && eps->which!=EPS_TARGET_MAGNITUDE)
70: SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
71: if (!eps->extraction) {
72: EPSSetExtraction(eps,EPS_RITZ);
73: } else if (eps->extraction!=EPS_RITZ) {
74: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type\n");
75: }
77: EPSAllocateSolution(eps);
78: VecDuplicateVecs(eps->t,eps->ncv,&ctx->AV);
79: PetscFree(eps->T);
80: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
81: EPSDefaultGetWork(eps,1);
83: /* dispatch solve method */
84: if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
85: eps->ops->solve = EPSSolve_Subspace;
86: return(0);
87: }
91: /*
92: EPSHessCond - Compute the inf-norm condition number of the upper
93: Hessenberg matrix H: cond(H) = norm(H)*norm(inv(H)).
94: This routine uses Gaussian elimination with partial pivoting to
95: compute the inverse explicitly.
96: */
97: static PetscErrorCode EPSHessCond(PetscInt n_,PetscScalar* H,PetscInt ldh_,PetscReal* cond)
98: {
99: #if defined(PETSC_MISSING_LAPACK_GETRF) || defined(SLEPC_MISSING_LAPACK_GETRI) || defined(SLEPC_MISSING_LAPACK_LANGE) || defined(SLEPC_MISSING_LAPACK_LANHS)
101: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF,GETRI - Lapack routines are unavailable.");
102: #else
104: PetscBLASInt *ipiv,lwork,info,n=n_,ldh=ldh_;
105: PetscScalar *work;
106: PetscReal hn,hin,*rwork;
107:
109: PetscLogEventBegin(EPS_Dense,0,0,0,0);
110: PetscMalloc(sizeof(PetscBLASInt)*n,&ipiv);
111: lwork = n*n;
112: PetscMalloc(sizeof(PetscScalar)*lwork,&work);
113: PetscMalloc(sizeof(PetscReal)*n,&rwork);
114: hn = LAPACKlanhs_("I",&n,H,&ldh,rwork);
115: LAPACKgetrf_(&n,&n,H,&ldh,ipiv,&info);
116: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGETRF %d",info);
117: LAPACKgetri_(&n,H,&ldh,ipiv,work,&lwork,&info);
118: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGETRI %d",info);
119: hin = LAPACKlange_("I",&n,&n,H,&ldh,rwork);
120: *cond = hn * hin;
121: PetscFree(ipiv);
122: PetscFree(work);
123: PetscFree(rwork);
124: PetscLogEventEnd(EPS_Dense,0,0,0,0);
125: return(0);
126: #endif
127: }
131: /*
132: EPSFindGroup - Find a group of nearly equimodular eigenvalues, provided
133: in arrays wr and wi, according to the tolerance grptol. Also the 2-norms
134: of the residuals must be passed-in (rsd). Arrays are processed from index
135: l to index m only. The output information is:
137: ngrp - number of entries of the group
138: ctr - (w(l)+w(l+ngrp-1))/2
139: ae - average of wr(l),...,wr(l+ngrp-1)
140: arsd - average of rsd(l),...,rsd(l+ngrp-1)
141: */
142: static PetscErrorCode EPSFindGroup(PetscInt l,PetscInt m,PetscScalar *wr,PetscScalar *wi,PetscReal *rsd,
143: PetscReal grptol,PetscInt *ngrp,PetscReal *ctr,PetscReal *ae,PetscReal *arsd)
144: {
145: PetscInt i;
146: PetscReal rmod,rmod1;
149: *ngrp = 0;
150: *ctr = 0;
151:
152: rmod = SlepcAbsEigenvalue(wr[l],wi[l]);
154: for (i=l;i<m;) {
155: rmod1 = SlepcAbsEigenvalue(wr[i],wi[i]);
156: if (PetscAbsReal(rmod-rmod1) > grptol*(rmod+rmod1)) break;
157: *ctr = (rmod+rmod1)/2.0;
158: if (wi[i] != 0.0) {
159: (*ngrp)+=2;
160: i+=2;
161: } else {
162: (*ngrp)++;
163: i++;
164: }
165: }
167: *ae = 0;
168: *arsd = 0;
170: if (*ngrp) {
171: for (i=l;i<l+*ngrp;i++) {
172: (*ae) += PetscRealPart(wr[i]);
173: (*arsd) += rsd[i]*rsd[i];
174: }
175: *ae = *ae / *ngrp;
176: *arsd = PetscSqrtScalar(*arsd / *ngrp);
177: }
178: return(0);
179: }
183: /*
184: EPSSchurResidualNorms - Computes the column norms of residual vectors
185: OP*V(1:n,l:m) - V*T(1:m,l:m), where, on entry, OP*V has been computed and
186: stored in AV. ldt is the leading dimension of T. On exit, rsd(l) to
187: rsd(m) contain the computed norms.
188: */
189: static PetscErrorCode EPSSchurResidualNorms(EPS eps,Vec *V,Vec *AV,PetscScalar *T,PetscInt l,PetscInt m,PetscInt ldt,PetscReal *rsd)
190: {
192: PetscInt i,k;
193: #if defined(PETSC_USE_COMPLEX)
194: PetscScalar t;
195: #endif
198: for (i=l;i<m;i++) {
199: if (i==m-1 || T[i+1+ldt*i]==0.0) k=i+1; else k=i+2;
200: VecCopy(AV[i],eps->work[0]);
201: SlepcVecMAXPBY(eps->work[0],1.0,-1.0,k,T+ldt*i,V);
202: #if !defined(PETSC_USE_COMPLEX)
203: VecDot(eps->work[0],eps->work[0],rsd+i);
204: #else
205: VecDot(eps->work[0],eps->work[0],&t);
206: rsd[i] = PetscRealPart(t);
207: #endif
208: }
210: for (i=l;i<m;i++) {
211: if (i == m-1) {
212: rsd[i] = PetscSqrtReal(rsd[i]);
213: } else if (T[i+1+(ldt*i)]==0.0) {
214: rsd[i] = PetscSqrtReal(rsd[i]);
215: } else {
216: rsd[i] = PetscSqrtReal((rsd[i]+rsd[i+1])/2.0);
217: rsd[i+1] = rsd[i];
218: i++;
219: }
220: }
221: return(0);
222: }
226: PetscErrorCode EPSSolve_Subspace(EPS eps)
227: {
229: EPS_SUBSPACE *ctx = (EPS_SUBSPACE *)eps->data;
230: PetscInt i,k,ngrp,nogrp,*itrsd,*itrsdold,
231: nxtsrr,idsrr,idort,nxtort,nv,ncv = eps->ncv,its;
232: PetscScalar *T=eps->T,*U;
233: PetscReal arsd,oarsd,ctr,octr,ae,oae,*rsd,norm,tcond=1.0;
234: PetscBool breakdown;
235: /* Parameters */
236: PetscInt init = 5; /* Number of initial iterations */
237: PetscReal stpfac = 1.5, /* Max num of iter before next SRR step */
238: alpha = 1.0, /* Used to predict convergence of next residual */
239: beta = 1.1, /* Used to predict convergence of next residual */
240: grptol = 1e-8, /* Tolerance for EPSFindGroup */
241: cnvtol = 1e-6; /* Convergence criterion for cnv */
242: PetscInt orttol = 2; /* Number of decimal digits whose loss
243: can be tolerated in orthogonalization */
246: its = 0;
247: PetscMalloc(sizeof(PetscScalar)*ncv*ncv,&U);
248: PetscMalloc(sizeof(PetscReal)*ncv,&rsd);
249: PetscMalloc(sizeof(PetscInt)*ncv,&itrsd);
250: PetscMalloc(sizeof(PetscInt)*ncv,&itrsdold);
252: for (i=0;i<ncv;i++) {
253: rsd[i] = 0.0;
254: itrsd[i] = -1;
255: }
256:
257: /* Complete the initial basis with random vectors and orthonormalize them */
258: k = eps->nini;
259: while (k<ncv) {
260: SlepcVecSetRandom(eps->V[k],eps->rand);
261: IPOrthogonalize(eps->ip,eps->nds,eps->DS,k,PETSC_NULL,eps->V,eps->V[k],PETSC_NULL,&norm,&breakdown);
262: if (norm>0.0 && !breakdown) {
263: VecScale(eps->V[k],1.0/norm);
264: k++;
265: }
266: }
268: while (eps->its<eps->max_it) {
269: eps->its++;
270: nv = PetscMin(eps->nconv+eps->mpd,ncv);
271:
272: /* Find group in previously computed eigenvalues */
273: EPSFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,rsd,grptol,&nogrp,&octr,&oae,&oarsd);
275: /* Compute a Rayleigh-Ritz projection step
276: on the active columns (idx) */
278: /* 1. AV(:,idx) = OP * V(:,idx) */
279: for (i=eps->nconv;i<nv;i++) {
280: STApply(eps->OP,eps->V[i],ctx->AV[i]);
281: }
283: /* 2. T(:,idx) = V' * AV(:,idx) */
284: for (i=eps->nconv;i<nv;i++) {
285: VecMDot(ctx->AV[i],nv,eps->V,T+i*ncv);
286: }
288: /* 3. Reduce projected matrix to Hessenberg form: [U,T] = hess(T) */
289: EPSDenseHessenberg(nv,eps->nconv,T,ncv,U);
290:
291: /* 4. Reduce T to quasi-triangular (Schur) form */
292: EPSDenseSchur(nv,eps->nconv,T,ncv,U,eps->eigr,eps->eigi);
294: /* 5. Sort diagonal elements in T and accumulate rotations on U */
295: EPSSortDenseSchur(eps,nv,eps->nconv,T,ncv,U,eps->eigr,eps->eigi);
296:
297: /* 6. AV(:,idx) = AV * U(:,idx) */
298: SlepcUpdateVectors(nv,ctx->AV,eps->nconv,nv,U,nv,PETSC_FALSE);
299:
300: /* 7. V(:,idx) = V * U(:,idx) */
301: SlepcUpdateVectors(nv,eps->V,eps->nconv,nv,U,nv,PETSC_FALSE);
302:
303: /* Convergence check */
304: EPSSchurResidualNorms(eps,eps->V,ctx->AV,T,eps->nconv,nv,ncv,rsd);
306: for (i=eps->nconv;i<nv;i++) {
307: itrsdold[i] = itrsd[i];
308: itrsd[i] = its;
309: eps->errest[i] = rsd[i];
310: }
311:
312: for (;;) {
313: /* Find group in currently computed eigenvalues */
314: EPSFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,eps->errest,grptol,&ngrp,&ctr,&ae,&arsd);
315: if (ngrp!=nogrp) break;
316: if (ngrp==0) break;
317: if (PetscAbsScalar(ae-oae)>ctr*cnvtol*(itrsd[eps->nconv]-itrsdold[eps->nconv])) break;
318: if (arsd>ctr*eps->tol) break;
319: eps->nconv = eps->nconv + ngrp;
320: if (eps->nconv>=nv) break;
321: }
322:
323: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
324: if (eps->nconv>=eps->nev) break;
325:
326: /* Compute nxtsrr (iteration of next projection step) */
327: nxtsrr = PetscMin(eps->max_it,PetscMax((PetscInt)floor(stpfac*its),init));
328:
329: if (ngrp!=nogrp || ngrp==0 || arsd>=oarsd) {
330: idsrr = nxtsrr - its;
331: } else {
332: idsrr = (PetscInt)floor(alpha+beta*(itrsdold[eps->nconv]-itrsd[eps->nconv])*log(arsd/eps->tol)/log(arsd/oarsd));
333: idsrr = PetscMax(1,idsrr);
334: }
335: nxtsrr = PetscMin(nxtsrr,its+idsrr);
337: /* Compute nxtort (iteration of next orthogonalization step) */
338: PetscMemcpy(U,T,sizeof(PetscScalar)*ncv*ncv);
339: EPSHessCond(nv,U,ncv,&tcond);
340: idort = PetscMax(1,(PetscInt)floor(orttol/PetscMax(1,log10(tcond))));
341: nxtort = PetscMin(its+idort,nxtsrr);
343: /* V(:,idx) = AV(:,idx) */
344: for (i=eps->nconv;i<nv;i++) {
345: VecCopy(ctx->AV[i],eps->V[i]);
346: }
347: its++;
349: /* Orthogonalization loop */
350: do {
351: while (its<nxtort) {
352:
353: /* AV(:,idx) = OP * V(:,idx) */
354: for (i=eps->nconv;i<nv;i++) {
355: STApply(eps->OP,eps->V[i],ctx->AV[i]);
356: }
357:
358: /* V(:,idx) = AV(:,idx) with normalization */
359: for (i=eps->nconv;i<nv;i++) {
360: VecCopy(ctx->AV[i],eps->V[i]);
361: VecNorm(eps->V[i],NORM_INFINITY,&norm);
362: VecScale(eps->V[i],1/norm);
363: }
364:
365: its++;
366: }
367: /* Orthonormalize vectors */
368: for (i=eps->nconv;i<nv;i++) {
369: IPOrthogonalize(eps->ip,eps->nds,eps->DS,i,PETSC_NULL,eps->V,eps->V[i],PETSC_NULL,&norm,&breakdown);
370: if (breakdown) {
371: SlepcVecSetRandom(eps->V[i],eps->rand);
372: IPOrthogonalize(eps->ip,eps->nds,eps->DS,i,PETSC_NULL,eps->V,eps->V[i],PETSC_NULL,&norm,&breakdown);
373: }
374: VecScale(eps->V[i],1/norm);
375: }
376: nxtort = PetscMin(its+idort,nxtsrr);
377: } while (its<nxtsrr);
378: }
380: PetscFree(U);
381: PetscFree(rsd);
382: PetscFree(itrsd);
383: PetscFree(itrsdold);
385: if (eps->nconv == eps->nev) eps->reason = EPS_CONVERGED_TOL;
386: else eps->reason = EPS_DIVERGED_ITS;
387: return(0);
388: }
392: PetscErrorCode EPSReset_Subspace(EPS eps)
393: {
395: EPS_SUBSPACE *ctx = (EPS_SUBSPACE *)eps->data;
398: VecDestroyVecs(eps->ncv,&ctx->AV);
399: PetscFree(eps->T);
400: EPSReset_Default(eps);
401: return(0);
402: }
406: PetscErrorCode EPSDestroy_Subspace(EPS eps)
407: {
411: PetscFree(eps->data);
412: return(0);
413: }
415: EXTERN_C_BEGIN
418: PetscErrorCode EPSCreate_Subspace(EPS eps)
419: {
423: PetscNewLog(eps,EPS_SUBSPACE,&eps->data);
424: eps->ops->setup = EPSSetUp_Subspace;
425: eps->ops->destroy = EPSDestroy_Subspace;
426: eps->ops->reset = EPSReset_Subspace;
427: eps->ops->backtransform = EPSBackTransform_Default;
428: eps->ops->computevectors = EPSComputeVectors_Schur;
429: return(0);
430: }
431: EXTERN_C_END