Actual source code: subspace.c

  1: /*                       

  3:    SLEPc eigensolver: "subspace"

  5:    Method: Subspace Iteration

  7:    Description:

  9:        This solver implements a version of the subspace iteration (or
 10:        simultaneous iteration) for computing an orthogonal basis of an
 11:        invariant subspace associated to the dominant eigenpairs. 

 13:    Algorithm:

 15:        The implemented algorithm is based on the SRRIT implementation (see
 16:        reference below).

 18:        The basic subspace iteration is a generalization of the power
 19:        method to m vectors, enforcing orthogonality between them to avoid
 20:        linear dependence. In addition, this implementation performs a
 21:        Rayleigh-Ritz projection procedure in order to improve convergence.
 22:        Deflation is handled by locking converged eigenvectors. For better
 23:        performance, orthogonalization and projection are performed only
 24:        when necessary.

 26:    References:

 28:        [1] G.W. Stewart and Z. Bai, "Algorithm 776. SRRIT - A Fortran 
 29:        Subroutine to Calculate the Dominant Invariant Subspace of a 
 30:        Nonsymmetric Matrix", ACM Transactions on Mathematical Software, 
 31:        23(4), pp. 494-513 (1997).

 33:    Last update: June 2004

 35: */
 36:  #include src/eps/epsimpl.h
 37:  #include slepcblaslapack.h

 41: PetscErrorCode EPSSetUp_SUBSPACE(EPS eps)
 42: {
 43:   PetscErrorCode ierr;
 44:   int            N;

 47:   VecGetSize(eps->vec_initial,&N);
 48:   if (eps->ncv) {
 49:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 50:   }
 51:   else eps->ncv = PetscMax(2*eps->nev,eps->nev+15);
 52:   eps->ncv = PetscMin(eps->ncv,N);
 53:   if (!eps->max_it) eps->max_it = PetscMax(100,N);
 54:   if (!eps->tol) eps->tol = 1.e-7;
 55:   if (eps->which!=EPS_LARGEST_MAGNITUDE)
 56:     SETERRQ(1,"Wrong value of eps->which");
 57:   EPSAllocateSolution(eps);
 58:   if (eps->T) { PetscFree(eps->T); }
 59:   PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
 60:   EPSDefaultGetWork(eps,eps->ncv);
 61:   return(0);
 62: }

 66: /*
 67:    EPSHessCond - Compute the inf-norm condition number of the upper 
 68:    Hessenberg matrix H: cond(H) = norm(H)*norm(inv(H)).
 69:    This routine uses Gaussian elimination with partial pivoting to 
 70:    compute the inverse explicitly. 
 71: */
 72: static PetscErrorCode EPSHessCond(PetscScalar* H,int n, PetscReal* cond)
 73: {
 74:   PetscErrorCode ierr;
 75:   int            *ipiv,lwork,info;
 76:   PetscScalar    *work;
 77:   PetscReal      hn,hin,*rwork;
 78: 
 80: #if defined(PETSC_MISSING_LAPACK_GETRF)
 81:   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
 82: #endif 
 83:   PetscMalloc(sizeof(int)*n,&ipiv);
 84:   lwork = n*n;
 85:   PetscMalloc(sizeof(PetscScalar)*lwork,&work);
 86:   PetscMalloc(sizeof(PetscReal)*n,&rwork);
 87:   hn = LAlanhs_("I",&n,H,&n,rwork,1);
 88:   LAgetrf_(&n,&n,H,&n,ipiv,&info);
 89:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGETRF %d",info);
 90:   LAgetri_(&n,H,&n,ipiv,work,&lwork,&info);
 91:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGETRI %d",info);
 92:   hin = LAlange_("I",&n,&n,H,&n,rwork,1);
 93:   *cond = hn * hin;
 94:   PetscFree(ipiv);
 95:   PetscFree(work);
 96:   PetscFree(rwork);
 97:   return(0);
 98: }

