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-2012, Universitat 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 <slepc-private/epsimpl.h>
 25: #include <slepc-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,PetscBool trans,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscBool *breakdown)
 44: {
 46:   PetscInt       j,m = *M;
 47:   PetscReal      norm;

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

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

 75:    Input Parameters:
 76:      eps   - the eigensolver; some error estimates are updated in eps->errest 
 77:      getall - whether all residuals must be computed
 78:      kini  - initial value of k (the loop variable)
 79:      nits  - number of iterations of the loop
 80:      V     - set of basis vectors (used only if trueresidual is activated)
 81:      nv    - number of vectors to process (dimension of Q, columns of V)
 82:      beta  - norm of f (the residual vector of the Arnoldi/Lanczos factorization)
 83:      corrf - correction factor for residual estimates (only in harmonic KS)

 85:    Output Parameters:
 86:      kout  - the first index where the convergence test failed
 87: */
 88: PetscErrorCode EPSKrylovConvergence(EPS eps,PetscBool getall,PetscInt kini,PetscInt nits,Vec *V,PetscInt nv,PetscReal beta,PetscReal corrf,PetscInt *kout)
 89: {
 91:   PetscInt       k,newk,marker,ld;
 92:   PetscScalar    re,im,*Zr,*Zi,*X;
 93:   PetscReal      resnorm;
 94:   PetscBool      isshift,refined;
 95:   Vec            x,y;

 98:   if (eps->trueres) {
 99:     VecDuplicate(eps->t,&x);
100:     VecDuplicate(eps->t,&y);
101:   }
102:   DSGetLeadingDimension(eps->ds,&ld);
103:   DSGetRefined(eps->ds,&refined);
104:   PetscObjectTypeCompare((PetscObject)eps->OP,STSHIFT,&isshift);
105:   marker = -1;
106:   if (eps->trackall) getall = PETSC_TRUE;
107:   for (k=kini;k<kini+nits;k++) {
108:     /* eigenvalue */
109:     re = eps->eigr[k];
110:     im = eps->eigi[k];
111:     if (eps->trueres || isshift) {
112:       STBackTransform(eps->OP,1,&re,&im);
113:     }
114:     newk = k;
115:     DSVectors(eps->ds,DS_MAT_X,&newk,&resnorm);
116:     if (eps->trueres) {
117:       DSGetArray(eps->ds,DS_MAT_X,&X);
118:       Zr = X+k*ld;
119:       if (newk==k+1) Zi = X+newk*ld;
120:       else Zi = PETSC_NULL;
121:       EPSComputeRitzVector(eps,Zr,Zi,V,nv,x,y);
122:       DSRestoreArray(eps->ds,DS_MAT_X,&X);
123:       EPSComputeResidualNorm_Private(eps,re,im,x,y,&resnorm);
124:     }
125:     else if (!refined) resnorm *= beta*corrf;
126:     /* error estimate */
127:     (*eps->conv_func)(eps,re,im,resnorm,&eps->errest[k],eps->conv_ctx);
128:     if (marker==-1 && eps->errest[k] >= eps->tol) marker = k;
129:     if (newk==k+1) { eps->errest[k+1] = eps->errest[k]; k++; }
130:     if (marker!=-1 && !getall) break;
131:   }
132:   if (marker!=-1) k = marker;
133:   *kout = k;
134:   if (eps->trueres) {
135:     VecDestroy(&x);
136:     VecDestroy(&y);
137:   }
138:   return(0);
139: }

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

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

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

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

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

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

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