Actual source code: krylov.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    Common subroutines for all Krylov-type solvers.

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

  8:    This file is part of SLEPc.

 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 = beta*v_m * e_m^T

 37:    where the columns of V are the Arnoldi vectors (which are B-orthonormal),
 38:    H is an upper Hessenberg matrix, e_m is the m-th vector of the canonical basis.
 39:    On exit, beta contains the B-norm of V[m] before normalization.
 40: */
 41: PetscErrorCode EPSBasicArnoldi(EPS eps,PetscBool trans,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,PetscReal *beta,PetscBool *breakdown)
 42: {
 44:   PetscInt       j,m = *M;
 45:   Vec            vj,vj1;

 48:   BVSetActiveColumns(eps->V,0,m);
 49:   for (j=k;j<m;j++) {
 50:     BVGetColumn(eps->V,j,&vj);
 51:     BVGetColumn(eps->V,j+1,&vj1);
 52:     if (trans) {
 53:       STApplyTranspose(eps->st,vj,vj1);
 54:     } else {
 55:       STApply(eps->st,vj,vj1);
 56:     }
 57:     BVRestoreColumn(eps->V,j,&vj);
 58:     BVRestoreColumn(eps->V,j+1,&vj1);
 59:     BVOrthogonalizeColumn(eps->V,j+1,H+ldh*j,beta,breakdown);
 60:     H[j+1+ldh*j] = *beta;
 61:     if (*breakdown) {
 62:       *M = j+1;
 63:       break;
 64:     } else {
 65:       BVScaleColumn(eps->V,j+1,1.0/(*beta));
 66:     }
 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:      getall - whether all residuals must be computed
 80:      kini  - initial value of k (the loop variable)
 81:      nits  - number of iterations of the loop
 82:      V     - set of basis vectors (used only if trueresidual is activated)
 83:      nv    - number of vectors to process (dimension of Q, columns of V)
 84:      beta  - norm of f (the residual vector of the Arnoldi/Lanczos factorization)
 85:      corrf - correction factor for residual estimates (only in harmonic KS)

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

100:   RGIsTrivial(eps->rg,&istrivial);
101:   if (eps->trueres) {
102:     BVGetVec(eps->V,&x);
103:     BVGetVec(eps->V,&y);
104:   }
105:   DSGetLeadingDimension(eps->ds,&ld);
106:   DSGetRefined(eps->ds,&refined);
107:   PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
108:   marker = -1;
109:   if (eps->trackall) getall = PETSC_TRUE;
110:   for (k=kini;k<kini+nits;k++) {
111:     /* eigenvalue */
112:     re = eps->eigr[k];
113:     im = eps->eigi[k];
114:     if (!istrivial || eps->trueres || isshift || eps->conv==EPS_CONV_NORM) {
115:       STBackTransform(eps->st,1,&re,&im);
116:     }
117:     if (!istrivial) {
118:       RGCheckInside(eps->rg,1,&re,&im,&inside);
119:       if (marker==-1 && inside<=0) marker = k;
120:       if (!(eps->trueres || isshift || eps->conv==EPS_CONV_NORM)) {  /* make sure eps->converged below uses the right value */
121:         re = eps->eigr[k];
122:         im = eps->eigi[k];
123:       }
124:     }
125:     newk = k;
126:     DSVectors(eps->ds,DS_MAT_X,&newk,&resnorm);
127:     if (eps->trueres) {
128:       DSGetArray(eps->ds,DS_MAT_X,&X);
129:       Zr = X+k*ld;
130:       if (newk==k+1) Zi = X+newk*ld;
131:       else Zi = NULL;
132:       EPSComputeRitzVector(eps,Zr,Zi,eps->V,x,y);
133:       DSRestoreArray(eps->ds,DS_MAT_X,&X);
134:       EPSComputeResidualNorm_Private(eps,re,im,x,y,&resnorm);
135:     }
136:     else if (!refined) resnorm *= beta*corrf;
137:     /* error estimate */
138:     (*eps->converged)(eps,re,im,resnorm,&eps->errest[k],eps->convergedctx);
139:     if (marker==-1 && eps->errest[k] >= eps->tol) marker = k;
140:     if (newk==k+1) {
141:       eps->errest[k+1] = eps->errest[k];
142:       k++;
143:     }
144:     if (marker!=-1 && !getall) break;
145:   }
146:   if (marker!=-1) k = marker;
147:   *kout = k;
148:   if (eps->trueres) {
149:     VecDestroy(&x);
150:     VecDestroy(&y);
151:   }
152:   return(0);
153: }

157: /*
158:    EPSFullLanczos - Computes an m-step Lanczos factorization with full
159:    reorthogonalization.  At each Lanczos step, the corresponding Lanczos
160:    vector is orthogonalized with respect to all previous Lanczos vectors.
161:    This is equivalent to computing an m-step Arnoldi factorization and
162:    exploting symmetry of the operator.

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

167:                     OP * V - V * T = beta_m*v_m * e_m^T

169:    where the columns of V are the Lanczos vectors (which are B-orthonormal),
170:    T is a real symmetric tridiagonal matrix, and e_m is the m-th vector of
171:    the canonical basis. The tridiagonal is stored as two arrays: alpha
172:    contains the diagonal elements, beta the off-diagonal. On exit, the last
173:    element of beta contains the B-norm of V[m] before normalization.
174: */
175: PetscErrorCode EPSFullLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown)
176: {
178:   PetscInt       j,m = *M;
179:   Vec            vj,vj1;
180:   PetscScalar    *hwork,lhwork[100];

183:   if (m > 100) {
184:     PetscMalloc1(m,&hwork);
185:   } else hwork = lhwork;

187:   BVSetActiveColumns(eps->V,0,m);
188:   for (j=k;j<m;j++) {
189:     BVGetColumn(eps->V,j,&vj);
190:     BVGetColumn(eps->V,j+1,&vj1);
191:     STApply(eps->st,vj,vj1);
192:     BVRestoreColumn(eps->V,j,&vj);
193:     BVRestoreColumn(eps->V,j+1,&vj1);
194:     BVOrthogonalizeColumn(eps->V,j+1,hwork,beta+j,breakdown);
195:     alpha[j] = PetscRealPart(hwork[j]);
196:     if (*breakdown) {
197:       *M = j+1;
198:       break;
199:     } else {
200:       BVScaleColumn(eps->V,j+1,1.0/beta[j]);
201:     }
202:   }
203:   if (m > 100) {
204:     PetscFree(hwork);
205:   }
206:   return(0);
207: }

