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-2013, 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 = 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) {
52: STApplyTranspose(eps->st,V[j],V[j+1]);
53: } else {
54: STApply(eps->st,V[j],V[j+1]);
55: }
56: IPOrthogonalize(eps->ip,eps->nds,eps->defl,j+1,NULL,V,V[j+1],H+ldh*j,&norm,breakdown);
57: H[j+1+ldh*j] = norm;
58: if (*breakdown) {
59: *M = j+1;
60: *beta = norm;
61: return(0);
62: } else {
63: VecScale(V[j+1],1/norm);
64: }
65: }
66: if (trans) {
67: STApplyTranspose(eps->st,V[m-1],f);
68: } else {
69: STApply(eps->st,V[m-1],f);
70: }
71: IPOrthogonalize(eps->ip,eps->nds,eps->defl,m,NULL,V,f,H+ldh*(m-1),beta,NULL);
72: return(0);
73: }
77: /*
78: EPSKrylovConvergence - Implements the loop that checks for convergence
79: in Krylov methods.
81: Input Parameters:
82: eps - the eigensolver; some error estimates are updated in eps->errest
83: getall - whether all residuals must be computed
84: kini - initial value of k (the loop variable)
85: nits - number of iterations of the loop
86: V - set of basis vectors (used only if trueresidual is activated)
87: nv - number of vectors to process (dimension of Q, columns of V)
88: beta - norm of f (the residual vector of the Arnoldi/Lanczos factorization)
89: corrf - correction factor for residual estimates (only in harmonic KS)
91: Output Parameters:
92: kout - the first index where the convergence test failed
93: */
94: PetscErrorCode EPSKrylovConvergence(EPS eps,PetscBool getall,PetscInt kini,PetscInt nits,Vec *V,PetscInt nv,PetscReal beta,PetscReal corrf,PetscInt *kout)
95: {
97: PetscInt k,newk,marker,ld;
98: PetscScalar re,im,*Zr,*Zi,*X;
99: PetscReal resnorm;
100: PetscBool isshift,refined;
101: Vec x,y;
104: if (eps->trueres) {
105: VecDuplicate(eps->t,&x);
106: VecDuplicate(eps->t,&y);
107: }
108: DSGetLeadingDimension(eps->ds,&ld);
109: DSGetRefined(eps->ds,&refined);
110: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
111: marker = -1;
112: if (eps->trackall) getall = PETSC_TRUE;
113: for (k=kini;k<kini+nits;k++) {
114: /* eigenvalue */
115: re = eps->eigr[k];
116: im = eps->eigi[k];
117: if (eps->trueres || isshift) {
118: STBackTransform(eps->st,1,&re,&im);
119: }
120: newk = k;
121: DSVectors(eps->ds,DS_MAT_X,&newk,&resnorm);
122: if (eps->trueres) {
123: DSGetArray(eps->ds,DS_MAT_X,&X);
124: Zr = X+k*ld;
125: if (newk==k+1) Zi = X+newk*ld;
126: else Zi = NULL;
127: EPSComputeRitzVector(eps,Zr,Zi,V,nv,x,y);
128: DSRestoreArray(eps->ds,DS_MAT_X,&X);
129: EPSComputeResidualNorm_Private(eps,re,im,x,y,&resnorm);
130: }
131: else if (!refined) resnorm *= beta*corrf;
132: /* error estimate */
133: (*eps->converged)(eps,re,im,resnorm,&eps->errest[k],eps->convergedctx);
134: if (marker==-1 && eps->errest[k] >= eps->tol) marker = k;
135: if (newk==k+1) {
136: eps->errest[k+1] = eps->errest[k];
137: k++;
138: }
139: if (marker!=-1 && !getall) break;
140: }
141: if (marker!=-1) k = marker;
142: *kout = k;
143: if (eps->trueres) {
144: VecDestroy(&x);
145: VecDestroy(&y);
146: }
147: return(0);
148: }
152: /*
153: EPSFullLanczos - Computes an m-step Lanczos factorization with full
154: reorthogonalization. At each Lanczos step, the corresponding Lanczos
155: vector is orthogonalized with respect to all previous Lanczos vectors.
156: This is equivalent to computing an m-step Arnoldi factorization and
157: exploting symmetry of the operator.
159: The first k columns are assumed to be locked and therefore they are
160: not modified. On exit, the following relation is satisfied:
162: OP * V - V * T = f * e_m^T
164: where the columns of V are the Lanczos vectors (which are B-orthonormal),
165: T is a real symmetric tridiagonal matrix, f is the residual vector and e_m
166: is the m-th vector of the canonical basis. The tridiagonal is stored as
167: two arrays: alpha contains the diagonal elements, beta the off-diagonal.
168: The vector f is B-orthogonal to the columns of V. On exit, the last element
169: of beta contains the B-norm of f and the next Lanczos vector can be
170: computed as v_{m+1} = f / beta(end).
172: */
173: PetscErrorCode EPSFullLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscBool *breakdown)
174: {
176: PetscInt j,m = *M;
177: PetscReal norm;
178: PetscScalar *hwork,lhwork[100];
181: if (m > 100) {
182: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&hwork);
183: } else hwork = lhwork;
185: for (j=k;j<m-1;j++) {
186: STApply(eps->st,V[j],V[j+1]);
187: IPOrthogonalize(eps->ip,eps->nds,eps->defl,j+1,NULL,V,V[j+1],hwork,&norm,breakdown);
188: alpha[j] = PetscRealPart(hwork[j]);
189: beta[j] = norm;
190: if (*breakdown) {
191: *M = j+1;
192: if (m > 100) {
193: PetscFree(hwork);
194: }
195: return(0);
196: } else {
197: VecScale(V[j+1],1.0/norm);
198: }
199: }
200: STApply(eps->st,V[m-1],f);
201: IPOrthogonalize(eps->ip,eps->nds,eps->defl,m,NULL,V,f,hwork,&norm,NULL);
202: alpha[m-1] = PetscRealPart(hwork[m-1]);
203: beta[m-1] = norm;
205: if (m > 100) {
206: PetscFree(hwork);
207: }
208: return(0);
209: }