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-2012, Universitat Politecnica de Valencia, Spain
11: This file is part of SLEPc.
12:
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)
33: {
35: PetscInt j,m = *M;
36: PetscReal norm,min=1.0;
37: PetscScalar *hwork,lhwork[100];
40: if (m > 100) {
41: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&hwork);
42: } else {
43: hwork = lhwork;
44: }
46: for (j=k;j<m-1;j++) {
47: STApply(eps->OP,V[j],V[j+1]);
48: IPPseudoOrthogonalize(eps->ip,j+1,V,omega,V[j+1],hwork,&norm,breakdown);
49: alpha[j] = PetscRealPart(hwork[j]);
50: beta[j] = PetscAbsReal(norm);
51: min = PetscMin(min,beta[j]);
52: omega[j+1] = (norm<0.0)?-1:1;
53:
54: if (breakdown && *breakdown ) {
55: *M = j+1;
56: if (m > 100) {
57: PetscFree(hwork);
58: }
59: return(0);
60: } else {
61: VecScale(V[j+1],1.0/norm);
62: }
63: }
64: STApply(eps->OP,V[m-1],f);
65: IPPseudoOrthogonalize(eps->ip,m,V,omega,f,hwork,&norm,PETSC_NULL);
66: alpha[m-1] = PetscRealPart(hwork[m-1]);
67: beta[m-1] =PetscAbsReal(norm);
68: omega[m] = (norm<0.0)?-1:1;
69: if (m > 100) {
70: PetscFree(hwork);
71: }
72: return(0);
73: }
77: PetscErrorCode EPSSolve_KrylovSchur_Indefinite(EPS eps)
78: {
79: PetscErrorCode ierr;
80: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
81: PetscInt i,k,l,ld,nv;
82: Vec u=eps->work[0];
83: PetscScalar *Q;
84: PetscReal *a,*b,*r,beta,norm,*omega;
85: PetscBool breakdown=PETSC_FALSE;
88: DSGetLeadingDimension(eps->ds,&ld);
89: /* Get the starting Lanczos vector */
90: SlepcVecSetRandom(eps->V[0],eps->rand);
91: IPNorm(eps->ip,eps->V[0],&norm);
92: if (norm==0.0) SETERRQ(((PetscObject)eps)->comm,1,"Initial vector is zero or belongs to the deflation space");
93: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
94: omega[0] = (norm > 0)?1.0:-1.0;
95: beta = PetscAbsReal(norm);
96: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
97: VecScale(eps->V[0],1.0/norm);
98: l = 0;
99:
100: /* Restart loop */
101: while (eps->reason == EPS_CONVERGED_ITERATING) {
102: eps->its++;
104: /* Compute an nv-step Lanczos factorization */
105: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
106: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
107: b = a + ld;
108: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
109: EPSFullLanczosIndef(eps,a,b,omega,eps->V,eps->nconv+l,&nv,u,PETSC_NULL);
110: beta = b[nv-1];
111: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
112: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
113: DSSetDimensions(eps->ds,nv,PETSC_IGNORE,eps->nconv,eps->nconv+l);
114: if (l==0) {
115: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
116: } else {
117: DSSetState(eps->ds,DS_STATE_RAW);
118: }
120: /* Solve projected problem */
121: DSSolve(eps->ds,eps->eigr,eps->eigi);
122: DSSort(eps->ds,eps->eigr,eps->eigi,PETSC_NULL,PETSC_NULL,PETSC_NULL);
124: /* Check convergence */
125: EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,eps->V,nv,beta,1.0,&k);
126: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
127: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
128:
129: /* Update l */
130: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
131: else {
132: l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
133: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
134: if(*(a+ld+k+l-1)!=0){
135: if (k+l<nv-1) l = l+1;
136: else l = l-1;
137: }
138: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
139: }
140:
141: if (eps->reason == EPS_CONVERGED_ITERATING) {
142: if (breakdown) {
143: SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Breakdown in Indefinite Krylov-Schur (beta=%g)",beta);
144: } else {
145: /* Prepare the Rayleigh quotient for restart */
146: DSGetArray(eps->ds,DS_MAT_Q,&Q);
147: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
148: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
149: r = a + 2*ld;
150: for (i=k;i<k+l;i++) {
151: r[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
152: }
153: omega[k+l] = omega[nv];
154: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
155: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
156: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
157: }
158: }
159: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
160: DSGetArray(eps->ds,DS_MAT_Q,&Q);
161: SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,ld,PETSC_FALSE);
162: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
164: /* Normalize u and append it to V */
165: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
166: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
167: VecAXPBY(eps->V[k+l],1.0/(beta*omega[k+l]),0.0,u);
168: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
169: }
171: EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv);
172: eps->nconv = k;
174: }
175: DSSetDimensions(eps->ds,eps->nconv,PETSC_IGNORE,0,0);
176: DSSetState(eps->ds,DS_STATE_RAW);
177: return(0);
178: }