102: /*
103:    EPSFindGroup - Find a group of nearly equimodular eigenvalues, provided 
104:    in arrays wr and wi, according to the tolerance grptol. Also the 2-norms
105:    of the residuals must be passed-in (rsd). Arrays are processed from index 
106:    l to index m only. The output information is:

108:    ngrp - number of entries of the group
109:    ctr  - (w(l)+w(L+ngrp-1))/2
110:    ae   - average of wr(l),...,wr(l+ngrp-1)
111:    arsd - average of rsd(l),...,rsd(l+ngrp-1)
112: */
113: static PetscErrorCode EPSFindGroup(int l,int m,PetscScalar *wr,PetscScalar *wi,PetscReal *rsd,
114:   PetscReal grptol,int *ngrp,PetscReal *ctr,PetscReal *ae,PetscReal *arsd)
115: {
116:   int       i;
117:   PetscReal rmod,rmod1;

120:   *ngrp = 0;
121:   *ctr = 0;
122: 
123:   rmod = SlepcAbsEigenvalue(wr[l],wi[l]);

125:   for (i=l;i<m;) {
126:     rmod1 = SlepcAbsEigenvalue(wr[i],wi[i]);
127:     if (PetscAbsReal(rmod-rmod1) > grptol*(rmod+rmod1)) break;
128:     *ctr = (rmod+rmod1)/2.0;
129:     if (wi[i] != 0.0) {
130:       (*ngrp)+=2;
131:       i+=2;
132:     } else {
133:       (*ngrp)++;
134:       i++;
135:     }
136:   }

138:   *ae = 0;
139:   *arsd = 0;

141:   if (*ngrp) {
142:     for (i=l;i<l+*ngrp;i++) {
143:       (*ae) += PetscRealPart(wr[i]);
144:       (*arsd) += rsd[i]*rsd[i];
145:     }
146:     *ae = *ae / *ngrp;
147:     *arsd = PetscSqrtScalar(*arsd / *ngrp);
148:   }
149:   return(0);
150: }

154: /*
155:    EPSSchurResidualNorms - Computes the column norms of residual vectors
156:    OP*V(1:n,l:m) - V*T(1:m,l:m) were on entry, OP*V has been computed and 
157:    stored in AV. ldt is the leading dimension of T. On exit, rsd(l) to
158:    rsd(m) contain the computed norms.
159: */
160: static PetscErrorCode EPSSchurResidualNorms(EPS eps,Vec *V,Vec *AV,PetscScalar *T,int l,int m,int ldt,PetscReal *rsd)
161: {
162:   PetscErrorCode ierr;
163:   int            i;
164:   PetscScalar    zero = 0.0,minus = -1.0;
165: #if defined(PETSC_USE_COMPLEX)
166:   PetscScalar    t;
167: #endif

170:   for (i=l;i<m;i++) {
171:     VecSet(&zero,eps->work[0]);
172:     VecMAXPY(m,T+ldt*i,eps->work[0],V);
173:     VecWAXPY(&minus,eps->work[0],AV[i],eps->work[1]);
174: #if !defined(PETSC_USE_COMPLEX)
175:     VecDot(eps->work[1],eps->work[1],rsd+i);
176: #else
177:     VecDot(eps->work[1],eps->work[1],&t);
178:     rsd[i] = PetscRealPart(t);
179: #endif    
180:   }

182:   for (i=l;i<m;i++) {
183:     if (i == m-1) {
184:       rsd[i] = sqrt(rsd[i]);
185:     } else if (T[i+1+(ldt*i)]==0.0) {
186:       rsd[i] = sqrt(rsd[i]);
187:     } else {
188:       rsd[i] = sqrt(rsd[i]+rsd[i+1])/2.0;
189:       rsd[i+1] = rsd[i];
190:       i++;
191:     }
192:   }
193:   return(0);
194: }

