Actual source code: epskrylov.c
slepc-3.6.1 2015-09-03
1: /*
2: Common subroutines for all Krylov-type solvers.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, 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: EPSDelayedArnoldi - This function is equivalent to EPSBasicArnoldi but
75: performs the computation in a different way. The main idea is that
76: reorthogonalization is delayed to the next Arnoldi step. This version is
77: more scalable but in some cases convergence may stagnate.
78: */
79: PetscErrorCode EPSDelayedArnoldi(EPS eps,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,PetscReal *beta,PetscBool *breakdown)
80: {
82: PetscInt i,j,m=*M;
83: Vec u,t;
84: PetscScalar shh[100],*lhh,dot,dot2;
85: PetscReal norm1=0.0,norm2;
86: Vec vj,vj1,vj2;
89: if (m<=100) lhh = shh;
90: else {
91: PetscMalloc1(m,&lhh);
92: }
93: BVCreateVec(eps->V,&u);
94: BVCreateVec(eps->V,&t);
96: BVSetActiveColumns(eps->V,0,m);
97: for (j=k;j<m;j++) {
98: BVGetColumn(eps->V,j,&vj);
99: BVGetColumn(eps->V,j+1,&vj1);
100: STApply(eps->st,vj,vj1);
101: BVRestoreColumn(eps->V,j,&vj);
102: BVRestoreColumn(eps->V,j+1,&vj1);
104: BVDotColumnBegin(eps->V,j+1,H+ldh*j);
105: if (j>k) {
106: BVDotColumnBegin(eps->V,j,lhh);
107: BVGetColumn(eps->V,j,&vj);
108: VecDotBegin(vj,vj,&dot);
109: }
110: if (j>k+1) {
111: BVNormVecBegin(eps->V,u,NORM_2,&norm2);
112: BVGetColumn(eps->V,j-2,&vj2);
113: VecDotBegin(u,vj2,&dot2);
114: }
116: BVDotColumnEnd(eps->V,j+1,H+ldh*j);
117: if (j>k) {
118: BVDotColumnEnd(eps->V,j,lhh);
119: VecDotEnd(vj,vj,&dot);
120: BVRestoreColumn(eps->V,j,&vj);
121: }
122: if (j>k+1) {
123: BVNormVecEnd(eps->V,u,NORM_2,&norm2);
124: VecDotEnd(u,vj2,&dot2);
125: BVRestoreColumn(eps->V,j-2,&vj2);
126: }
128: if (j>k) {
129: norm1 = PetscSqrtReal(PetscRealPart(dot));
130: for (i=0;i<j;i++)
131: H[ldh*j+i] = H[ldh*j+i]/norm1;
132: H[ldh*j+j] = H[ldh*j+j]/dot;
134: BVCopyVec(eps->V,j,t);
135: BVScaleColumn(eps->V,j,1.0/norm1);
136: BVScaleColumn(eps->V,j+1,1.0/norm1);
137: }
139: BVMultColumn(eps->V,-1.0,1.0,j+1,H+ldh*j);
141: if (j>k) {
142: BVSetActiveColumns(eps->V,0,j);
143: BVMultVec(eps->V,-1.0,1.0,t,lhh);
144: BVSetActiveColumns(eps->V,0,m);
145: for (i=0;i<j;i++)
146: H[ldh*(j-1)+i] += lhh[i];
147: }
149: if (j>k+1) {
150: BVGetColumn(eps->V,j-1,&vj1);
151: VecCopy(u,vj1);
152: BVRestoreColumn(eps->V,j-1,&vj1);
153: BVScaleColumn(eps->V,j-1,1.0/norm2);
154: H[ldh*(j-2)+j-1] = norm2;
155: }
157: if (j<m-1) {
158: VecCopy(t,u);
159: }
160: }
162: BVNormVec(eps->V,t,NORM_2,&norm2);
163: VecScale(t,1.0/norm2);
164: BVGetColumn(eps->V,m-1,&vj1);
165: VecCopy(t,vj1);
166: BVRestoreColumn(eps->V,m-1,&vj1);
167: H[ldh*(m-2)+m-1] = norm2;
169: BVDotColumn(eps->V,m,lhh);
171: BVMultColumn(eps->V,-1.0,1.0,m,lhh);
172: for (i=0;i<m;i++)
173: H[ldh*(m-1)+i] += lhh[i];
175: BVNormColumn(eps->V,m,NORM_2,beta);
176: BVScaleColumn(eps->V,m,1.0 / *beta);
177: *breakdown = PETSC_FALSE;
179: if (m>100) { PetscFree(lhh); }
180: VecDestroy(&u);
181: VecDestroy(&t);
182: return(0);
183: }
187: /*
188: EPSDelayedArnoldi1 - This function is similar to EPSDelayedArnoldi,
189: but without reorthogonalization (only delayed normalization).
190: */
191: PetscErrorCode EPSDelayedArnoldi1(EPS eps,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,PetscReal *beta,PetscBool *breakdown)
192: {
194: PetscInt i,j,m=*M;
195: PetscScalar dot;
196: PetscReal norm=0.0;
197: Vec vj,vj1;
200: BVSetActiveColumns(eps->V,0,m);
201: for (j=k;j<m;j++) {
202: BVGetColumn(eps->V,j,&vj);
203: BVGetColumn(eps->V,j+1,&vj1);
204: STApply(eps->st,vj,vj1);
205: BVRestoreColumn(eps->V,j+1,&vj1);
206: BVDotColumnBegin(eps->V,j+1,H+ldh*j);
207: if (j>k) {
208: VecDotBegin(vj,vj,&dot);
209: }
210: BVDotColumnEnd(eps->V,j+1,H+ldh*j);
211: if (j>k) {
212: VecDotEnd(vj,vj,&dot);
213: }
214: BVRestoreColumn(eps->V,j,&vj);
216: if (j>k) {
217: norm = PetscSqrtReal(PetscRealPart(dot));
218: BVScaleColumn(eps->V,j,1.0/norm);
219: H[ldh*(j-1)+j] = norm;
221: for (i=0;i<j;i++)
222: H[ldh*j+i] = H[ldh*j+i]/norm;
223: H[ldh*j+j] = H[ldh*j+j]/dot;
224: BVScaleColumn(eps->V,j+1,1.0/norm);
225: *beta = norm;
226: }
227: BVMultColumn(eps->V,-1.0,1.0,j+1,H+ldh*j);
228: }
230: *breakdown = PETSC_FALSE;
231: return(0);
232: }
236: /*
237: EPSKrylovConvergence - Implements the loop that checks for convergence
238: in Krylov methods.
240: Input Parameters:
241: eps - the eigensolver; some error estimates are updated in eps->errest
242: getall - whether all residuals must be computed
243: kini - initial value of k (the loop variable)
244: nits - number of iterations of the loop
245: V - set of basis vectors (used only if trueresidual is activated)
246: nv - number of vectors to process (dimension of Q, columns of V)
247: beta - norm of f (the residual vector of the Arnoldi/Lanczos factorization)
248: corrf - correction factor for residual estimates (only in harmonic KS)
250: Output Parameters:
251: kout - the first index where the convergence test failed
252: */
253: PetscErrorCode EPSKrylovConvergence(EPS eps,PetscBool getall,PetscInt kini,PetscInt nits,PetscReal beta,PetscReal corrf,PetscInt *kout)
254: {
256: PetscInt k,newk,marker,ld,inside;
257: PetscScalar re,im,*Zr,*Zi,*X;
258: PetscReal resnorm;
259: PetscBool isshift,refined,istrivial;
260: Vec x,y,w[3];
263: RGIsTrivial(eps->rg,&istrivial);
264: if (eps->trueres) {
265: BVCreateVec(eps->V,&x);
266: BVCreateVec(eps->V,&y);
267: BVCreateVec(eps->V,&w[0]);
268: BVCreateVec(eps->V,&w[2]);
269: #if !defined(PETSC_USE_COMPLEX)
270: BVCreateVec(eps->V,&w[1]);
271: #else
272: w[1] = NULL;
273: #endif
274: }
275: DSGetLeadingDimension(eps->ds,&ld);
276: DSGetRefined(eps->ds,&refined);
277: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
278: marker = -1;
279: if (eps->trackall) getall = PETSC_TRUE;
280: for (k=kini;k<kini+nits;k++) {
281: /* eigenvalue */
282: re = eps->eigr[k];
283: im = eps->eigi[k];
284: if (!istrivial || eps->trueres || isshift || eps->conv==EPS_CONV_NORM) {
285: STBackTransform(eps->st,1,&re,&im);
286: }
287: if (!istrivial) {
288: RGCheckInside(eps->rg,1,&re,&im,&inside);
289: if (marker==-1 && inside<=0) marker = k;
290: if (!(eps->trueres || isshift || eps->conv==EPS_CONV_NORM)) { /* make sure eps->converged below uses the right value */
291: re = eps->eigr[k];
292: im = eps->eigi[k];
293: }
294: }
295: newk = k;
296: DSVectors(eps->ds,DS_MAT_X,&newk,&resnorm);
297: if (eps->trueres) {
298: DSGetArray(eps->ds,DS_MAT_X,&X);
299: Zr = X+k*ld;
300: if (newk==k+1) Zi = X+newk*ld;
301: else Zi = NULL;
302: EPSComputeRitzVector(eps,Zr,Zi,eps->V,x,y);
303: DSRestoreArray(eps->ds,DS_MAT_X,&X);
304: EPSComputeResidualNorm_Private(eps,re,im,x,y,w,&resnorm);
305: }
306: else if (!refined) resnorm *= beta*corrf;
307: /* error estimate */
308: (*eps->converged)(eps,re,im,resnorm,&eps->errest[k],eps->convergedctx);
309: if (marker==-1 && eps->errest[k] >= eps->tol) marker = k;
310: if (newk==k+1) {
311: eps->errest[k+1] = eps->errest[k];
312: k++;
313: }
314: if (marker!=-1 && !getall) break;
315: }
316: if (marker!=-1) k = marker;
317: *kout = k;
318: if (eps->trueres) {
319: VecDestroy(&x);
320: VecDestroy(&y);
321: VecDestroy(&w[0]);
322: VecDestroy(&w[2]);
323: #if !defined(PETSC_USE_COMPLEX)
324: VecDestroy(&w[1]);
325: #endif
326: }
327: return(0);
328: }
332: /*
333: EPSFullLanczos - Computes an m-step Lanczos factorization with full
334: reorthogonalization. At each Lanczos step, the corresponding Lanczos
335: vector is orthogonalized with respect to all previous Lanczos vectors.
336: This is equivalent to computing an m-step Arnoldi factorization and
337: exploting symmetry of the operator.
339: The first k columns are assumed to be locked and therefore they are
340: not modified. On exit, the following relation is satisfied:
342: OP * V - V * T = beta_m*v_m * e_m^T
344: where the columns of V are the Lanczos vectors (which are B-orthonormal),
345: T is a real symmetric tridiagonal matrix, and e_m is the m-th vector of
346: the canonical basis. The tridiagonal is stored as two arrays: alpha
347: contains the diagonal elements, beta the off-diagonal. On exit, the last
348: element of beta contains the B-norm of V[m] before normalization.
349: */
350: PetscErrorCode EPSFullLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown)
351: {
353: PetscInt j,m = *M;
354: Vec vj,vj1;
355: PetscScalar *hwork,lhwork[100];
358: if (m > 100) {
359: PetscMalloc1(m,&hwork);
360: } else hwork = lhwork;
362: BVSetActiveColumns(eps->V,0,m);
363: for (j=k;j<m;j++) {
364: BVGetColumn(eps->V,j,&vj);
365: BVGetColumn(eps->V,j+1,&vj1);
366: STApply(eps->st,vj,vj1);
367: BVRestoreColumn(eps->V,j,&vj);
368: BVRestoreColumn(eps->V,j+1,&vj1);
369: BVOrthogonalizeColumn(eps->V,j+1,hwork,beta+j,breakdown);
370: alpha[j] = PetscRealPart(hwork[j]);
371: if (*breakdown) {
372: *M = j+1;
373: break;
374: } else {
375: BVScaleColumn(eps->V,j+1,1.0/beta[j]);
376: }
377: }
378: if (m > 100) {
379: PetscFree(hwork);
380: }
381: return(0);
382: }
386: PetscErrorCode EPSPseudoLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,PetscReal *cos,Vec w)
387: {
389: PetscInt j,m = *M,i,ld,l;
390: Vec vj,vj1;
391: PetscScalar *hwork,lhwork[100];
392: PetscReal norm,norm1,norm2,t,*f,sym=0.0,fro=0.0;
393: PetscBLASInt j_,one=1;
396: DSGetLeadingDimension(eps->ds,&ld);
397: DSGetDimensions(eps->ds,NULL,NULL,&l,NULL,NULL);
398: if (cos) *cos = 1.0;
399: if (m > 100) {
400: PetscMalloc1(m,&hwork);
401: } else hwork = lhwork;
403: BVSetActiveColumns(eps->V,0,m);
404: for (j=k;j<m;j++) {
405: BVGetColumn(eps->V,j,&vj);
406: BVGetColumn(eps->V,j+1,&vj1);
407: STApply(eps->st,vj,vj1);
408: BVRestoreColumn(eps->V,j,&vj);
409: BVRestoreColumn(eps->V,j+1,&vj1);
410: BVOrthogonalizeColumn(eps->V,j+1,hwork,&norm,breakdown);
411: alpha[j] = PetscRealPart(hwork[j]);
412: beta[j] = PetscAbsReal(norm);
413: DSGetArrayReal(eps->ds,DS_MAT_T,&f);
414: if (j==k) {
415: for (i=l;i<j-1;i++) hwork[i]-= f[2*ld+i];
416: for (i=0;i<l;i++) hwork[i] = 0.0;
417: }
418: DSRestoreArrayReal(eps->ds,DS_MAT_T,&f);
419: hwork[j-1] -= beta[j-1];
420: PetscBLASIntCast(j,&j_);
421: sym = SlepcAbs(BLASnrm2_(&j_,hwork,&one),sym);
422: fro = SlepcAbs(fro,SlepcAbs(alpha[j],beta[j]));
423: if (j>0) fro = SlepcAbs(fro,beta[j-1]);
424: if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*eps->tol)) { *symmlost = PETSC_TRUE; *M=j+1; break; }
425: omega[j+1] = (norm<0.0)? -1.0: 1.0;
426: BVScaleColumn(eps->V,j+1,1.0/norm);
427: /* */
428: if (cos) {
429: BVGetColumn(eps->V,j+1,&vj1);
430: VecNorm(vj1,NORM_2,&norm1);
431: BVApplyMatrix(eps->V,vj1,w);
432: BVRestoreColumn(eps->V,j+1,&vj1);
433: VecNorm(w,NORM_2,&norm2);
434: t = 1.0/(norm1*norm2);
435: if (*cos>t) *cos = t;
436: }
437: }
438: if (m > 100) {
439: PetscFree(hwork);
440: }
441: return(0);
442: }