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:  #include src/eps/epsimpl.h
 21:  #include slepcblaslapack.h

 25: PetscErrorCode EPSSetUp_SUBSPACE(EPS eps)
 26: {
 28:   PetscInt       N;

 31:   VecGetSize(eps->vec_initial,&N);
 32:   if (eps->nev > N) eps->nev = N;
 33:   if (eps->ncv) {
 34:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 35:   }
 36:   else eps->ncv = PetscMin(N,PetscMax(2*eps->nev,eps->nev+15));
 37:   if (!eps->max_it) eps->max_it = PetscMax(100,N);
 38:   if (!eps->tol) eps->tol = 1.e-7;
 39:   if (eps->which!=EPS_LARGEST_MAGNITUDE)
 40:     SETERRQ(1,"Wrong value of eps->which");
 41:   EPSAllocateSolution(eps);
 42:   PetscFree(eps->T);
 43:   PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
 44:   EPSDefaultGetWork(eps,eps->ncv);
 45:   return(0);
 46: }

 50: /*
 51:    EPSHessCond - Compute the inf-norm condition number of the upper 
 52:    Hessenberg matrix H: cond(H) = norm(H)*norm(inv(H)).
 53:    This routine uses Gaussian elimination with partial pivoting to 
 54:    compute the inverse explicitly. 
 55: */
 56: static PetscErrorCode EPSHessCond(PetscScalar* H,int n, PetscReal* cond)
 57: {
 58: #if defined(PETSC_MISSING_LAPACK_GETRF) || defined(SLEPC_MISSING_LAPACK_GETRI) || defined(SLEPC_MISSING_LAPACK_LANGE) || defined(SLEPC_MISSING_LAPACK_LANHS)
 60:   SETERRQ(PETSC_ERR_SUP,"GETRF,GETRI - Lapack routines are unavailable.");
 61: #else
 63:   int            *ipiv,lwork,info;
 64:   PetscScalar    *work;
 65:   PetscReal      hn,hin,*rwork;
 66: 
 68:   PetscMalloc(sizeof(int)*n,&ipiv);
 69:   lwork = n*n;
 70:   PetscMalloc(sizeof(PetscScalar)*lwork,&work);
 71:   PetscMalloc(sizeof(PetscReal)*n,&rwork);
 72:   hn = LAPACKlanhs_("I",&n,H,&n,rwork,1);
 73:   LAPACKgetrf_(&n,&n,H,&n,ipiv,&info);
 74:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGETRF %d",info);
 75:   LAPACKgetri_(&n,H,&n,ipiv,work,&lwork,&info);
 76:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGETRI %d",info);
 77:   hin = LAPACKlange_("I",&n,&n,H,&n,rwork,1);
 78:   *cond = hn * hin;
 79:   PetscFree(ipiv);
 80:   PetscFree(work);
 81:   PetscFree(rwork);
 82:   return(0);
 83: #endif
 84: }

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

 94:    ngrp - number of entries of the group
 95:    ctr  - (w(l)+w(l+ngrp-1))/2
 96:    ae   - average of wr(l),...,wr(l+ngrp-1)
 97:    arsd - average of rsd(l),...,rsd(l+ngrp-1)
 98: */
 99: static PetscErrorCode EPSFindGroup(int l,int m,PetscScalar *wr,PetscScalar *wi,PetscReal *rsd,
100:   PetscReal grptol,int *ngrp,PetscReal *ctr,PetscReal *ae,PetscReal *arsd)
101: {
102:   int       i;
103:   PetscReal rmod,rmod1;

106:   *ngrp = 0;
107:   *ctr = 0;
108: 
109:   rmod = SlepcAbsEigenvalue(wr[l],wi[l]);

111:   for (i=l;i<m;) {
112:     rmod1 = SlepcAbsEigenvalue(wr[i],wi[i]);
113:     if (PetscAbsReal(rmod-rmod1) > grptol*(rmod+rmod1)) break;
114:     *ctr = (rmod+rmod1)/2.0;
115:     if (wi[i] != 0.0) {
116:       (*ngrp)+=2;
117:       i+=2;
118:     } else {
119:       (*ngrp)++;
120:       i++;
121:     }
122:   }

124:   *ae = 0;
125:   *arsd = 0;

127:   if (*ngrp) {
128:     for (i=l;i<l+*ngrp;i++) {
129:       (*ae) += PetscRealPart(wr[i]);
130:       (*arsd) += rsd[i]*rsd[i];
131:     }
132:     *ae = *ae / *ngrp;
133:     *arsd = PetscSqrtScalar(*arsd / *ngrp);
134:   }
135:   return(0);
136: }

