Actual source code: krylovschur.c

  1: /*                       

  3:    SLEPc eigensolver: "krylovschur"

  5:    Method: Krylov-Schur

  7:    Algorithm:

  9:        Single-vector Krylov-Schur method for both symmetric and non-symmetric
 10:        problems.

 12:    References:

 14:        [1] "Krylov-Schur Methods in SLEPc", SLEPc Technical Report STR-7, 
 15:            available at http://www.grycap.upv.es/slepc.

 17:        [2] G.W. Stewart, "A Krylov-Schur Algorithm for Large Eigenproblems",
 18:            SIAM J. Matrix Analysis and App., 23(3), pp. 601-614, 2001. 

 20:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 21:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 22:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

 24:    This file is part of SLEPc.
 25:       
 26:    SLEPc is free software: you can redistribute it and/or modify it under  the
 27:    terms of version 3 of the GNU Lesser General Public License as published by
 28:    the Free Software Foundation.

 30:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 31:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 32:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 33:    more details.

 35:    You  should have received a copy of the GNU Lesser General  Public  License
 36:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 37:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 38: */

 40: #include <private/epsimpl.h>                /*I "slepceps.h" I*/
 41: #include <slepcblaslapack.h>

 43: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS);

 50: PetscErrorCode EPSSetUp_KrylovSchur(EPS eps)
 51: {
 53:   PetscBool      issinv;

 56:   /* spectrum slicing requires special treatment of default values */
 57:   if (eps->which==EPS_ALL) {
 58:     if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(((PetscObject)eps)->comm,1,"Must define a computational interval when using EPS_ALL");
 59:     if (!eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,1,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
 60:     if (!((PetscObject)(eps->OP))->type_name) { /* default to shift-and-invert */
 61:       STSetType(eps->OP,STSINVERT);
 62:     }
 63:     PetscTypeCompareAny((PetscObject)eps->OP,&issinv,STSINVERT,STCAYLEY,"");
 64:     if (!issinv) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
 65: #if defined(PETSC_USE_REAL_DOUBLE)
 66:     if (eps->tol==PETSC_DEFAULT) eps->tol = 1e-10;  /* use tighter tolerance */
 67: #endif
 68:     if (eps->intb >= PETSC_MAX_REAL) { /* right-open interval */
 69:       if (eps->inta <= PETSC_MIN_REAL) SETERRQ(((PetscObject)eps)->comm,1,"The defined computational interval should have at least one of their sides bounded");
 70:       STSetDefaultShift(eps->OP,eps->inta);
 71:     }
 72:     else { STSetDefaultShift(eps->OP,eps->intb); }

 74:     if (eps->nev==1) eps->nev = 40;  /* nev not set, use default value */
 75:     if (eps->nev<10) SETERRQ(((PetscObject)eps)->comm,1,"nev cannot be less than 10 in spectrum slicing runs");
 76:     eps->ops->backtransform = PETSC_NULL;
 77:   }

 79:   /* proceed with the general case */
 80:   if (eps->ncv) { /* ncv set */
 81:     if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
 82:   } else if (eps->mpd) { /* mpd set */
 83:     eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
 84:   } else { /* neither set: defaults depend on nev being small or large */
 85:     if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
 86:     else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
 87:   }
 88:   if (!eps->mpd) eps->mpd = eps->ncv;
 89:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must not be larger than nev+mpd");
 90:   if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 91:   if (!eps->which) { EPSDefaultSetWhich(eps); }
 92:   if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
 93:     SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");

 95:   if (!eps->extraction) {
 96:     EPSSetExtraction(eps,EPS_RITZ);
 97:   } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC) {
 98:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
 99:   }

101:   EPSAllocateSolution(eps);
102:   PetscFree(eps->T);
103:   if (!eps->ishermitian || eps->extraction==EPS_HARMONIC) {
104:     PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
105:   }
106:   EPSDefaultGetWork(eps,1);

108:   /* dispatch solve method */
109:   if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
110:   if (eps->ishermitian) {
111:     if (eps->which==EPS_ALL) eps->ops->solve = EPSSolve_KrylovSchur_Slice;
112:     else {
113:       switch (eps->extraction) {
114:         case EPS_RITZ:     eps->ops->solve = EPSSolve_KrylovSchur_Symm; break;
115:         case EPS_HARMONIC: eps->ops->solve = EPSSolve_KrylovSchur_Harmonic; break;
116:         default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
117:       }
118:     }
119:   } else {
120:     switch (eps->extraction) {
121:       case EPS_RITZ:     eps->ops->solve = EPSSolve_KrylovSchur_Default; break;
122:       case EPS_HARMONIC: eps->ops->solve = EPSSolve_KrylovSchur_Harmonic; break;
123:       default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
124:     }
125:   }
126:   return(0);
127: }

