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-2012, Universitat 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 <slepc-private/epsimpl.h>                /*I "slepceps.h" I*/
 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:   EPS_SUBSPACE   *ctx = (EPS_SUBSPACE*)eps->data;

 56:   if (eps->ncv) { /* ncv set */
 57:     if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
 58:   }
 59:   else if (eps->mpd) { /* mpd set */
 60:     eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
 61:   }
 62:   else { /* neither set: defaults depend on nev being small or large */
 63:     if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
 64:     else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
 65:   }
 66:   if (!eps->mpd) eps->mpd = eps->ncv;
 67:   if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 68:   if (!eps->which) { EPSDefaultSetWhich(eps); }
 69:   if (eps->which!=EPS_LARGEST_MAGNITUDE && eps->which!=EPS_TARGET_MAGNITUDE) SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
 70:   if (!eps->extraction) {
 71:     EPSSetExtraction(eps,EPS_RITZ);
 72:   } else if (eps->extraction!=EPS_RITZ) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
 73:   if (eps->arbit_func) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");

 75:   EPSAllocateSolution(eps);
 76:   VecDuplicateVecs(eps->t,eps->ncv,&ctx->AV);
 77:   if (eps->ishermitian) {
 78:     DSSetType(eps->ds,DSHEP);
 79:   } else {
 80:     DSSetType(eps->ds,DSNHEP);
 81:   }
 82:   DSAllocate(eps->ds,eps->ncv);
 83:   EPSDefaultGetWork(eps,1);

 85:   /* dispatch solve method */
 86:   if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
 87:   eps->ops->solve = EPSSolve_Subspace;
 88:   return(0);
 89: }

 93: /*
 94:    EPSSubspaceFindGroup - Find a group of nearly equimodular eigenvalues, provided 
 95:    in arrays wr and wi, according to the tolerance grptol. Also the 2-norms
 96:    of the residuals must be passed in (rsd). Arrays are processed from index 
 97:    l to index m only. The output information is:

 99:    ngrp - number of entries of the group
100:    ctr  - (w(l)+w(l+ngrp-1))/2
101:    ae   - average of wr(l),...,wr(l+ngrp-1)
102:    arsd - average of rsd(l),...,rsd(l+ngrp-1)
103: */
104: static PetscErrorCode EPSSubspaceFindGroup(PetscInt l,PetscInt m,PetscScalar *wr,PetscScalar *wi,PetscReal *rsd,PetscReal grptol,PetscInt *ngrp,PetscReal *ctr,PetscReal *ae,PetscReal *arsd)
105: {
106:   PetscInt  i;
107:   PetscReal rmod,rmod1;

110:   *ngrp = 0;
111:   *ctr = 0;
112:   rmod = SlepcAbsEigenvalue(wr[l],wi[l]);

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

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

142: /*
143:    EPSSubspaceResidualNorms - Computes the column norms of residual vectors
144:    OP*V(1:n,l:m) - V*T(1:m,l:m), where, on entry, OP*V has been computed and 
145:    stored in AV. ldt is the leading dimension of T. On exit, rsd(l) to
146:    rsd(m) contain the computed norms.
147: */
148: static PetscErrorCode EPSSubspaceResidualNorms(Vec *V,Vec *AV,PetscScalar *T,PetscInt l,PetscInt m,PetscInt ldt,Vec w,PetscReal *rsd)
149: {
151:   PetscInt       i,k;
152:   PetscScalar    t;

155:   for (i=l;i<m;i++) {
156:     if (i==m-1 || T[i+1+ldt*i]==0.0) k=i+1;
157:     else k=i+2;
158:     VecCopy(AV[i],w);
159:     SlepcVecMAXPBY(w,1.0,-1.0,k,T+ldt*i,V);
160:     VecDot(w,w,&t);
161:     rsd[i] = PetscRealPart(t);
162:   }
163:   for (i=l;i<m;i++) {
164:     if (i == m-1) {
165:       rsd[i] = PetscSqrtReal(rsd[i]);
166:     } else if (T[i+1+(ldt*i)]==0.0) {
167:       rsd[i] = PetscSqrtReal(rsd[i]);
168:     } else {
169:       rsd[i] = PetscSqrtReal((rsd[i]+rsd[i+1])/2.0);
170:       rsd[i+1] = rsd[i];
171:       i++;
172:     }
173:   }
174:   return(0);
175: }

179: PetscErrorCode EPSSolve_Subspace(EPS eps)
180: {
182:   EPS_SUBSPACE   *ctx = (EPS_SUBSPACE*)eps->data;
183:   PetscInt       i,k,ld,ngrp,nogrp,*itrsd,*itrsdold,
184:                  nxtsrr,idsrr,idort,nxtort,nv,ncv = eps->ncv,its;
185:   PetscScalar    *T,*U;
186:   PetscReal      arsd,oarsd,ctr,octr,ae,oae,*rsd,norm,tcond=1.0;
187:   PetscBool      breakdown;
188:   /* Parameters */
189:   PetscInt       init = 5;        /* Number of initial iterations */
190:   PetscReal      stpfac = 1.5,    /* Max num of iter before next SRR step */
191:                  alpha = 1.0,     /* Used to predict convergence of next residual */
192:                  beta = 1.1,      /* Used to predict convergence of next residual */
193:                  grptol = 1e-8,   /* Tolerance for EPSSubspaceFindGroup */
194:                  cnvtol = 1e-6;   /* Convergence criterion for cnv */
195:   PetscInt       orttol = 2;      /* Number of decimal digits whose loss
196:                                      can be tolerated in orthogonalization */

199:   its = 0;
200:   PetscMalloc(sizeof(PetscReal)*ncv,&rsd);
201:   PetscMalloc(sizeof(PetscInt)*ncv,&itrsd);
202:   PetscMalloc(sizeof(PetscInt)*ncv,&itrsdold);
203:   DSGetLeadingDimension(eps->ds,&ld);

205:   for (i=0;i<ncv;i++) {
206:     rsd[i] = 0.0;
207:     itrsd[i] = -1;
208:   }
209: 
210:   /* Complete the initial basis with random vectors and orthonormalize them */
211:   k = eps->nini;
212:   while (k<ncv) {
213:     SlepcVecSetRandom(eps->V[k],eps->rand);
214:     IPOrthogonalize(eps->ip,eps->nds,eps->defl,k,PETSC_NULL,eps->V,eps->V[k],PETSC_NULL,&norm,&breakdown);
215:     if (norm>0.0 && !breakdown) {
216:       VecScale(eps->V[k],1.0/norm);
217:       k++;
218:     }
219:   }

221:   while (eps->its<eps->max_it) {
222:     eps->its++;
223:     nv = PetscMin(eps->nconv+eps->mpd,ncv);
224:     DSSetDimensions(eps->ds,nv,PETSC_IGNORE,eps->nconv,0);
225: 
226:     /* Find group in previously computed eigenvalues */
227:     EPSSubspaceFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,rsd,grptol,&nogrp,&octr,&oae,&oarsd);

229:     /* AV(:,idx) = OP * V(:,idx) */
230:     for (i=eps->nconv;i<nv;i++) {
231:       STApply(eps->OP,eps->V[i],ctx->AV[i]);
232:     }

234:     /* T(:,idx) = V' * AV(:,idx) */
235:     DSGetArray(eps->ds,DS_MAT_A,&T);
236:     for (i=eps->nconv;i<nv;i++) {
237:       VecMDot(ctx->AV[i],nv,eps->V,T+i*ld);
238:     }
239:     DSRestoreArray(eps->ds,DS_MAT_A,&T);
240:     DSSetState(eps->ds,DS_STATE_RAW);

242:     /* Solve projected problem */
243:     DSSolve(eps->ds,eps->eigr,eps->eigi);
244:     DSSort(eps->ds,eps->eigr,eps->eigi,PETSC_NULL,PETSC_NULL,PETSC_NULL);
245: 
246:     /* Update vectors V(:,idx) = V * U(:,idx) */
247:     DSGetArray(eps->ds,DS_MAT_Q,&U);
248:     SlepcUpdateVectors(nv,ctx->AV,eps->nconv,nv,U,ld,PETSC_FALSE);
249:     SlepcUpdateVectors(nv,eps->V,eps->nconv,nv,U,ld,PETSC_FALSE);
250:     DSRestoreArray(eps->ds,DS_MAT_Q,&U);
251: 
252:     /* Convergence check */
253:     DSGetArray(eps->ds,DS_MAT_A,&T);
254:     EPSSubspaceResidualNorms(eps->V,ctx->AV,T,eps->nconv,nv,ld,eps->work[0],rsd);
255:     DSRestoreArray(eps->ds,DS_MAT_A,&T);

257:     for (i=eps->nconv;i<nv;i++) {
258:       itrsdold[i] = itrsd[i];
259:       itrsd[i] = its;
260:       eps->errest[i] = rsd[i];
261:     }
262: 
263:     for (;;) {
264:       /* Find group in currently computed eigenvalues */
265:       EPSSubspaceFindGroup(eps->nconv,nv,eps->eigr,eps->eigi,eps->errest,grptol,&ngrp,&ctr,&ae,&arsd);
266:       if (ngrp!=nogrp) break;
267:       if (ngrp==0) break;
268:       if (PetscAbsScalar(ae-oae)>ctr*cnvtol*(itrsd[eps->nconv]-itrsdold[eps->nconv])) break;
269:       if (arsd>ctr*eps->tol) break;
270:       eps->nconv = eps->nconv + ngrp;
271:       if (eps->nconv>=nv) break;
272:     }
273: 
274:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
275:     if (eps->nconv>=eps->nev) break;
276: 
277:     /* Compute nxtsrr (iteration of next projection step) */
278:     nxtsrr = PetscMin(eps->max_it,PetscMax((PetscInt)floor(stpfac*its),init));
279: 
280:     if (ngrp!=nogrp || ngrp==0 || arsd>=oarsd) {
281:       idsrr = nxtsrr - its;
282:     } else {
283:       idsrr = (PetscInt)floor(alpha+beta*(itrsdold[eps->nconv]-itrsd[eps->nconv])*log(arsd/eps->tol)/log(arsd/oarsd));
284:       idsrr = PetscMax(1,idsrr);
285:     }
286:     nxtsrr = PetscMin(nxtsrr,its+idsrr);

288:     /* Compute nxtort (iteration of next orthogonalization step) */
289:     DSCond(eps->ds,&tcond);
290:     idort = PetscMax(1,(PetscInt)floor(orttol/PetscMax(1,log10(tcond))));
291:     nxtort = PetscMin(its+idort,nxtsrr);

293:     /* V(:,idx) = AV(:,idx) */
294:     for (i=eps->nconv;i<nv;i++) {
295:       VecCopy(ctx->AV[i],eps->V[i]);
296:     }
297:     its++;

299:     /* Orthogonalization loop */
300:     do {
301:       while (its<nxtort) {
302: 
303:         /* AV(:,idx) = OP * V(:,idx) */
304:         for (i=eps->nconv;i<nv;i++) {
305:           STApply(eps->OP,eps->V[i],ctx->AV[i]);
306:         }
307: 
308:         /* V(:,idx) = AV(:,idx) with normalization */
309:         for (i=eps->nconv;i<nv;i++) {
310:           VecCopy(ctx->AV[i],eps->V[i]);
311:           VecNorm(eps->V[i],NORM_INFINITY,&norm);
312:           VecScale(eps->V[i],1/norm);
313:         }
314:         its++;
315:       }
316:       /* Orthonormalize vectors */
317:       for (i=eps->nconv;i<nv;i++) {
318:         IPOrthogonalize(eps->ip,eps->nds,eps->defl,i,PETSC_NULL,eps->V,eps->V[i],PETSC_NULL,&norm,&breakdown);
319:         if (breakdown) {
320:           SlepcVecSetRandom(eps->V[i],eps->rand);
321:           IPOrthogonalize(eps->ip,eps->nds,eps->defl,i,PETSC_NULL,eps->V,eps->V[i],PETSC_NULL,&norm,&breakdown);
322:         }
323:         VecScale(eps->V[i],1/norm);
324:       }
325:       nxtort = PetscMin(its+idort,nxtsrr);
326:     } while (its<nxtsrr);
327:   }

329:   PetscFree(rsd);
330:   PetscFree(itrsd);
331:   PetscFree(itrsdold);

333:   if (eps->nconv == eps->nev) eps->reason = EPS_CONVERGED_TOL;
334:   else eps->reason = EPS_DIVERGED_ITS;
335:   /* truncate Schur decomposition and change the state to raw so that
336:      PSVectors() computes eigenvectors from scratch */
337:   DSSetDimensions(eps->ds,eps->nconv,PETSC_IGNORE,0,0);
338:   DSSetState(eps->ds,DS_STATE_RAW);
339:   return(0);
340: }

344: PetscErrorCode EPSReset_Subspace(EPS eps)
345: {
347:   EPS_SUBSPACE   *ctx = (EPS_SUBSPACE*)eps->data;

350:   VecDestroyVecs(eps->ncv,&ctx->AV);
351:   EPSReset_Default(eps);
352:   return(0);
353: }

357: PetscErrorCode EPSDestroy_Subspace(EPS eps)
358: {

362:   PetscFree(eps->data);
363:   return(0);
364: }

366: EXTERN_C_BEGIN
369: PetscErrorCode EPSCreate_Subspace(EPS eps)
370: {

374:   PetscNewLog(eps,EPS_SUBSPACE,&eps->data);
375:   eps->ops->setup                = EPSSetUp_Subspace;
376:   eps->ops->destroy              = EPSDestroy_Subspace;
377:   eps->ops->reset                = EPSReset_Subspace;
378:   eps->ops->backtransform        = EPSBackTransform_Default;
379:   eps->ops->computevectors       = EPSComputeVectors_Schur;
380:   return(0);
381: }
382: EXTERN_C_END