211: PetscErrorCode EPSPseudoLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal *cos,Vec w)
212: {
214:   PetscInt       j,m = *M;
215:   Vec            vj,vj1;
216:   PetscScalar    *hwork,lhwork[100];
217:   PetscReal      norm,norm1,norm2,t;

220:   if (cos) *cos = 1.0;
221:   if (m > 100) {
222:     PetscMalloc1(m,&hwork);
223:   } else hwork = lhwork;

225:   BVSetActiveColumns(eps->V,0,m);
226:   for (j=k;j<m;j++) {
227:     BVGetColumn(eps->V,j,&vj);
228:     BVGetColumn(eps->V,j+1,&vj1);
229:     STApply(eps->st,vj,vj1);
230:     BVRestoreColumn(eps->V,j,&vj);
231:     BVRestoreColumn(eps->V,j+1,&vj1);
232:     BVOrthogonalizeColumn(eps->V,j+1,hwork,&norm,breakdown);
233:     alpha[j] = PetscRealPart(hwork[j]);
234:     beta[j] = PetscAbsReal(norm);
235:     omega[j+1] = (norm<0.0)? -1.0: 1.0;
236:     BVScaleColumn(eps->V,j+1,1.0/norm);
237:     /* */
238:     BVGetColumn(eps->V,j+1,&vj1);
239:     VecNorm(vj1,NORM_2,&norm1);
240:     BVApplyMatrix(eps->V,vj1,w);
241:     BVRestoreColumn(eps->V,j+1,&vj1);
242:     VecNorm(w,NORM_2,&norm2);
243:     t = 1.0/(norm1*norm2);
244:     if (cos && *cos>t) *cos = t;
245:   }
246:   if (m > 100) {
247:     PetscFree(hwork);
248:   }
249:   return(0);
250: }