Actual source code: krylov.c

  1: /*                       
  2:    Common subroutines for all Krylov-type solvers.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.
  9:       
 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

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

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

 24:  #include private/epsimpl.h
 25:  #include private/slepcimpl.h
 26:  #include slepcblaslapack.h

 30: /*
 31:    EPSBasicArnoldi - Computes an m-step Arnoldi factorization. The first k
 32:    columns are assumed to be locked and therefore they are not modified. On
 33:    exit, the following relation is satisfied:

 35:                     OP * V - V * H = f * e_m^T

 37:    where the columns of V are the Arnoldi vectors (which are B-orthonormal),
 38:    H is an upper Hessenberg matrix, f is the residual vector and e_m is
 39:    the m-th vector of the canonical basis. The vector f is B-orthogonal to
 40:    the columns of V. On exit, beta contains the B-norm of f and the next 
 41:    Arnoldi vector can be computed as v_{m+1} = f / beta. 
 42: */
 43: PetscErrorCode EPSBasicArnoldi(EPS eps,PetscTruth trans,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscTruth *breakdown)
 44: {
 46:   PetscInt       j,m = *M;
 47:   PetscReal      norm;

 50: 
 51:   for (j=k;j<m-1;j++) {
 52:     if (trans) { STApplyTranspose(eps->OP,V[j],V[j+1]); }
 53:     else { STApply(eps->OP,V[j],V[j+1]); }
 54:     IPOrthogonalize(eps->ip,eps->nds,eps->DS,j+1,PETSC_NULL,V,V[j+1],H+ldh*j,&norm,breakdown);
 55:     H[j+1+ldh*j] = norm;
 56:     if (*breakdown) {
 57:       *M = j+1;
 58:       *beta = norm;
 59:       return(0);
 60:     } else {
 61:       VecScale(V[j+1],1/norm);
 62:     }
 63:   }
 64:   if (trans) { STApplyTranspose(eps->OP,V[m-1],f); }
 65:   else { STApply(eps->OP,V[m-1],f); }
 66:   IPOrthogonalize(eps->ip,eps->nds,eps->DS,m,PETSC_NULL,V,f,H+ldh*(m-1),beta,PETSC_NULL);
 67: 
 68:   return(0);
 69: }

 73: /*
 74:    EPSKrylovConvergence - Implements the loop that checks for convergence
 75:    in Krylov methods.

 77:    Input Parameters:
 78:      eps   - the eigensolver; some error estimates are updated in eps->errest 
 79:      issym - whether the projected problem is symmetric or not
 80:      kini  - initial value of k (the loop variable)
 81:      nits  - number of iterations of the loop
 82:      S     - Schur form of projected matrix (not referenced if issym)
 83:      lds   - leading dimension of S
 84:      Q     - Schur vectors of projected matrix (eigenvectors if issym)
 85:      V     - set of basis vectors (used only if trueresidual is activated)
 86:      nv    - number of vectors to process (dimension of Q, columns of V)
 87:      beta  - norm of f (the residual vector of the Arnoldi/Lanczos factorization)
 88:      corrf - correction factor for residual estimates (only in harmonic KS)

 90:    Output Parameters:
 91:      kout  - the first index where the convergence test failed

 93:    Workspace:
 94:      work is workspace to store 5*nv scalars, nv booleans and nv reals (only if !issym)
 95: */
 96: PetscErrorCode EPSKrylovConvergence(EPS eps,PetscTruth issym,PetscInt kini,PetscInt nits,PetscScalar *S,PetscInt lds,PetscScalar *Q,Vec *V,PetscInt nv,PetscReal beta,PetscReal corrf,PetscInt *kout,PetscScalar *work)
 97: {
 99:   PetscInt       k,marker;
100:   PetscScalar    re,im,*Z = work,*work2 = work;
101:   PetscReal      resnorm;
102:   PetscTruth     iscomplex,isshift;

105:   if (!issym) { Z = work; work2 = work+2*nv; }
106:   PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isshift);
107:   marker = -1;
108:   for (k=kini;k<kini+nits;k++) {
109:     /* eigenvalue */
110:     re = eps->eigr[k];
111:     im = eps->eigi[k];
112:     if (eps->trueres || isshift) {
113:       STBackTransform(eps->OP,1,&re,&im);
114:     }
115:     iscomplex = PETSC_FALSE;
116:     if (!issym && k<nv-1 && S[k+1+k*lds] != 0.0) iscomplex = PETSC_TRUE;
117:     /* residual norm */
118:     if (issym) {
119:       resnorm = beta*PetscAbsScalar(Q[(k-kini+1)*nv-1]);
120:     } else {
121:       DenseSelectedEvec(S,lds,Q,Z,k,iscomplex,nv,work2);
122:       if (iscomplex) resnorm = beta*SlepcAbsEigenvalue(Z[nv-1],Z[2*nv-1]);
123:       else resnorm = beta*PetscAbsScalar(Z[nv-1]);
124:     }
125:     if (eps->trueres) {
126:       if (issym) Z = Q+(k-kini)*nv;
127:       EPSComputeTrueResidual(eps,re,im,Z,V,nv,&resnorm);
128:     }
129:     else resnorm *= corrf;
130:     /* error estimate */
131:     (*eps->conv_func)(eps,re,im,resnorm,&eps->errest[k],eps->conv_ctx);
132:     if (marker==-1 && eps->errest[k] >= eps->tol) marker = k;
133:     if (iscomplex) { eps->errest[k+1] = eps->errest[k]; k++; }
134:     if (marker!=-1 && !eps->trackall) break;
135:   }
136:   if (marker!=-1) k = marker;
137:   *kout = k;