131: /*
132:    EPSProjectedKSNonsym - Solves the projected eigenproblem in the Krylov-Schur
133:    method (non-symmetric case).

135:    On input:
136:      l is the number of vectors kept in previous restart (0 means first restart)
137:      S is the projected matrix (leading dimension is lds)

139:    On output:
140:      S has (real) Schur form with diagonal blocks sorted appropriately
141:      Q contains the corresponding Schur vectors (order n, leading dimension n)
142: */
143: PetscErrorCode EPSProjectedKSNonsym(EPS eps,PetscInt l,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt n)
144: {
146:   PetscInt       i;

149:   if (l==0) {
150:     PetscMemzero(Q,n*n*sizeof(PetscScalar));
151:     for (i=0;i<n;i++)
152:       Q[i*(n+1)] = 1.0;
153:   } else {
154:     /* Reduce S to Hessenberg form, S <- Q S Q' */
155:     EPSDenseHessenberg(n,eps->nconv,S,lds,Q);
156:   }
157:   /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
158:   EPSDenseSchur(n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);
159:   /* Sort the remaining columns of the Schur form */
160:   EPSSortDenseSchur(eps,n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);
161:   return(0);
162: }

166: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
167: {
169:   PetscInt       i,k,l,lwork,nv;
170:   Vec            u=eps->work[0];
171:   PetscScalar    *S=eps->T,*Q,*work;
172:   PetscReal      beta;
173:   PetscBool      breakdown;

176:   PetscMemzero(S,eps->ncv*eps->ncv*sizeof(PetscScalar));
177:   PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&Q);
178:   lwork = 7*eps->ncv;
179:   PetscMalloc(lwork*sizeof(PetscScalar),&work);

181:   /* Get the starting Arnoldi vector */
182:   EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
183:   l = 0;
184: 
185:   /* Restart loop */
186:   while (eps->reason == EPS_CONVERGED_ITERATING) {
187:     eps->its++;

189:     /* Compute an nv-step Arnoldi factorization */
190:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
191:     EPSBasicArnoldi(eps,PETSC_FALSE,S,eps->ncv,eps->V,eps->nconv+l,&nv,u,&beta,&breakdown);
192:     VecScale(u,1.0/beta);

194:     /* Solve projected problem */
195:     EPSProjectedKSNonsym(eps,l,S,eps->ncv,Q,nv);

197:     /* Check convergence */
198:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->trackall,eps->nconv,nv-eps->nconv,S,eps->ncv,Q,eps->V,nv,beta,1.0,&k,work);
199:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
200:     if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
201: 
202:     /* Update l */
203:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
204:     else {
205:       l = (nv-k)/2;
206: #if !defined(PETSC_USE_COMPLEX)
207:       if (S[(k+l-1)*(eps->ncv+1)+1] != 0.0) {
208:         if (k+l<nv-1) l = l+1;
209:         else l = l-1;
210:       }
211: #endif
212:     }

214:     if (eps->reason == EPS_CONVERGED_ITERATING) {
215:       if (breakdown) {
216:         /* Start a new Arnoldi factorization */
217:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);
218:         EPSGetStartVector(eps,k,eps->V[k],&breakdown);
219:         if (breakdown) {
220:           eps->reason = EPS_DIVERGED_BREAKDOWN;
221:           PetscInfo(eps,"Unable to generate more start vectors\n");
222:         }
223:       } else {
224:         /* Prepare the Rayleigh quotient for restart */
225:         for (i=k;i<k+l;i++) {
226:           S[i*eps->ncv+k+l] = Q[(i+1)*nv-1]*beta;
227:         }
228:       }
229:     }
230:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
231:     SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,nv,PETSC_FALSE);

233:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
234:       VecCopy(u,eps->V[k+l]);
235:     }
236:     eps->nconv = k;

238:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
239:   }
240:   PetscFree(Q);
241:   PetscFree(work);
242:   return(0);
243: }

247: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
248: {

252:   PetscFree(eps->T);
253:   EPSReset_Default(eps);
254:   return(0);
255: }

260: PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
261: {
263:   eps->ops->setup          = EPSSetUp_KrylovSchur;
264:   eps->ops->reset          = EPSReset_KrylovSchur;
265:   eps->ops->backtransform  = EPSBackTransform_Default;
266:   eps->ops->computevectors = EPSComputeVectors_Schur;
267:   return(0);
268: }