140: /*
141:    EPSSchurResidualNorms - Computes the column norms of residual vectors
142:    OP*V(1:n,l:m) - V*T(1:m,l:m) were on entry, OP*V has been computed and 
143:    stored in AV. ldt is the leading dimension of T. On exit, rsd(l) to
144:    rsd(m) contain the computed norms.
145: */
146: static PetscErrorCode EPSSchurResidualNorms(EPS eps,Vec *V,Vec *AV,PetscScalar *T,int l,int m,int ldt,PetscReal *rsd)
147: {
149:   int            i;
150: #if defined(PETSC_USE_COMPLEX)
151:   PetscScalar    t;
152: #endif

155:   for (i=l;i<m;i++) {
156:     VecSet(eps->work[0],0.0);
157:     VecMAXPY(eps->work[0],m,T+ldt*i,V);
158:     VecWAXPY(eps->work[1],-1.0,eps->work[0],AV[i]);
159: #if !defined(PETSC_USE_COMPLEX)
160:     VecDot(eps->work[1],eps->work[1],rsd+i);
161: #else
162:     VecDot(eps->work[1],eps->work[1],&t);
163:     rsd[i] = PetscRealPart(t);
164: #endif    
165:   }

167:   for (i=l;i<m;i++) {
168:     if (i == m-1) {
169:       rsd[i] = sqrt(rsd[i]);
170:     } else if (T[i+1+(ldt*i)]==0.0) {
171:       rsd[i] = sqrt(rsd[i]);
172:     } else {
173:       rsd[i] = sqrt(rsd[i]+rsd[i+1])/2.0;
174:       rsd[i+1] = rsd[i];
175:       i++;
176:     }
177:   }
178:   return(0);
179: }