139:   return(0);
140: }

144: /*
145:    EPSFullLanczos - Computes an m-step Lanczos factorization with full
146:    reorthogonalization.  At each Lanczos step, the corresponding Lanczos
147:    vector is orthogonalized with respect to all previous Lanczos vectors.
148:    This is equivalent to computing an m-step Arnoldi factorization and
149:    exploting symmetry of the operator.

151:    The first k columns are assumed to be locked and therefore they are 
152:    not modified. On exit, the following relation is satisfied:

154:                     OP * V - V * T = f * e_m^T

156:    where the columns of V are the Lanczos vectors (which are B-orthonormal),
157:    T is a real symmetric tridiagonal matrix, f is the residual vector and e_m
158:    is the m-th vector of the canonical basis. The tridiagonal is stored as
159:    two arrays: alpha contains the diagonal elements, beta the off-diagonal.
160:    The vector f is B-orthogonal to the columns of V. On exit, the last element
161:    of beta contains the B-norm of f and the next Lanczos vector can be 
162:    computed as v_{m+1} = f / beta(end). 

164: */
165: PetscErrorCode EPSFullLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscTruth *breakdown)
166: {
168:   PetscInt       j,m = *M;
169:   PetscReal      norm;
170:   PetscScalar    *hwork,lhwork[100];

173:   if (m > 100) {
174:     PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&hwork);
175:   } else {
176:     hwork = lhwork;
177:   }

179:   for (j=k;j<m-1;j++) {
180:     STApply(eps->OP,V[j],V[j+1]);
181:     IPOrthogonalize(eps->ip,eps->nds,eps->DS,j+1,PETSC_NULL,V,V[j+1],hwork,&norm,breakdown);
182:     alpha[j-k] = PetscRealPart(hwork[j]);
183:     beta[j-k] = norm;
184:     if (*breakdown) {
185:       *M = j+1;
186:       if (m > 100) {
187:         PetscFree(hwork);
188:       }
189:       return(0);
190:     } else {
191:       VecScale(V[j+1],1.0/norm);
192:     }
193:   }
194:   STApply(eps->OP,V[m-1],f);
195:   IPOrthogonalize(eps->ip,eps->nds,eps->DS,m,PETSC_NULL,V,f,hwork,&norm,PETSC_NULL);
196:   alpha[m-1-k] = PetscRealPart(hwork[m-1]);
197:   beta[m-1-k] = norm;
198: 
199:   if (m > 100) {
200:     PetscFree(hwork);
201:   }
202:   return(0);
203: }