Actual source code: ks-indef.c
1: /*
3: SLEPc eigensolver: "krylovschur"
5: Method: Krylov-Schur for symmetric-indefinite eigenproblems
7: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
8: SLEPc - Scalable Library for Eigenvalue Problem Computations
9: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
11: This file is part of SLEPc.
13: SLEPc is free software: you can redistribute it and/or modify it under the
14: terms of version 3 of the GNU Lesser General Public License as published by
15: the Free Software Foundation.
17: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
18: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
20: more details.
22: You should have received a copy of the GNU Lesser General Public License
23: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
24: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
25: */
26: #include <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
27: #include <slepcblaslapack.h>
28: #include krylovschur.h
32: static PetscErrorCode EPSFullLanczosIndef(EPS eps,PetscReal *alpha,PetscReal *beta,PetscReal *omega,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscBool *breakdown,PetscReal *cos,Vec w)
33: {
35: PetscInt j,m = *M;
36: PetscScalar *hwork,lhwork[100];
37: PetscReal norm,norm1,norm2,t;
40: if (cos) *cos=1.0;
41: if (m > 100) {
42: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&hwork);
43: } else hwork = lhwork;
45: for (j=k;j<m-1;j++) {
46: STApply(eps->st,V[j],V[j+1]);
47: IPPseudoOrthogonalize(eps->ip,j+1,V,omega,V[j+1],hwork,&norm,breakdown);
48: VecScale(V[j+1],1.0/norm);
49: alpha[j] = PetscRealPart(hwork[j]);
50: beta[j] = PetscAbsReal(norm);
51: omega[j+1] = (norm<0.0)?-1.0:1.0;
52: /* */
53: VecNorm(V[j+1],NORM_2,&norm1);
54: IPApplyMatrix(eps->ip,V[j+1],w);
55: VecNorm(w,NORM_2,&norm2);
56: t=1/(norm1*norm2);
57: if (cos && *cos>t) *cos = t;
58: }
59: STApply(eps->st,V[m-1],f);
60: IPPseudoOrthogonalize(eps->ip,m,V,omega,f,hwork,&norm,NULL);
61: VecScale(f,1.0/norm);
62: alpha[m-1] = PetscRealPart(hwork[m-1]);
63: beta[m-1] =PetscAbsReal(norm);
64: omega[m] = (norm<0.0)?-1:1;
65: if (m > 100) {
66: PetscFree(hwork);
67: }
68: return(0);
69: }
73: PetscErrorCode EPSSolve_KrylovSchur_Indefinite(EPS eps)
74: {
75: PetscErrorCode ierr;
76: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
77: PetscInt i,k,l,ld,nv,t;
78: Vec u=eps->work[0],w=eps->work[1];
79: PetscScalar *Q;
80: PetscReal *a,*b,*r,beta,beta1,beta2,norm,*omega;
81: PetscBool breakdown=PETSC_FALSE;
84: DSGetLeadingDimension(eps->ds,&ld);
86: /* Get the starting Lanczos vector */
87: SlepcVecSetRandom(eps->V[0],eps->rand);
88: IPNorm(eps->ip,eps->V[0],&norm);
89: if (norm==0.0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Initial vector is zero or belongs to the deflation space");
90: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
91: omega[0] = (norm > 0)?1.0:-1.0;
92: beta = PetscAbsReal(norm);
93: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
94: VecScale(eps->V[0],1.0/norm);
95: l = 0;
97: /* Restart loop */
98: while (eps->reason == EPS_CONVERGED_ITERATING) {
99: eps->its++;
101: /* Compute an nv-step Lanczos factorization */
102: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
103: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
104: b = a + ld;
105: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
106: EPSFullLanczosIndef(eps,a,b,omega,eps->V,eps->nconv+l,&nv,u,&breakdown,NULL,w);
107: beta = b[nv-1];
108: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
109: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
110: DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
111: if (l==0) {
112: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
113: } else {
114: DSSetState(eps->ds,DS_STATE_RAW);
115: }
116: /* Solve projected problem */
117: DSSolve(eps->ds,eps->eigr,eps->eigi);
118: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
120: /* Check convergence */
121: DSGetDimensions(eps->ds,NULL,NULL,NULL,NULL,&t);
122: VecNorm(u,NORM_2,&beta1);
123: IPApplyMatrix(eps->ip,u,w);
124: VecNorm(w,NORM_2,&beta2);
125: beta1 = PetscMax(beta1,beta2);
126: EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,t-eps->nconv,eps->V,nv,beta*beta1,1.0,&k);
127: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
128: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
130: /* Update l */
131: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
132: else {
133: l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
134: l = PetscMin(l,t);
135: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
136: if (*(a+ld+k+l-1)!=0) {
137: if (k+l<t-1) l = l+1;
138: else l = l-1;
139: }
140: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
141: }
143: if (eps->reason == EPS_CONVERGED_ITERATING) {
144: if (breakdown) {
145: SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Breakdown in Indefinite Krylov-Schur (beta=%g)",beta);
146: } else {
147: /* Prepare the Rayleigh quotient for restart */
148: DSGetArray(eps->ds,DS_MAT_Q,&Q);
149: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
150: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
151: b = a + ld;
152: r = a + 2*ld;
153: for (i=k;i<k+l;i++) {
154: r[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
155: }
156: b[k+l-1] = r[k+l-1];
157: omega[k+l] = omega[nv];
158: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
159: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
160: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
161: }
162: }
163: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
164: DSGetArray(eps->ds,DS_MAT_Q,&Q);
165: SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,ld,PETSC_FALSE);
166: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
168: /* Append u to V */
169: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
170: VecCopy(u,eps->V[k+l]);
171: }
173: EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv);
174: eps->nconv = k;
175: }
176: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
177: DSSetState(eps->ds,DS_STATE_RAW);
178: return(0);
179: }