183: PetscErrorCode EPSSolve_SUBSPACE(EPS eps)
184: {
185: #if defined(SLEPC_MISSING_LAPACK_GEHRD) || defined(SLEPC_MISSING_LAPACK_ORGHR) || defined(SLEPC_MISSING_LAPACK_UNGHR)
186:   SETERRQ(PETSC_ERR_SUP,"GEHRD,ORGHR/UNGHR - Lapack routines are unavailable.");
187: #else
189:   int            i,j,ilo,lwork,info,ngrp,nogrp,*itrsd,*itrsdold,
190:                  nxtsrr,idsrr,*iwork,idort,nxtort,ncv = eps->ncv;
191:   PetscScalar    *T=eps->T,*U,*tau,*work;
192:   PetscReal      arsd,oarsd,ctr,octr,ae,oae,*rsd,*rsdold,norm,tcond;
193:   PetscTruth     breakdown;
194:   /* Parameters */
195:   int            init = 5;        /* Number of initial iterations */
196:   PetscReal      stpfac = 1.5,    /* Max num of iter before next SRR step */
197:                  alpha = 1.0,     /* Used to predict convergence of next residual */
198:                  beta = 1.1,      /* Used to predict convergence of next residual */
199:                  grptol = 1e-8,   /* Tolerance for EPSFindGroup */
200:                  cnvtol = 1e-6;   /* Convergence criterion for cnv */
201:   int            orttol = 2;      /* Number of decimal digits whose loss
202:                                      can be tolerated in orthogonalization */

205:   eps->its = 0;
206:   eps->nconv = 0;
207:   PetscMalloc(sizeof(PetscScalar)*ncv*ncv,&U);
208:   PetscMalloc(sizeof(PetscReal)*ncv,&rsd);
209:   PetscMalloc(sizeof(PetscReal)*ncv,&rsdold);
210:   PetscMalloc(sizeof(PetscScalar)*ncv,&tau);
211:   lwork = ncv*ncv;
212:   PetscMalloc(sizeof(PetscScalar)*lwork,&work);
213:   PetscMalloc(sizeof(int)*ncv,&itrsd);
214:   PetscMalloc(sizeof(int)*ncv,&itrsdold);
215:   PetscMalloc(sizeof(int)*ncv,&iwork);

217:   /* Generate a set of random initial vectors and orthonormalize them */
218:   for (i=0;i<ncv;i++) {
219:     SlepcVecSetRandom(eps->V[i]);
220:     eps->eigr[i] = 0.0;
221:     eps->eigi[i] = 0.0;
222:     rsd[i] = 0.0;
223:     itrsd[i] = -1;
224:   }
225:   EPSQRDecomposition(eps,eps->V,0,ncv,PETSC_NULL,0);
226: 
227:   while (eps->its<eps->max_it) {

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

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

235:     /* 1. AV(:,idx) = OP * V(:,idx) */
236:     for (i=eps->nconv;i<ncv;i++) {
237:       STApply(eps->OP,eps->V[i],eps->AV[i]);
238:     }

240:     /* 2. T(:,idx) = V' * AV(:,idx) */
241:     for (i=eps->nconv;i<ncv;i++) {
242:       VecMDot(eps->AV[i],ncv,eps->V,T+i*ncv);
243:     }

245:     /* 3. Reduce projected matrix to Hessenberg form: [U,T] = hess(T) */
246:     ilo = eps->nconv + 1;
247:     LAPACKgehrd_(&ncv,&ilo,&ncv,T,&ncv,tau,work,&lwork,&info);
248:     if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGEHRD %d",info);
249:     for (j=0;j<ncv-1;j++) {
250:       for (i=j+2;i<ncv;i++) {
251:         U[i+j*ncv] = T[i+j*ncv];
252:         T[i+j*ncv] = 0.0;
253:       }
254:     }
255:     LAPACKorghr_(&ncv,&ilo,&ncv,U,&ncv,tau,work,&lwork,&info);
256:     if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORGHR %d",info);
257: 
258:     /* 4. Reduce T to quasi-triangular (Schur) form */
259:     EPSDenseSchur(ncv,eps->nconv,T,ncv,U,eps->eigr,eps->eigi);

261:     /* 5. Sort diagonal elements in T and accumulate rotations on U */
262:     EPSSortDenseSchur(ncv,eps->nconv,T,ncv,U,eps->eigr,eps->eigi);
263: 
264:     /* 6. AV(:,idx) = AV * U(:,idx) */
265:     for (i=eps->nconv;i<ncv;i++) {
266:       VecSet(eps->work[i],0.0);
267:       VecMAXPY(eps->work[i],ncv,U+ncv*i,eps->AV);
268:     }
269:     for (i=eps->nconv;i<ncv;i++) {
270:       VecCopy(eps->work[i],eps->AV[i]);
271:     }
272: 
273:     /* 7. V(:,idx) = V * U(:,idx) */
274:     for (i=eps->nconv;i<ncv;i++) {
275:       VecSet(eps->work[i],0.0);
276:       VecMAXPY(eps->work[i],ncv,U+ncv*i,eps->V);
277:     }
278:     for (i=eps->nconv;i<ncv;i++) {
279:       VecCopy(eps->work[i],eps->V[i]);
280:     }
281: 
282:     /* Compute residuals */
283:     for (i=0;i<ncv;i++) { rsdold[i] = rsd[i]; }

285:     EPSSchurResidualNorms(eps,eps->V,eps->AV,T,eps->nconv,ncv,ncv,rsd);

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

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

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

332:     /* Orthogonalization loop */
333:     do {
334:       while (eps->its<nxtort) {
335: 
336:         /* AV(:,idx) = OP * V(:,idx) */
337:         for (i=eps->nconv;i<ncv;i++) {
338:           STApply(eps->OP,eps->V[i],eps->AV[i]);
339:         }
340: 
341:         /* V(:,idx) = AV(:,idx) with normalization */
342:         for (i=eps->nconv;i<ncv;i++) {
343:           VecCopy(eps->AV[i],eps->V[i]);
344:           VecNorm(eps->V[i],NORM_INFINITY,&norm);
345:           VecScale(eps->V[i],1/norm);
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:         VecScale(eps->V[i],1/norm);
358:       }
359:       nxtort = PetscMin(eps->its+idort,nxtsrr);
360:     } while (eps->its<nxtsrr);
361:   }

363:   PetscFree(U);
364:   PetscFree(rsd);
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: #endif 
377: }

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