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-2009, Universidad 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
 40:  #include slepcblaslapack.h

 44: PetscErrorCode EPSSetUp_SUBSPACE(EPS eps)
 45: {
 47:   PetscInt       i,N,nloc;
 48:   PetscScalar    *pAV;

 51:   VecGetSize(eps->vec_initial,&N);
 52:   if (eps->ncv) { /* ncv set */
 53:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 54:   }
 55:   else if (eps->mpd) { /* mpd set */
 56:     eps->ncv = PetscMin(N,eps->nev+eps->mpd);
 57:   }
 58:   else { /* neither set: defaults depend on nev being small or large */
 59:     if (eps->nev<500) eps->ncv = PetscMin(N,PetscMax(2*eps->nev,eps->nev+15));
 60:     else { eps->mpd = 500; eps->ncv = PetscMin(N,eps->nev+eps->mpd); }
 61:   }
 62:   if (!eps->mpd) eps->mpd = eps->ncv;
 63:   if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);
 64:   if (eps->which!=EPS_LARGEST_MAGNITUDE)
 65:     SETERRQ(1,"Wrong value of eps->which");
 66:   if (!eps->extraction) {
 67:     EPSSetExtraction(eps,EPS_RITZ);
 68:   } else if (eps->extraction!=EPS_RITZ) {
 69:     SETERRQ(PETSC_ERR_SUP,"Unsupported extraction type\n");
 70:   }

 72:   EPSAllocateSolution(eps);
 73:   VecGetLocalSize(eps->vec_initial,&nloc);
 74:   PetscMalloc(eps->ncv*sizeof(Vec),&eps->AV);
 75:   PetscMalloc(eps->ncv*nloc*sizeof(PetscScalar),&pAV);
 76:   for (i=0;i<eps->ncv;i++) {
 77:     VecCreateMPIWithArray(((PetscObject)eps)->comm,nloc,PETSC_DECIDE,pAV+i*nloc,&eps->AV[i]);
 78:   }
 79:   PetscFree(eps->T);
 80:   PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
 81:   EPSDefaultGetWork(eps,2);
 82:   return(0);
 83: }

 87: /*
 88:    EPSHessCond - Compute the inf-norm condition number of the upper 
 89:    Hessenberg matrix H: cond(H) = norm(H)*norm(inv(H)).
 90:    This routine uses Gaussian elimination with partial pivoting to 
 91:    compute the inverse explicitly. 
 92: */
 93: static PetscErrorCode EPSHessCond(PetscInt n_,PetscScalar* H,PetscInt ldh_,PetscReal* cond)
 94: {
 95: #if defined(PETSC_MISSING_LAPACK_GETRF) || defined(SLEPC_MISSING_LAPACK_GETRI) || defined(SLEPC_MISSING_LAPACK_LANGE) || defined(SLEPC_MISSING_LAPACK_LANHS)
 97:   SETERRQ(PETSC_ERR_SUP,"GETRF,GETRI - Lapack routines are unavailable.");
 98: #else
100:   PetscBLASInt   *ipiv,lwork,info,n=n_,ldh=ldh_;
101:   PetscScalar    *work;
102:   PetscReal      hn,hin,*rwork;
103: 
105:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
106:   PetscMalloc(sizeof(PetscBLASInt)*n,&ipiv);
107:   lwork = n*n;
108:   PetscMalloc(sizeof(PetscScalar)*lwork,&work);
109:   PetscMalloc(sizeof(PetscReal)*n,&rwork);
110:   hn = LAPACKlanhs_("I",&n,H,&ldh,rwork);
111:   LAPACKgetrf_(&n,&n,H,&ldh,ipiv,&info);
112:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGETRF %d",info);
113:   LAPACKgetri_(&n,H,&ldh,ipiv,work,&lwork,&info);
114:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGETRI %d",info);
115:   hin = LAPACKlange_("I",&n,&n,H,&ldh,rwork);
116:   *cond = hn * hin;
117:   PetscFree(ipiv);
118:   PetscFree(work);
119:   PetscFree(rwork);
120:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
121:   return(0);
122: #endif
123: }

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

133:    ngrp - number of entries of the group
134:    ctr  - (w(l)+w(l+ngrp-1))/2
135:    ae   - average of wr(l),...,wr(l+ngrp-1)
136:    arsd - average of rsd(l),...,rsd(l+ngrp-1)
137: */
138: static PetscErrorCode EPSFindGroup(PetscInt l,PetscInt m,PetscScalar *wr,PetscScalar *wi,PetscReal *rsd,
139:   PetscReal grptol,PetscInt *ngrp,PetscReal *ctr,PetscReal *ae,PetscReal *arsd)
140: {
141:   PetscInt  i;
142:   PetscReal rmod,rmod1;

145:   *ngrp = 0;
146:   *ctr = 0;
147: 
148:   rmod = SlepcAbsEigenvalue(wr[l],wi[l]);

150:   for (i=l;i<m;) {
151:     rmod1 = SlepcAbsEigenvalue(wr[i],wi[i]);
152:     if (PetscAbsReal(rmod-rmod1) > grptol*(rmod+rmod1)) break;
153:     *ctr = (rmod+rmod1)/2.0;
154:     if (wi[i] != 0.0) {
155:       (*ngrp)+=2;
156:       i+=2;
157:     } else {
158:       (*ngrp)++;
159:       i++;
160:     }
161:   }

163:   *ae = 0;
164:   *arsd = 0;

166:   if (*ngrp) {
167:     for (i=l;i<l+*ngrp;i++) {
168:       (*ae) += PetscRealPart(wr[i]);
169:       (*arsd) += rsd[i]*rsd[i];
170:     }
171:     *ae = *ae / *ngrp;
172:     *arsd = PetscSqrtScalar(*arsd / *ngrp);
173:   }
174:   return(0);
175: }