198: PetscErrorCode EPSSolve_SUBSPACE(EPS eps)
199: {
200:   PetscErrorCode ierr;
201:   int            i,j,ilo,lwork,info,ngrp,nogrp,*itrsd,*itrsdold,
202:                  nxtsrr,idsrr,*iwork,idort,nxtort,ncv = eps->ncv;
203:   PetscScalar    *T=eps->T,*U,*tau,*work,t;
204:   PetscReal      arsd,oarsd,ctr,octr,ae,oae,*rsdold,norm,tcond;
205:   PetscTruth     breakdown;
206:   /* Parameters */
207:   int            init = 5;        /* Number of initial iterations */
208:   PetscReal      stpfac = 1.5,    /* Max num of iter before next SRR step */
209:                  alpha = 1.0,     /* Used to predict when  */
210:                  beta = 1.1,      /* Used to predict when  */
211:                  grptol = 1e-8,   /* Tolerance for EPSFindGroup */
212:                  cnvtol = 1e-6;   /* Convergence criterion for cnv */
213:   int            orttol = 2;      /* Number of decimal digits whose loss
214:                                      can be tolerated in orthogonalization */

217: #if defined(PETSC_BLASLAPACK_ESSL_ONLY)
218:   SETERRQ(PETSC_ERR_SUP,"GEHRD,ORGHR - Lapack routines are unavailable.");
219: #endif 
220:   eps->its = 0;
221:   eps->nconv = 0;
222:   PetscMalloc(sizeof(PetscScalar)*ncv*ncv,&U);
223:   PetscMalloc(sizeof(PetscReal)*ncv,&rsdold);
224:   PetscMalloc(sizeof(PetscScalar)*ncv,&tau);
225:   lwork = ncv*ncv;
226:   PetscMalloc(sizeof(PetscScalar)*lwork,&work);
227:   PetscMalloc(sizeof(int)*ncv,&itrsd);
228:   PetscMalloc(sizeof(int)*ncv,&itrsdold);
229:   PetscMalloc(sizeof(int)*ncv,&iwork);

231:   /* Generate a set of random initial vectors and orthonormalize them */
232:   for (i=0;i<ncv;i++) {
233:     SlepcVecSetRandom(eps->V[i]);
234:     eps->eigr[i] = 0.0;
235:     eps->eigi[i] = 0.0;
236:     eps->errest[i] = 0.0;
237:     itrsd[i] = -1;
238:   }
239:   EPSQRDecomposition(eps,eps->V,0,ncv,PETSC_NULL,0);
240: 
241:   while (eps->its<eps->max_it) {

243:     /* Find group in previously computed eigenvalues */
244:     EPSFindGroup(eps->nconv,ncv,eps->eigr,eps->eigi,eps->errest,grptol,&nogrp,&octr,&oae,&oarsd);

246:     /* Compute a Rayleigh-Ritz projection step 
247:        on the active columns (idx) */

249:     /* 1. AV(:,idx) = OP * V(:,idx) */
250:     for (i=eps->nconv;i<ncv;i++) {
251:       STApply(eps->OP,eps->V[i],eps->AV[i]);
252:     }

254:     /* 2. T(:,idx) = V' * AV(:,idx) */
255:     for (i=eps->nconv;i<ncv;i++) {
256:       VecMDot(ncv,eps->AV[i],eps->V,T+i*ncv);
257:     }

259:     /* 3. Reduce projected matrix to Hessenberg form: [U,T] = hess(T) */
260:     ilo = eps->nconv + 1;
261:     LAgehrd_(&ncv,&ilo,&ncv,T,&ncv,tau,work,&lwork,&info);
262:     if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGEHRD %d",info);
263:     for (j=0;j<ncv-1;j++) {
264:       for (i=j+2;i<ncv;i++) {
265:         U[i+j*ncv] = T[i+j*ncv];
266:         T[i+j*ncv] = 0.0;
267:       }
268:     }
269:     LAorghr_(&ncv,&ilo,&ncv,U,&ncv,tau,work,&lwork,&info);
270:     if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORGHR %d",info);
271: 
272:     /* 4. Reduce T to quasi-triangular (Schur) form */
273:     EPSDenseSchur(ncv,eps->nconv,T,U,eps->eigr,eps->eigi);

275:     /* 5. Sort diagonal elements in T and accumulate rotations on U */
276:     EPSSortDenseSchur(ncv,eps->nconv,T,U,eps->eigr,eps->eigi);
277: 
278:     /* 6. AV(:,idx) = AV * U(:,idx) */
279:     EPSReverseProjection(eps,eps->AV,U,eps->nconv,ncv,eps->work);
280: 
281:     /* 7. V(:,idx) = V * U(:,idx) */
282:     EPSReverseProjection(eps,eps->V,U,eps->nconv,ncv,eps->work);
283: 
284:     /* Compute residuals */
285:     for (i=0;i<ncv;i++) { rsdold[i] = eps->errest[i]; }

287:     EPSSchurResidualNorms(eps,eps->V,eps->AV,T,eps->nconv,ncv,ncv,eps->errest);

289:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,ncv);
290: 
291:     /* Convergence check */
292:     for (i=0;i<ncv;i++) { itrsdold[i] = itrsd[i]; }
293:     for (i=eps->nconv;i<ncv;i++) { itrsd[i] = eps->its; }
294: 
295:     for (;;) {
296:       /* Find group in currently computed eigenvalues */
297:       EPSFindGroup(eps->nconv,ncv,eps->eigr,eps->eigi,eps->errest,grptol,&ngrp,&ctr,&ae,&arsd);
298:       if (ngrp!=nogrp) break;
299:       if (ngrp==0) break;
300:       if (PetscAbsScalar(ae-oae)>ctr*cnvtol*(itrsd[eps->nconv]-itrsdold[eps->nconv])) break;
301:       if (arsd>ctr*eps->tol) break;
302:       eps->nconv = eps->nconv + ngrp;
303:       if (eps->nconv>=ncv) break;
304:     }
305: 
306:     if (eps->nconv>=eps->nev) break;
307: 
308:     /* Compute nxtsrr (iteration of next projection step) */
309:     nxtsrr = PetscMin(eps->max_it,PetscMax((int)floor(stpfac*eps->its), init));
310: 
311:     if (ngrp!=nogrp || ngrp==0 || arsd>=oarsd) {
312:       idsrr = nxtsrr - eps->its;
313:     } else {
314:       idsrr = (int)floor(alpha+beta*(itrsdold[eps->nconv]-itrsd[eps->nconv])*log(arsd/eps->tol)/log(arsd/oarsd));
315:       idsrr = PetscMax(1,idsrr);
316:     }
317:     nxtsrr = PetscMin(nxtsrr,eps->its+idsrr);

319:     /* Compute nxtort (iteration of next orthogonalization step) */
320:     PetscMemcpy(U,T,sizeof(PetscScalar)*ncv);
321:     EPSHessCond(U,ncv,&tcond);
322:     idort = PetscMax(1,(int)floor(orttol/PetscMax(1,log10(tcond))));
323:     nxtort = PetscMin(eps->its+idort, nxtsrr);

325:     /* V(:,idx) = AV(:,idx) */
326:     for (i=eps->nconv;i<ncv;i++) {
327:       VecCopy(eps->AV[i],eps->V[i]);
328:     }
329:     eps->its++;

331:     /* Orthogonalization loop */
332:     do {
333:       while (eps->its<nxtort) {
334: 
335:         /* AV(:,idx) = OP * V(:,idx) */
336:         for (i=eps->nconv;i<ncv;i++) {
337:           STApply(eps->OP,eps->V[i],eps->AV[i]);
338:         }
339: 
340:         /* V(:,idx) = AV(:,idx) with normalization */
341:         for (i=eps->nconv;i<ncv;i++) {
342:           VecCopy(eps->AV[i],eps->V[i]);
343:           VecNorm(eps->V[i],NORM_INFINITY,&norm);
344:           t = 1 / norm;
345:           VecScale(&t,eps->V[i]);
346:         }
347: 
348:         eps->its++;
349:       }
350:       /* Orthonormalize vectors */
351:       for (i=eps->nconv;i<ncv;i++) {
352:         EPSOrthogonalize(eps,i+eps->nds,eps->DSV,eps->V[i],PETSC_NULL,&norm,&breakdown);
353:         if (breakdown) {
354:           SlepcVecSetRandom(eps->V[i]);
355:           EPSOrthogonalize(eps,i+eps->nds,eps->DSV,eps->V[i],PETSC_NULL,&norm,&breakdown);
356:         }
357:         t = 1 / norm;
358:         VecScale(&t,eps->V[i]);
359:       }
360:       nxtort = PetscMin(eps->its+idort,nxtsrr);
361:     } while (eps->its<nxtsrr);
362:   }

364:   PetscFree(U);
365:   PetscFree(rsdold);
366:   PetscFree(tau);
367:   PetscFree(work);
368:   PetscFree(itrsd);
369:   PetscFree(itrsdold);
370:   PetscFree(iwork);

372:   if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
373:   else eps->reason = EPS_DIVERGED_ITS;

375:   return(0);
376: }

378: EXTERN_C_BEGIN
381: PetscErrorCode EPSCreate_SUBSPACE(EPS eps)
382: {
384:   eps->ops->solve                = EPSSolve_SUBSPACE;
385:   eps->ops->setup                = EPSSetUp_SUBSPACE;
386:   eps->ops->destroy              = EPSDestroy_Default;
387:   eps->ops->backtransform        = EPSBackTransform_Default;
388:   eps->ops->computevectors       = EPSComputeVectors_Schur;
389:   return(0);
390: }
391: EXTERN_C_END