Actual source code: ks-indef.c
slepc-3.5.2 2014-10-10
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-2014, 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>
27: #include krylovschur.h
31: PetscErrorCode EPSSolve_KrylovSchur_Indefinite(EPS eps)
32: {
33: PetscErrorCode ierr;
34: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
35: PetscInt i,k,l,ld,nv,t;
36: Mat U;
37: Vec vomega,u,w=eps->work[0];
38: PetscScalar *Q,*aux;
39: PetscReal *a,*b,*r,beta,beta1,beta2,*omega;
40: PetscBool breakdown=PETSC_FALSE;
43: DSGetLeadingDimension(eps->ds,&ld);
45: /* Get the starting Lanczos vector */
46: EPSGetStartVector(eps,0,NULL);
48: /* Extract sigma[0] from BV, computed during normalization */
49: BVSetActiveColumns(eps->V,0,1);
50: VecCreateSeq(PETSC_COMM_SELF,1,&vomega);
51: BVGetSignature(eps->V,vomega);
52: VecGetArray(vomega,&aux);
53: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
54: omega[0] = PetscRealPart(aux[0]);
55: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
56: VecRestoreArray(vomega,&aux);
57: VecDestroy(&vomega);
58: l = 0;
60: /* Restart loop */
61: while (eps->reason == EPS_CONVERGED_ITERATING) {
62: eps->its++;
64: /* Compute an nv-step Lanczos factorization */
65: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
66: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
67: b = a + ld;
68: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
69: EPSPseudoLanczos(eps,a,b,omega,eps->nconv+l,&nv,&breakdown,NULL,w);
70: beta = b[nv-1];
71: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
72: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
73: DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
74: if (l==0) {
75: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
76: } else {
77: DSSetState(eps->ds,DS_STATE_RAW);
78: }
79: BVSetActiveColumns(eps->V,eps->nconv,nv);
81: /* Solve projected problem */
82: DSSolve(eps->ds,eps->eigr,eps->eigi);
83: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
85: /* Check convergence */
86: DSGetDimensions(eps->ds,NULL,NULL,NULL,NULL,&t);
87: BVGetColumn(eps->V,nv,&u);
88: VecNorm(u,NORM_2,&beta1);
89: BVRestoreColumn(eps->V,nv,&u);
90: VecNorm(w,NORM_2,&beta2); /* w contains B*V[nv] */
91: beta1 = PetscMax(beta1,beta2);
92: EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,t-eps->nconv,beta*beta1,1.0,&k);
93: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
94: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
96: /* Update l */
97: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
98: else {
99: l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
100: l = PetscMin(l,t);
101: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
102: if (*(a+ld+k+l-1)!=0) {
103: if (k+l<t-1) l = l+1;
104: else l = l-1;
105: }
106: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
107: }
109: if (eps->reason == EPS_CONVERGED_ITERATING) {
110: if (breakdown) {
111: SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Breakdown in Indefinite Krylov-Schur (beta=%g)",beta);
112: } else {
113: /* Prepare the Rayleigh quotient for restart */
114: DSGetArray(eps->ds,DS_MAT_Q,&Q);
115: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
116: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
117: b = a + ld;
118: r = a + 2*ld;
119: for (i=k;i<k+l;i++) {
120: r[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
121: }
122: b[k+l-1] = r[k+l-1];
123: omega[k+l] = omega[nv];
124: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
125: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
126: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
127: }
128: }
129: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
130: DSGetMat(eps->ds,DS_MAT_Q,&U);
131: BVMultInPlace(eps->V,U,eps->nconv,k+l);
132: MatDestroy(&U);
134: /* Move restart vector and update signature */
135: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
136: BVCopyColumn(eps->V,nv,k+l);
137: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
138: VecCreateSeq(PETSC_COMM_SELF,k+l,&vomega);
139: VecGetArray(vomega,&aux);
140: for (i=0;i<k+l;i++) aux[i] = omega[i];
141: VecRestoreArray(vomega,&aux);
142: BVSetActiveColumns(eps->V,0,k+l);
143: BVSetSignature(eps->V,vomega);
144: VecDestroy(&vomega);
145: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
146: }
148: EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv);
149: eps->nconv = k;
150: }
151: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
152: DSSetState(eps->ds,DS_STATE_RAW);
153: return(0);
154: }