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-2010, 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

 42: PetscErrorCode EPSSolve_SUBSPACE(EPS);

 44: typedef struct {
 45:   Vec *AV;
 46: } EPS_SUBSPACE;

 50: PetscErrorCode EPSSetUp_SUBSPACE(EPS eps)
 51: {
 53:   PetscInt       i;
 54:   EPS_SUBSPACE   *ctx = (EPS_SUBSPACE *)eps->data;
 55:   PetscScalar    *pAV;

 58:   if (eps->ncv) { /* ncv set */
 59:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 60:   }
 61:   else if (eps->mpd) { /* mpd set */
 62:     eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
 63:   }
 64:   else { /* neither set: defaults depend on nev being small or large */
 65:     if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
 66:     else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
 67:   }
 68:   if (!eps->mpd) eps->mpd = eps->ncv;
 69:   if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 70:   if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
 71:   if (eps->which!=EPS_LARGEST_MAGNITUDE)
 72:     SETERRQ(1,"Wrong value of eps->which");
 73:   if (!eps->extraction) {
 74:     EPSSetExtraction(eps,EPS_RITZ);
 75:   } else if (eps->extraction!=EPS_RITZ) {
 76:     SETERRQ(PETSC_ERR_SUP,"Unsupported extraction type\n");
 77:   }

 79:   EPSAllocateSolution(eps);
 80:   PetscMalloc(eps->ncv*sizeof(Vec),&ctx->AV);
 81:   PetscMalloc(eps->ncv*eps->nloc*sizeof(PetscScalar),&pAV);
 82:   for (i=0;i<eps->ncv;i++) {
 83:     VecCreateMPIWithArray(((PetscObject)eps)->comm,eps->nloc,PETSC_DECIDE,pAV+i*eps->nloc,&ctx->AV[i]);
 84:   }
 85:   PetscFree(eps->T);
 86:   PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
 87:   EPSDefaultGetWork(eps,1);

 89:   /* dispatch solve method */
 90:   if (eps->leftvecs) SETERRQ(PETSC_ERR_SUP,"Left vectors not supported in this solver");
 91:   eps->ops->solve = EPSSolve_SUBSPACE;
 92:   return(0);
 93: }

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

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

143:    ngrp - number of entries of the group
144:    ctr  - (w(l)+w(l+ngrp-1))/2
145:    ae   - average of wr(l),...,wr(l+ngrp-1)
146:    arsd - average of rsd(l),...,rsd(l+ngrp-1)
147: */
148: static PetscErrorCode EPSFindGroup(PetscInt l,PetscInt m,PetscScalar *wr,PetscScalar *wi,PetscReal *rsd,
149:   PetscReal grptol,PetscInt *ngrp,PetscReal *ctr,PetscReal *ae,PetscReal *arsd)
150: {
151:   PetscInt  i;
152:   PetscReal rmod,rmod1;

155:   *ngrp = 0;
156:   *ctr = 0;
157: 
158:   rmod = SlepcAbsEigenvalue(wr[l],wi[l]);

160:   for (i=l;i<m;) {
161:     rmod1 = SlepcAbsEigenvalue(wr[i],wi[i]);
162:     if (PetscAbsReal(rmod-rmod1) > grptol*(rmod+rmod1)) break;
163:     *ctr = (rmod+rmod1)/2.0;
164:     if (wi[i] != 0.0) {
165:       (*ngrp)+=2;
166:       i+=2;
167:     } else {
168:       (*ngrp)++;
169:       i++;
170:     }
171:   }

173:   *ae = 0;
174:   *arsd = 0;

176:   if (*ngrp) {
177:     for (i=l;i<l+*ngrp;i++) {
178:       (*ae) += PetscRealPart(wr[i]);
179:       (*arsd) += rsd[i]*rsd[i];
180:     }
181:     *ae = *ae / *ngrp;
182:     *arsd = PetscSqrtScalar(*arsd / *ngrp);
183:   }
184:   return(0);
185: }

