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:    Last update: Feb 2009

 22:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 24:    Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain

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

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

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

 42:  #include private/epsimpl.h
 43:  #include slepcblaslapack.h

 45: PetscErrorCode EPSSolve_KRYLOVSCHUR_DEFAULT(EPS eps);

 51: PetscErrorCode EPSSetUp_KRYLOVSCHUR(EPS eps)
 52: {
 54:   PetscInt       N;

 57:   VecGetSize(eps->vec_initial,&N);
 58:   if (eps->ncv) { /* ncv set */
 59:     if (eps->ncv<eps->nev+1) SETERRQ(1,"The value of ncv must be at least nev+1");
 60:   }
 61:   else if (eps->mpd) { /* mpd set */
 62:     eps->ncv = PetscMin(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(N,PetscMax(2*eps->nev,eps->nev+15));
 66:     else { eps->mpd = 500; eps->ncv = PetscMin(N,eps->nev+eps->mpd); }
 67:   }
 68:   if (!eps->mpd) eps->mpd = eps->ncv;
 69:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(1,"The value of ncv must not be larger than nev+mpd");
 70:   if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);
 71:   if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
 72:     SETERRQ(1,"Wrong value of eps->which");

 74:   if (!eps->extraction) {
 75:     EPSSetExtraction(eps,EPS_RITZ);
 76:   } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC) {
 77:     SETERRQ(PETSC_ERR_SUP,"Unsupported extraction type");
 78:   }

 80:   EPSAllocateSolution(eps);
 81:   PetscFree(eps->T);
 82:   if (!eps->ishermitian || eps->extraction==EPS_HARMONIC) {
 83:     PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
 84:   }
 85:   EPSDefaultGetWork(eps,2);
 86:   return(0);
 87: }

 91: /*
 92:    EPSProjectedKSNonsym - Solves the projected eigenproblem in the Krylov-Schur
 93:    method (non-symmetric case).

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

 99:    On output:
100:      S has (real) Schur form with diagonal blocks sorted appropriately
101:      Q contains the corresponding Schur vectors (order n, leading dimension n)
102: */
103: PetscErrorCode EPSProjectedKSNonsym(EPS eps,PetscInt l,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt n)
104: {
106:   PetscInt       i;

109:   if (l==0) {
110:     PetscMemzero(Q,n*n*sizeof(PetscScalar));
111:     for (i=0;i<n;i++)
112:       Q[i*(n+1)] = 1.0;
113:   } else {
114:     /* Reduce S to Hessenberg form, S <- Q S Q' */
115:     EPSDenseHessenberg(n,eps->nconv,S,lds,Q);
116:   }
117:   /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
118:   EPSDenseSchur(n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);
119:   /* Sort the remaining columns of the Schur form */
120:   EPSSortDenseSchur(n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi,eps->which);
121:   return(0);
122: }

126: PetscErrorCode EPSSolve_KRYLOVSCHUR(EPS eps)
127: {

131:   if (eps->ishermitian) {
132:     switch (eps->extraction) {
133:       case EPS_RITZ:
134:         EPSSolve_KRYLOVSCHUR_SYMM(eps);
135:         break;
136:       case EPS_HARMONIC:
137:         EPSSolve_KRYLOVSCHUR_HARMONIC(eps);
138:         break;
139:       default: SETERRQ(PETSC_ERR_SUP,"Unsupported extraction type");
140:     }
141:   } else {
142:     switch (eps->extraction) {
143:       case EPS_RITZ:
144:         EPSSolve_KRYLOVSCHUR_DEFAULT(eps);
145:         break;
146:       case EPS_HARMONIC:
147:         EPSSolve_KRYLOVSCHUR_HARMONIC(eps);
148:         break;
149:       default: SETERRQ(PETSC_ERR_SUP,"Unsupported extraction type");
150:     }
151:   }
152:   return(0);
153: }

157: PetscErrorCode EPSSolve_KRYLOVSCHUR_DEFAULT(EPS eps)
158: {
160:   PetscInt       i,k,l,lwork,nv;
161:   Vec            u=eps->work[1];
162:   PetscScalar    *S=eps->T,*Q,*work;
163:   PetscReal      beta;
164:   PetscTruth     breakdown;

167:   PetscMemzero(S,eps->ncv*eps->ncv*sizeof(PetscScalar));
168:   PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&Q);
169:   lwork = (eps->ncv+4)*eps->ncv;
170:   PetscMalloc(lwork*sizeof(PetscScalar),&work);

172:   /* Get the starting Arnoldi vector */
173:   EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
174:   l = 0;
175: 
176:   /* Restart loop */
177:   while (eps->reason == EPS_CONVERGED_ITERATING) {
178:     eps->its++;

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

185:     /* Solve projected problem and compute residual norm estimates */
186:     EPSProjectedKSNonsym(eps,l,S,eps->ncv,Q,nv);
187:     ArnoldiResiduals(S,eps->ncv,Q,beta,eps->nconv,nv,eps->eigr,eps->eigi,eps->errest,work);

189:     /* Check convergence */
190:     k = eps->nconv;
191:     while (k<nv && eps->errest[k]<eps->tol) k++;
192:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
193:     if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
194: 
195:     /* Update l */
196:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
197:     else {
198:       l = (nv-k)/2;
199: #if !defined(PETSC_USE_COMPLEX)
200:       if (S[(k+l-1)*(eps->ncv+1)+1] != 0.0) {
201:         if (k+l<nv-1) l = l+1;
202:         else l = l-1;
203:       }
204: #endif
205:     }
206: 
207:     if (eps->reason == EPS_CONVERGED_ITERATING) {
208:       if (breakdown) {
209:         /* Start a new Arnoldi factorization */
210:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%i norm=%g)\n",eps->its,beta);
211:         EPSGetStartVector(eps,k,eps->V[k],&breakdown);
212:         if (breakdown) {
213:           eps->reason = EPS_DIVERGED_BREAKDOWN;
214:           PetscInfo(eps,"Unable to generate more start vectors\n");
215:         }
216:       } else {
217:         /* Prepare the Rayleigh quotient for restart */
218:         for (i=k;i<k+l;i++) {
219:           S[i*eps->ncv+k+l] = Q[(i+1)*nv-1]*beta;
220:         }
221:       }
222:     }
223:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
224:     SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,nv,PETSC_FALSE);

226:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
227:       VecCopy(u,eps->V[k+l]);
228:     }
229:     eps->nconv = k;

231:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
232: 
233:   }

235:   PetscFree(Q);
236:   PetscFree(work);
237:   return(0);
238: }

243: PetscErrorCode EPSCreate_KRYLOVSCHUR(EPS eps)
244: {
246:   eps->data                      = PETSC_NULL;
247:   eps->ops->solve                = EPSSolve_KRYLOVSCHUR;
248:   eps->ops->solvets              = PETSC_NULL;
249:   eps->ops->setup                = EPSSetUp_KRYLOVSCHUR;
250:   eps->ops->setfromoptions       = PETSC_NULL;
251:   eps->ops->destroy              = EPSDestroy_Default;
252:   eps->ops->view                 = PETSC_NULL;
253:   eps->ops->backtransform        = EPSBackTransform_Default;
254:   eps->ops->computevectors       = EPSComputeVectors_Schur;
255:   return(0);
256: }