179: /*
180:    EPSSchurResidualNorms - Computes the column norms of residual vectors
181:    OP*V(1:n,l:m) - V*T(1:m,l:m), where, on entry, OP*V has been computed and 
182:    stored in AV. ldt is the leading dimension of T. On exit, rsd(l) to
183:    rsd(m) contain the computed norms.
184: */
185: static PetscErrorCode EPSSchurResidualNorms(EPS eps,Vec *V,Vec *AV,PetscScalar *T,PetscInt l,PetscInt m,PetscInt ldt,PetscReal *rsd)
186: {
188:   PetscInt       i,k;
189: #if defined(PETSC_USE_COMPLEX)
190:   PetscScalar    t;
191: #endif

194:   for (i=l;i<m;i++) {
195:     VecSet(eps->work[0],0.0);
196:     if (i==m-1 || T[i+1+ldt*i]==0.0) k=i+1; else k=i+2;
197:     VecMAXPY(eps->work[0],k,T+ldt*i,V);
198:     VecWAXPY(eps->work[1],-1.0,eps->work[0],AV[i]);
199: #if !defined(PETSC_USE_COMPLEX)
200:     VecDot(eps->work[1],eps->work[1],rsd+i);
201: #else
202:     VecDot(eps->work[1],eps->work[1],&t);
203:     rsd[i] = PetscRealPart(t);
204: #endif    
205:   }

207:   for (i=l;i<m;i++) {
208:     if (i == m-1) {
209:       rsd[i] = sqrt(rsd[i]);
210:     } else if (T[i+1+(ldt*i)]==0.0) {
211:       rsd[i] = sqrt(rsd[i]);
212:     } else {
213:       rsd[i] = sqrt((rsd[i]+rsd[i+1])/2.0);
214:       rsd[i+1] = rsd[i];
215:       i++;
216:     }
217:   }
218:   return(0);
219: }