189: /*
190:    EPSSchurResidualNorms - Computes the column norms of residual vectors
191:    OP*V(1:n,l:m) - V*T(1:m,l:m), where, on entry, OP*V has been computed and 
192:    stored in AV. ldt is the leading dimension of T. On exit, rsd(l) to
193:    rsd(m) contain the computed norms.
194: */
195: static PetscErrorCode EPSSchurResidualNorms(EPS eps,Vec *V,Vec *AV,PetscScalar *T,PetscInt l,PetscInt m,PetscInt ldt,PetscReal *rsd)
196: {
198:   PetscInt       i,k;
199: #if defined(PETSC_USE_COMPLEX)
200:   PetscScalar    t;
201: #endif

204:   for (i=l;i<m;i++) {
205:     if (i==m-1 || T[i+1+ldt*i]==0.0) k=i+1; else k=i+2;
206:     VecCopy(AV[i],eps->work[0]);
207:     SlepcVecMAXPBY(eps->work[0],1.0,-1.0,k,T+ldt*i,V);
208: #if !defined(PETSC_USE_COMPLEX)
209:     VecDot(eps->work[0],eps->work[0],rsd+i);
210: #else
211:     VecDot(eps->work[0],eps->work[0],&t);
212:     rsd[i] = PetscRealPart(t);
213: #endif    
214:   }

216:   for (i=l;i<m;i++) {
217:     if (i == m-1) {
218:       rsd[i] = sqrt(rsd[i]);
219:     } else if (T[i+1+(ldt*i)]==0.0) {
220:       rsd[i] = sqrt(rsd[i]);
221:     } else {
222:       rsd[i] = sqrt((rsd[i]+rsd[i+1])/2.0);
223:       rsd[i+1] = rsd[i];
224:       i++;
225:     }
226:   }
227:   return(0);
228: }

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

252:   its = 0;
253:   PetscMalloc(sizeof(PetscScalar)*ncv*ncv,&U);
254:   PetscMalloc(sizeof(PetscReal)*ncv,&rsd);
255:   PetscMalloc(sizeof(PetscInt)*ncv,&itrsd);
256:   PetscMalloc(sizeof(PetscInt)*ncv,&itrsdold);

258:   for (i=0;i<ncv;i++) {
259:     rsd[i] = 0.0;
260:     itrsd[i] = -1;
261:   }
262: 
263:   /* Complete the initial basis with random vectors and orthonormalize them */
264:   k = eps->nini;
265:   while (k<ncv) {
266:     SlepcVecSetRandom(eps->V[k],eps->rand);
267:     IPOrthogonalize(eps->ip,eps->nds,eps->DS,k,PETSC_NULL,eps->V,eps->V[k],PETSC_NULL,&norm,&breakdown);
268:     if (norm>0.0 && !breakdown) {
269:       VecScale(eps->V[k],1.0/norm);
270:       k++;
271:     }
272:   }

274:   while (eps->its<eps->max_it) {
275:     eps->its++;
276:     nv = PetscMin(eps->nconv+eps->mpd,ncv);
277: 
278:     /* Find group in previously computed eigenvalues */
279:     EPSFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,rsd,grptol,&nogrp,&octr,&oae,&oarsd);

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

284:     /* 1. AV(:,idx) = OP * V(:,idx) */
285:     for (i=eps->nconv;i<nv;i++) {
286:       STApply(eps->OP,eps->V[i],ctx->AV[i]);
287:     }

289:     /* 2. T(:,idx) = V' * AV(:,idx) */
290:     for (i=eps->nconv;i<nv;i++) {
291:       VecMDot(ctx->AV[i],nv,eps->V,T+i*ncv);
292:     }

294:     /* 3. Reduce projected matrix to Hessenberg form: [U,T] = hess(T) */
295:     EPSDenseHessenberg(nv,eps->nconv,T,ncv,U);
296: 
297:     /* 4. Reduce T to quasi-triangular (Schur) form */
298:     EPSDenseSchur(nv,eps->nconv,T,ncv,U,eps->eigr,eps->eigi);

300:     /* 5. Sort diagonal elements in T and accumulate rotations on U */
301:     EPSSortDenseSchur(eps,nv,eps->nconv,T,ncv,U,eps->eigr,eps->eigi);
302: 
303:     /* 6. AV(:,idx) = AV * U(:,idx) */
304:     SlepcUpdateVectors(nv,ctx->AV,eps->nconv,nv,U,nv,PETSC_FALSE);
305: 
306:     /* 7. V(:,idx) = V * U(:,idx) */
307:     SlepcUpdateVectors(nv,eps->V,eps->nconv,nv,U,nv,PETSC_FALSE);
308: 
309:     /* Convergence check */
310:     EPSSchurResidualNorms(eps,eps->V,ctx->AV,T,eps->nconv,nv,ncv,rsd);

