Actual source code: ks-harm.c
1: /*
3: SLEPc eigensolver: "krylovschur"
5: Method: Krylov-Schur with harmonic extraction
7: References:
9: [1] "Practical Implementation of Harmonic Krylov-Schur", SLEPc Technical
10: Report STR-9, available at http://www.grycap.upv.es/slepc.
12: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
13: SLEPc - Scalable Library for Eigenvalue Problem Computations
14: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
16: This file is part of SLEPc.
17:
18: SLEPc is free software: you can redistribute it and/or modify it under the
19: terms of version 3 of the GNU Lesser General Public License as published by
20: the Free Software Foundation.
22: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
23: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
24: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
25: more details.
27: You should have received a copy of the GNU Lesser General Public License
28: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
29: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30: */
32: #include <private/epsimpl.h> /*I "slepceps.h" I*/
33: #include <slepcblaslapack.h>
37: /*
38: EPSTranslateHarmonic - Computes a translation of the Krylov decomposition
39: in order to perform a harmonic extraction.
41: On input:
42: S is the Rayleigh quotient (order m, leading dimension is lds)
43: tau is the translation amount
44: b is assumed to be beta*e_m^T
46: On output:
47: g = (B-sigma*eye(m))'\b
48: S is updated as S + g*b'
50: Workspace:
51: work is workspace to store a working copy of S and the pivots (int
52: of length m)
53: */
54: PetscErrorCode EPSTranslateHarmonic(PetscInt m_,PetscScalar *S,PetscInt lds,PetscScalar tau,PetscScalar beta,PetscScalar *g,PetscScalar *work)
55: {
56: #if defined(PETSC_MISSING_LAPACK_GETRF) || defined(PETSC_MISSING_LAPACK_GETRS)
58: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF,GETRS - Lapack routines are unavailable.");
59: #else
61: PetscInt i,j;
62: PetscBLASInt info,m,one = 1;
63: PetscScalar *B = work;
64: PetscBLASInt *ipiv = (PetscBLASInt*)(work+m_*m_);
67: m = PetscBLASIntCast(m_);
68: /* Copy S to workspace B */
69: for (i=0;i<m;i++)
70: for (j=0;j<m;j++)
71: B[i+j*m] = S[i+j*lds];
72: /* Vector g initialy stores b */
73: PetscMemzero(g,m*sizeof(PetscScalar));
74: g[m-1] = beta;
75:
76: /* g = (B-sigma*eye(m))'\b */
77: for (i=0;i<m;i++)
78: B[i+i*m] -= tau;
79: LAPACKgetrf_(&m,&m,B,&m,ipiv,&info);
80: if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
81: if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
82: PetscLogFlops(2.0*m*m*m/3.0);
83: LAPACKgetrs_("C",&m,&one,B,&m,ipiv,g,&m,&info);
84: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
85: PetscLogFlops(2.0*m*m-m);
87: /* S = S + g*b' */
88: for (i=0;i<m;i++)
89: S[i+(m-1)*lds] = S[i+(m-1)*lds] + g[i]*beta;
90: return(0);
91: #endif
92: }
96: /*
97: EPSRecoverHarmonic - Computes a translation of the truncated Krylov
98: decomposition in order to recover the original non-translated state
100: On input:
101: S is the truncated Rayleigh quotient (size n, leading dimension m)
102: k and l indicate the active columns of S
103: [U, u] is the basis of the Krylov subspace
104: g is the vector computed in the original translation
105: Q is the similarity transformation used to reduce to sorted Schur form
107: On output:
108: S is updated as S + g*b'
109: u is re-orthonormalized with respect to U
110: b is re-scaled
111: g is destroyed
113: Workspace:
114: ghat is workspace to store a vector of length n
115: */
116: PetscErrorCode EPSRecoverHarmonic(PetscScalar *S,PetscInt n_,PetscInt k,PetscInt l,PetscInt m_,PetscScalar *g,PetscScalar *Q,Vec *U,Vec u,PetscScalar *ghat)
117: {
119: PetscBLASInt one=1,ncol=k+l,n,m;
120: PetscScalar done=1.0,dmone=-1.0,dzero=0.0;
121: PetscReal gamma,gnorm;
122: PetscBLASInt i,j;
125: n = PetscBLASIntCast(n_);
126: m = PetscBLASIntCast(m_);
128: /* g^ = -Q(:,idx)'*g */
129: BLASgemv_("C",&n,&ncol,&dmone,Q,&n,g,&one,&dzero,ghat,&one);
131: /* S = S + g^*b' */
132: for (i=0;i<k+l;i++) {
133: for (j=k;j<k+l;j++) {
134: S[i+j*m] += ghat[i]*S[k+l+j*m];
135: }
136: }
138: /* g~ = (I-Q(:,idx)*Q(:,idx)')*g = g+Q(:,idx)*g^ */
139: BLASgemv_("N",&n,&ncol,&done,Q,&n,ghat,&one,&done,g,&one);
141: /* gamma u^ = u - U*g~ */
142: SlepcVecMAXPBY(u,1.0,-1.0,m,g,U);
144: /* Renormalize u */
145: gnorm = 0.0;
146: for (i=0;i<n;i++)
147: gnorm = gnorm + PetscRealPart(g[i]*PetscConj(g[i]));
148: gamma = PetscSqrtReal(1.0+gnorm);
149: VecScale(u,1.0/gamma);
151: /* b = gamma*b */
152: for (i=k;i<k+l;i++) {
153: S[i*m+k+l] *= gamma;
154: }
155: return(0);
156: }
160: PetscErrorCode EPSSolve_KrylovSchur_Harmonic(EPS eps)
161: {
163: PetscInt i,k,l,lwork,nv;
164: Vec u=eps->work[0];
165: PetscScalar *S=eps->T,*Q,*g,*work;
166: PetscReal beta,gnorm;
167: PetscBool breakdown;
170: PetscMemzero(S,eps->ncv*eps->ncv*sizeof(PetscScalar));
171: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&Q);
172: lwork = PetscMax((eps->ncv+1)*eps->ncv,7*eps->ncv);
173: PetscMalloc(lwork*sizeof(PetscScalar),&work);
174: PetscMalloc(eps->ncv*sizeof(PetscScalar),&g);
176: /* Get the starting Arnoldi vector */
177: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
178: l = 0;
179:
180: /* Restart loop */
181: while (eps->reason == EPS_CONVERGED_ITERATING) {
182: eps->its++;
184: /* Compute an nv-step Arnoldi factorization */
185: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
186: EPSBasicArnoldi(eps,PETSC_FALSE,S,eps->ncv,eps->V,eps->nconv+l,&nv,u,&beta,&breakdown);
187: VecScale(u,1.0/beta);
189: /* Compute translation of Krylov decomposition */
190: EPSTranslateHarmonic(nv,S,eps->ncv,eps->target,(PetscScalar)beta,g,work);
191: gnorm = 0.0;
192: for (i=0;i<nv;i++)
193: gnorm = gnorm + PetscRealPart(g[i]*PetscConj(g[i]));
195: /* Solve projected problem and compute residual norm estimates */
196: EPSProjectedKSNonsym(eps,l,S,eps->ncv,Q,nv);
198: /* Check convergence */
199: EPSKrylovConvergence(eps,PETSC_FALSE,eps->trackall,eps->nconv,nv-eps->nconv,S,eps->ncv,Q,eps->V,nv,beta,PetscSqrtReal(1.0+gnorm),&k,work);
200: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
201: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
202:
203: /* Update l */
204: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
205: else {
206: l = (nv-k)/2;
207: #if !defined(PETSC_USE_COMPLEX)
208: if (S[(k+l-1)*(eps->ncv+1)+1] != 0.0) {
209: if (k+l<nv-1) l = l+1;
210: else l = l-1;
211: }
212: #endif
213: }
214:
215: if (eps->reason == EPS_CONVERGED_ITERATING) {
216: if (breakdown) {
217: /* Start a new Arnoldi factorization */
218: PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);
219: EPSGetStartVector(eps,k,eps->V[k],&breakdown);
220: if (breakdown) {
221: eps->reason = EPS_DIVERGED_BREAKDOWN;
222: PetscInfo(eps,"Unable to generate more start vectors\n");
223: }
224: } else {
225: /* Prepare the Rayleigh quotient for restart */
226: for (i=k;i<k+l;i++) {
227: S[i*eps->ncv+k+l] = Q[(i+1)*nv-1]*beta;
228: }
229: EPSRecoverHarmonic(S,nv,k,l,eps->ncv,g,Q,eps->V,u,work);
230: }
231: }
232: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
233: SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,nv,PETSC_FALSE);
234:
235: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
236: VecCopy(u,eps->V[k+l]);
237: }
238: eps->nconv = k;
239: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
240: }
242: PetscFree(Q);
243: PetscFree(work);
244: PetscFree(g);
245: return(0);
246: }