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: June 2004
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: SLEPc - Scalable Library for Eigenvalue Problem Computations
21: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
23: This file is part of SLEPc. See the README file for conditions of use
24: and additional information.
25: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26: */
28: #include src/eps/epsimpl.h
29: #include slepcblaslapack.h
33: PetscErrorCode EPSSetUp_SUBSPACE(EPS eps)
34: {
36: PetscInt N;
39: VecGetSize(eps->vec_initial,&N);
40: if (eps->ncv) {
41: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
42: }
43: else eps->ncv = PetscMin(N,PetscMax(2*eps->nev,eps->nev+15));
44: if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);
45: if (eps->which!=EPS_LARGEST_MAGNITUDE)
46: SETERRQ(1,"Wrong value of eps->which");
47: EPSAllocateSolution(eps);
48: PetscFree(eps->T);
49: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
50: EPSDefaultGetWork(eps,eps->ncv);
51: return(0);
52: }
56: /*
57: EPSHessCond - Compute the inf-norm condition number of the upper
58: Hessenberg matrix H: cond(H) = norm(H)*norm(inv(H)).
59: This routine uses Gaussian elimination with partial pivoting to
60: compute the inverse explicitly.
61: */
62: static PetscErrorCode EPSHessCond(PetscScalar* H,int n, PetscReal* cond)
63: {
64: #if defined(PETSC_MISSING_LAPACK_GETRF) || defined(SLEPC_MISSING_LAPACK_GETRI) || defined(SLEPC_MISSING_LAPACK_LANGE) || defined(SLEPC_MISSING_LAPACK_LANHS)
66: SETERRQ(PETSC_ERR_SUP,"GETRF,GETRI - Lapack routines are unavailable.");
67: #else
69: int *ipiv,lwork,info;
70: PetscScalar *work;
71: PetscReal hn,hin,*rwork;
72:
74: PetscLogEventBegin(EPS_Dense,0,0,0,0);
75: PetscMalloc(sizeof(int)*n,&ipiv);
76: lwork = n*n;
77: PetscMalloc(sizeof(PetscScalar)*lwork,&work);
78: PetscMalloc(sizeof(PetscReal)*n,&rwork);
79: hn = LAPACKlanhs_("I",&n,H,&n,rwork,1);
80: LAPACKgetrf_(&n,&n,H,&n,ipiv,&info);
81: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGETRF %d",info);
82: LAPACKgetri_(&n,H,&n,ipiv,work,&lwork,&info);
83: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGETRI %d",info);
84: hin = LAPACKlange_("I",&n,&n,H,&n,rwork,1);
85: *cond = hn * hin;
86: PetscFree(ipiv);
87: PetscFree(work);
88: PetscFree(rwork);
89: PetscLogEventEnd(EPS_Dense,0,0,0,0);
90: return(0);
91: #endif
92: }
96: /*
97: EPSFindGroup - Find a group of nearly equimodular eigenvalues, provided
98: in arrays wr and wi, according to the tolerance grptol. Also the 2-norms
99: of the residuals must be passed-in (rsd). Arrays are processed from index
100: l to index m only. The output information is:
102: ngrp - number of entries of the group
103: ctr - (w(l)+w(l+ngrp-1))/2
104: ae - average of wr(l),...,wr(l+ngrp-1)
105: arsd - average of rsd(l),...,rsd(l+ngrp-1)
106: */
107: static PetscErrorCode EPSFindGroup(int l,int m,PetscScalar *wr,PetscScalar *wi,PetscReal *rsd,
108: PetscReal grptol,int *ngrp,PetscReal *ctr,PetscReal *ae,PetscReal *arsd)
109: {
110: int i;
111: PetscReal rmod,rmod1;
114: *ngrp = 0;
115: *ctr = 0;
116:
117: rmod = SlepcAbsEigenvalue(wr[l],wi[l]);
119: for (i=l;i<m;) {
120: rmod1 = SlepcAbsEigenvalue(wr[i],wi[i]);
121: if (PetscAbsReal(rmod-rmod1) > grptol*(rmod+rmod1)) break;
122: *ctr = (rmod+rmod1)/2.0;
123: if (wi[i] != 0.0) {
124: (*ngrp)+=2;
125: i+=2;
126: } else {
127: (*ngrp)++;
128: i++;
129: }
130: }
132: *ae = 0;
133: *arsd = 0;
135: if (*ngrp) {
136: for (i=l;i<l+*ngrp;i++) {
137: (*ae) += PetscRealPart(wr[i]);
138: (*arsd) += rsd[i]*rsd[i];
139: }
140: *ae = *ae / *ngrp;
141: *arsd = PetscSqrtScalar(*arsd / *ngrp);
142: }
143: return(0);
144: }
148: /*
149: EPSSchurResidualNorms - Computes the column norms of residual vectors
150: OP*V(1:n,l:m) - V*T(1:m,l:m) were on entry, OP*V has been computed and
151: stored in AV. ldt is the leading dimension of T. On exit, rsd(l) to
152: rsd(m) contain the computed norms.
153: */
154: static PetscErrorCode EPSSchurResidualNorms(EPS eps,Vec *V,Vec *AV,PetscScalar *T,int l,int m,int ldt,PetscReal *rsd)
155: {
157: int i;
158: #if defined(PETSC_USE_COMPLEX)
159: PetscScalar t;
160: #endif
163: for (i=l;i<m;i++) {
164: VecSet(eps->work[0],0.0);
165: VecMAXPY(eps->work[0],m,T+ldt*i,V);
166: VecWAXPY(eps->work[1],-1.0,eps->work[0],AV[i]);
167: #if !defined(PETSC_USE_COMPLEX)
168: VecDot(eps->work[1],eps->work[1],rsd+i);
169: #else
170: VecDot(eps->work[1],eps->work[1],&t);
171: rsd[i] = PetscRealPart(t);
172: #endif
173: }
175: for (i=l;i<m;i++) {
176: if (i == m-1) {
177: rsd[i] = sqrt(rsd[i]);
178: } else if (T[i+1+(ldt*i)]==0.0) {
179: rsd[i] = sqrt(rsd[i]);
180: } else {
181: rsd[i] = sqrt(rsd[i]+rsd[i+1])/2.0;
182: rsd[i+1] = rsd[i];
183: i++;
184: }
185: }
186: return(0);
187: }
191: PetscErrorCode EPSSolve_SUBSPACE(EPS eps)
192: {
194: int i,ngrp,nogrp,*itrsd,*itrsdold,
195: nxtsrr,idsrr,idort,nxtort,ncv = eps->ncv,its;
196: PetscScalar *T=eps->T,*U;
197: PetscReal arsd,oarsd,ctr,octr,ae,oae,*rsd,*rsdold,norm,tcond;
198: PetscTruth breakdown;
199: /* Parameters */
200: int init = 5; /* Number of initial iterations */
201: PetscReal stpfac = 1.5, /* Max num of iter before next SRR step */
202: alpha = 1.0, /* Used to predict convergence of next residual */
203: beta = 1.1, /* Used to predict convergence of next residual */
204: grptol = 1e-8, /* Tolerance for EPSFindGroup */
205: cnvtol = 1e-6; /* Convergence criterion for cnv */
206: int orttol = 2; /* Number of decimal digits whose loss
207: can be tolerated in orthogonalization */
210: its = 0;
211: PetscMalloc(sizeof(PetscScalar)*ncv*ncv,&U);
212: PetscMalloc(sizeof(PetscReal)*ncv,&rsd);
213: PetscMalloc(sizeof(PetscReal)*ncv,&rsdold);
214: PetscMalloc(sizeof(int)*ncv,&itrsd);
215: PetscMalloc(sizeof(int)*ncv,&itrsdold);
217: /* Generate a set of random initial vectors and orthonormalize them */
218: for (i=0;i<ncv;i++) {
219: SlepcVecSetRandom(eps->V[i]);
220: rsd[i] = 0.0;
221: itrsd[i] = -1;
222: }
223: IPQRDecomposition(eps->ip,eps->V,0,ncv,PETSC_NULL,0,eps->work[0]);
224:
225: while (eps->its<eps->max_it) {
226: eps->its++;
227:
228: /* Find group in previously computed eigenvalues */
229: EPSFindGroup(eps->nconv,ncv,eps->eigr,eps->eigi,rsd,grptol,&nogrp,&octr,&oae,&oarsd);
231: /* Compute a Rayleigh-Ritz projection step
232: on the active columns (idx) */
234: /* 1. AV(:,idx) = OP * V(:,idx) */
235: for (i=eps->nconv;i<ncv;i++) {
236: STApply(eps->OP,eps->V[i],eps->AV[i]);
237: }
239: /* 2. T(:,idx) = V' * AV(:,idx) */
240: for (i=eps->nconv;i<ncv;i++) {
241: VecMDot(eps->AV[i],ncv,eps->V,T+i*ncv);
242: }
244: /* 3. Reduce projected matrix to Hessenberg form: [U,T] = hess(T) */
245: EPSDenseHessenberg(ncv,eps->nconv,T,ncv,U);
246:
247: /* 4. Reduce T to quasi-triangular (Schur) form */
248: EPSDenseSchur(ncv,eps->nconv,T,ncv,U,eps->eigr,eps->eigi);
250: /* 5. Sort diagonal elements in T and accumulate rotations on U */
251: EPSSortDenseSchur(ncv,eps->nconv,T,ncv,U,eps->eigr,eps->eigi,eps->which);
252:
253: /* 6. AV(:,idx) = AV * U(:,idx) */
254: for (i=eps->nconv;i<ncv;i++) {
255: VecSet(eps->work[i],0.0);
256: VecMAXPY(eps->work[i],ncv,U+ncv*i,eps->AV);
257: }
258: for (i=eps->nconv;i<ncv;i++) {
259: VecCopy(eps->work[i],eps->AV[i]);
260: }
261:
262: /* 7. V(:,idx) = V * U(:,idx) */
263: for (i=eps->nconv;i<ncv;i++) {
264: VecSet(eps->work[i],0.0);
265: VecMAXPY(eps->work[i],ncv,U+ncv*i,eps->V);
266: }
267: for (i=eps->nconv;i<ncv;i++) {
268: VecCopy(eps->work[i],eps->V[i]);
269: }
270:
271: /* Compute residuals */
272: for (i=0;i<ncv;i++) { rsdold[i] = rsd[i]; }
274: EPSSchurResidualNorms(eps,eps->V,eps->AV,T,eps->nconv,ncv,ncv,rsd);
276: for (i=0;i<ncv;i++) {
277: eps->errest[i] = rsd[i] / SlepcAbsEigenvalue(eps->eigr[i],eps->eigi[i]);
278: }
279: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,ncv);
280:
281: /* Convergence check */
282: for (i=0;i<ncv;i++) { itrsdold[i] = itrsd[i]; }
283: for (i=eps->nconv;i<ncv;i++) { itrsd[i] = its; }
284:
285: for (;;) {
286: /* Find group in currently computed eigenvalues */
287: EPSFindGroup(eps->nconv,ncv,eps->eigr,eps->eigi,rsd,grptol,&ngrp,&ctr,&ae,&arsd);
288: if (ngrp!=nogrp) break;
289: if (ngrp==0) break;
290: if (PetscAbsScalar(ae-oae)>ctr*cnvtol*(itrsd[eps->nconv]-itrsdold[eps->nconv])) break;
291: if (arsd>ctr*eps->tol) break;
292: eps->nconv = eps->nconv + ngrp;
293: if (eps->nconv>=ncv) break;
294: }
295:
296: if (eps->nconv>=eps->nev) break;
297:
298: /* Compute nxtsrr (iteration of next projection step) */
299: nxtsrr = PetscMin(eps->max_it,PetscMax((int)floor(stpfac*its), init));
300:
301: if (ngrp!=nogrp || ngrp==0 || arsd>=oarsd) {
302: idsrr = nxtsrr - its;
303: } else {
304: idsrr = (int)floor(alpha+beta*(itrsdold[eps->nconv]-itrsd[eps->nconv])*log(arsd/eps->tol)/log(arsd/oarsd));
305: idsrr = PetscMax(1,idsrr);
306: }
307: nxtsrr = PetscMin(nxtsrr,its+idsrr);
309: /* Compute nxtort (iteration of next orthogonalization step) */
310: PetscMemcpy(U,T,sizeof(PetscScalar)*ncv);
311: EPSHessCond(U,ncv,&tcond);
312: idort = PetscMax(1,(int)floor(orttol/PetscMax(1,log10(tcond))));
313: nxtort = PetscMin(its+idort, nxtsrr);
315: /* V(:,idx) = AV(:,idx) */
316: for (i=eps->nconv;i<ncv;i++) {
317: VecCopy(eps->AV[i],eps->V[i]);
318: }
319: its++;
321: /* Orthogonalization loop */
322: do {
323: while (its<nxtort) {
324:
325: /* AV(:,idx) = OP * V(:,idx) */
326: for (i=eps->nconv;i<ncv;i++) {
327: STApply(eps->OP,eps->V[i],eps->AV[i]);
328: }
329:
330: /* V(:,idx) = AV(:,idx) with normalization */
331: for (i=eps->nconv;i<ncv;i++) {
332: VecCopy(eps->AV[i],eps->V[i]);
333: VecNorm(eps->V[i],NORM_INFINITY,&norm);
334: VecScale(eps->V[i],1/norm);
335: }
336:
337: its++;
338: }
339: /* Orthonormalize vectors */
340: for (i=eps->nconv;i<ncv;i++) {
341: IPOrthogonalize(eps->ip,i+eps->nds,PETSC_NULL,eps->DSV,eps->V[i],PETSC_NULL,&norm,&breakdown,eps->work[0]);
342: if (breakdown) {
343: SlepcVecSetRandom(eps->V[i]);
344: IPOrthogonalize(eps->ip,i+eps->nds,PETSC_NULL,eps->DSV,eps->V[i],PETSC_NULL,&norm,&breakdown,eps->work[0]);
345: }
346: VecScale(eps->V[i],1/norm);
347: }
348: nxtort = PetscMin(its+idort,nxtsrr);
349: } while (its<nxtsrr);
350: }
352: PetscFree(U);
353: PetscFree(rsd);
354: PetscFree(rsdold);
355: PetscFree(itrsd);
356: PetscFree(itrsdold);
358: if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
359: else eps->reason = EPS_DIVERGED_ITS;
361: return(0);
362: }
367: PetscErrorCode EPSCreate_SUBSPACE(EPS eps)
368: {
370: eps->ops->solve = EPSSolve_SUBSPACE;
371: eps->ops->setup = EPSSetUp_SUBSPACE;
372: eps->ops->destroy = EPSDestroy_Default;
373: eps->ops->backtransform = EPSBackTransform_Default;
374: eps->ops->computevectors = EPSComputeVectors_Schur;
375: return(0);
376: }