223: PetscErrorCode EPSSolve_SUBSPACE(EPS eps)
224: {
226:   PetscInt       i,ngrp,nogrp,*itrsd,*itrsdold,
227:                  nxtsrr,idsrr,idort,nxtort,nv,ncv = eps->ncv,its;
228:   PetscScalar    *T=eps->T,*U;
229:   PetscReal      arsd,oarsd,ctr,octr,ae,oae,*rsd,*rsdold,norm,tcond=1.0;
230:   PetscTruth     breakdown;
231:   /* Parameters */
232:   PetscInt       init = 5;        /* Number of initial iterations */
233:   PetscReal      stpfac = 1.5,    /* Max num of iter before next SRR step */
234:                  alpha = 1.0,     /* Used to predict convergence of next residual */
235:                  beta = 1.1,      /* Used to predict convergence of next residual */
236:                  grptol = 1e-8,   /* Tolerance for EPSFindGroup */
237:                  cnvtol = 1e-6;   /* Convergence criterion for cnv */
238:   PetscInt       orttol = 2;      /* Number of decimal digits whose loss
239:                                      can be tolerated in orthogonalization */

242:   its = 0;
243:   PetscMalloc(sizeof(PetscScalar)*ncv*ncv,&U);
244:   PetscMalloc(sizeof(PetscReal)*ncv,&rsd);
245:   PetscMalloc(sizeof(PetscReal)*ncv,&rsdold);
246:   PetscMalloc(sizeof(PetscInt)*ncv,&itrsd);
247:   PetscMalloc(sizeof(PetscInt)*ncv,&itrsdold);

249:   /* Generate a set of random initial vectors and orthonormalize them */
250:   for (i=0;i<ncv;i++) {
251:     SlepcVecSetRandom(eps->V[i]);
252:     rsd[i] = 0.0;
253:     itrsd[i] = -1;
254:   }
255:   IPQRDecomposition(eps->ip,eps->V,0,ncv,PETSC_NULL,0,eps->work[0]);
256: 
257:   while (eps->its<eps->max_it) {
258:     eps->its++;
259:     nv = PetscMin(eps->nconv+eps->mpd,ncv);
260: 
261:     /* Find group in previously computed eigenvalues */
262:     EPSFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,rsd,grptol,&nogrp,&octr,&oae,&oarsd);

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

267:     /* 1. AV(:,idx) = OP * V(:,idx) */
268:     for (i=eps->nconv;i<nv;i++) {
269:       STApply(eps->OP,eps->V[i],eps->AV[i]);
270:     }

272:     /* 2. T(:,idx) = V' * AV(:,idx) */
273:     for (i=eps->nconv;i<nv;i++) {
274:       VecMDot(eps->AV[i],nv,eps->V,T+i*ncv);
275:     }

277:     /* 3. Reduce projected matrix to Hessenberg form: [U,T] = hess(T) */
278:     EPSDenseHessenberg(nv,eps->nconv,T,ncv,U);
279: 
280:     /* 4. Reduce T to quasi-triangular (Schur) form */
281:     EPSDenseSchur(nv,eps->nconv,T,ncv,U,eps->eigr,eps->eigi);

283:     /* 5. Sort diagonal elements in T and accumulate rotations on U */
284:     EPSSortDenseSchur(nv,eps->nconv,T,ncv,U,eps->eigr,eps->eigi,eps->which);
285: 
286:     /* 6. AV(:,idx) = AV * U(:,idx) */
287:     SlepcUpdateVectors(nv,eps->AV,eps->nconv,nv,U,nv,PETSC_FALSE);
288: 
289:     /* 7. V(:,idx) = V * U(:,idx) */
290:     SlepcUpdateVectors(nv,eps->V,eps->nconv,nv,U,nv,PETSC_FALSE);
291: 
292:     /* Compute residuals */
293:     for (i=0;i<nv;i++) { rsdold[i] = rsd[i]; }

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

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

330:     /* Compute nxtort (iteration of next orthogonalization step) */
331:     PetscMemcpy(U,T,sizeof(PetscScalar)*ncv*ncv);
332:     EPSHessCond(nv,U,ncv,&tcond);
333:     idort = PetscMax(1,(PetscInt)floor(orttol/PetscMax(1,log10(tcond))));
334:     nxtort = PetscMin(its+idort, nxtsrr);

336:     /* V(:,idx) = AV(:,idx) */
337:     for (i=eps->nconv;i<nv;i++) {
338:       VecCopy(eps->AV[i],eps->V[i]);
339:     }
340:     its++;

342:     /* Orthogonalization loop */
343:     do {
344:       while (its<nxtort) {
345: 
346:         /* AV(:,idx) = OP * V(:,idx) */
347:         for (i=eps->nconv;i<nv;i++) {
348:           STApply(eps->OP,eps->V[i],eps->AV[i]);
349:         }
350: 
351:         /* V(:,idx) = AV(:,idx) with normalization */
352:         for (i=eps->nconv;i<nv;i++) {
353:           VecCopy(eps->AV[i],eps->V[i]);
354:           VecNorm(eps->V[i],NORM_INFINITY,&norm);
355:           VecScale(eps->V[i],1/norm);
356:         }
357: 
358:         its++;
359:       }
360:       /* Orthonormalize vectors */
361:       for (i=eps->nconv;i<nv;i++) {
362:         IPOrthogonalize(eps->ip,i+eps->nds,PETSC_NULL,eps->DSV,eps->V[i],PETSC_NULL,&norm,&breakdown,eps->work[0],PETSC_NULL);
363:         if (breakdown) {
364:           SlepcVecSetRandom(eps->V[i]);
365:           IPOrthogonalize(eps->ip,i+eps->nds,PETSC_NULL,eps->DSV,eps->V[i],PETSC_NULL,&norm,&breakdown,eps->work[0],PETSC_NULL);
366:         }
367:         VecScale(eps->V[i],1/norm);
368:       }
369:       nxtort = PetscMin(its+idort,nxtsrr);
370:     } while (its<nxtsrr);
371:   }

373:   PetscFree(U);
374:   PetscFree(rsd);
375:   PetscFree(rsdold);
376:   PetscFree(itrsd);
377:   PetscFree(itrsdold);

379:   if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
380:   else eps->reason = EPS_DIVERGED_ITS;

382:   return(0);
383: }

387: PetscErrorCode EPSDestroy_SUBSPACE(EPS eps)
388: {
390:   PetscInt       i;
391:   PetscScalar    *pAV;

394:   VecGetArray(eps->AV[0],&pAV);
395:   for (i=0;i<eps->ncv;i++) {
396:     VecDestroy(eps->AV[i]);
397:   }
398:   PetscFree(pAV);
399:   PetscFree(eps->AV);
400:   EPSDestroy_Default(eps);
401:   return(0);
402: }

407: PetscErrorCode EPSCreate_SUBSPACE(EPS eps)
408: {
410:   eps->ops->solve                = EPSSolve_SUBSPACE;
411:   eps->ops->setup                = EPSSetUp_SUBSPACE;
412:   eps->ops->destroy              = EPSDestroy_SUBSPACE;
413:   eps->ops->backtransform        = EPSBackTransform_Default;
414:   eps->ops->computevectors       = EPSComputeVectors_Schur;
415:   return(0);
416: }