312:     for (i=eps->nconv;i<nv;i++) {
313:       itrsdold[i] = itrsd[i];
314:       itrsd[i] = its;
315:       eps->errest[i] = rsd[i];
316:     }
317: 
318:     for (;;) {
319:       /* Find group in currently computed eigenvalues */
320:       EPSFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,eps->errest,grptol,&ngrp,&ctr,&ae,&arsd);
321:       if (ngrp!=nogrp) break;
322:       if (ngrp==0) break;
323:       if (PetscAbsScalar(ae-oae)>ctr*cnvtol*(itrsd[eps->nconv]-itrsdold[eps->nconv])) break;
324:       if (arsd>ctr*eps->tol) break;
325:       eps->nconv = eps->nconv + ngrp;
326:       if (eps->nconv>=nv) break;
327:     }
328: 
329:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
330:     if (eps->nconv>=eps->nev) break;
331: 
332:     /* Compute nxtsrr (iteration of next projection step) */
333:     nxtsrr = PetscMin(eps->max_it,PetscMax((PetscInt)floor(stpfac*its), init));
334: 
335:     if (ngrp!=nogrp || ngrp==0 || arsd>=oarsd) {
336:       idsrr = nxtsrr - its;
337:     } else {
338:       idsrr = (PetscInt)floor(alpha+beta*(itrsdold[eps->nconv]-itrsd[eps->nconv])*log(arsd/eps->tol)/log(arsd/oarsd));
339:       idsrr = PetscMax(1,idsrr);
340:     }
341:     nxtsrr = PetscMin(nxtsrr,its+idsrr);

343:     /* Compute nxtort (iteration of next orthogonalization step) */
344:     PetscMemcpy(U,T,sizeof(PetscScalar)*ncv*ncv);
345:     EPSHessCond(nv,U,ncv,&tcond);
346:     idort = PetscMax(1,(PetscInt)floor(orttol/PetscMax(1,log10(tcond))));
347:     nxtort = PetscMin(its+idort, nxtsrr);

349:     /* V(:,idx) = AV(:,idx) */
350:     for (i=eps->nconv;i<nv;i++) {
351:       VecCopy(ctx->AV[i],eps->V[i]);
352:     }
353:     its++;

355:     /* Orthogonalization loop */
356:     do {
357:       while (its<nxtort) {
358: 
359:         /* AV(:,idx) = OP * V(:,idx) */
360:         for (i=eps->nconv;i<nv;i++) {
361:           STApply(eps->OP,eps->V[i],ctx->AV[i]);
362:         }
363: 
364:         /* V(:,idx) = AV(:,idx) with normalization */
365:         for (i=eps->nconv;i<nv;i++) {
366:           VecCopy(ctx->AV[i],eps->V[i]);
367:           VecNorm(eps->V[i],NORM_INFINITY,&norm);
368:           VecScale(eps->V[i],1/norm);
369:         }
370: 
371:         its++;
372:       }
373:       /* Orthonormalize vectors */
374:       for (i=eps->nconv;i<nv;i++) {
375:         IPOrthogonalize(eps->ip,eps->nds,eps->DS,i,PETSC_NULL,eps->V,eps->V[i],PETSC_NULL,&norm,&breakdown);
376:         if (breakdown) {
377:           SlepcVecSetRandom(eps->V[i],eps->rand);
378:           IPOrthogonalize(eps->ip,eps->nds,eps->DS,i,PETSC_NULL,eps->V,eps->V[i],PETSC_NULL,&norm,&breakdown);
379:         }
380:         VecScale(eps->V[i],1/norm);
381:       }
382:       nxtort = PetscMin(its+idort,nxtsrr);
383:     } while (its<nxtsrr);
384:   }

386:   PetscFree(U);
387:   PetscFree(rsd);
388:   PetscFree(itrsd);
389:   PetscFree(itrsdold);

391:   if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
392:   else eps->reason = EPS_DIVERGED_ITS;

394:   return(0);
395: }

399: PetscErrorCode EPSDestroy_SUBSPACE(EPS eps)
400: {
402:   PetscInt       i;
403:   PetscScalar    *pAV;
404:   EPS_SUBSPACE   *ctx = (EPS_SUBSPACE *)eps->data;

407:   VecGetArray(ctx->AV[0],&pAV);
408:   for (i=0;i<eps->ncv;i++) {
409:     VecDestroy(ctx->AV[i]);
410:   }
411:   PetscFree(pAV);
412:   PetscFree(ctx->AV);
413:   EPSDestroy_Default(eps);
414:   return(0);
415: }

420: PetscErrorCode EPSCreate_SUBSPACE(EPS eps)
421: {
423:   EPS_SUBSPACE   *ctx;

426:   PetscNew(EPS_SUBSPACE,&ctx);
427:   PetscLogObjectMemory(eps,sizeof(EPS_SUBSPACE));
428:   eps->data                      = (void *) ctx;
429:   eps->ops->setup                = EPSSetUp_SUBSPACE;
430:   eps->ops->destroy              = EPSDestroy_SUBSPACE;
431:   eps->ops->backtransform        = EPSBackTransform_Default;
432:   eps->ops->computevectors       = EPSComputeVectors_Schur;
433:   return(0);
434: }