Actual source code: ks-indef.c
slepc-3.7.3 2016-09-29
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-2016, 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,nconv=0;
36: Mat U;
37: Vec vomega,w=eps->work[0];
38: PetscScalar *Q,*aux;
39: PetscReal *a,*b,*r,beta,beta1=1.0,*omega;
40: PetscBool breakdown=PETSC_FALSE,symmlost=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,&symmlost,NULL,w);
70: if (symmlost) {
71: eps->reason = EPS_DIVERGED_SYMMETRY_LOST;
72: if (nv==eps->nconv+l+1) { eps->nconv = nconv; break; }
73: }
74: beta = b[nv-1];
75: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
76: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
77: DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
78: if (l==0) {
79: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
80: } else {
81: DSSetState(eps->ds,DS_STATE_RAW);
82: }
83: BVSetActiveColumns(eps->V,eps->nconv,nv);
85: /* Solve projected problem */
86: DSSolve(eps->ds,eps->eigr,eps->eigi);
87: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
89: /* Check convergence */
90: DSGetDimensions(eps->ds,NULL,NULL,NULL,NULL,&t);
91: #if 0
92: /* take into account also left residual */
93: BVGetColumn(eps->V,nv,&u);
94: VecNorm(u,NORM_2,&beta1);
95: BVRestoreColumn(eps->V,nv,&u);
96: VecNorm(w,NORM_2,&beta2); /* w contains B*V[nv] */
97: beta1 = PetscMax(beta1,beta2);
98: #endif
99: EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,t-eps->nconv,beta*beta1,1.0,&k);
100: (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
101: nconv = k;
103: /* Update l */
104: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
105: else {
106: l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
107: l = PetscMin(l,t);
108: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
109: if (*(a+ld+k+l-1)!=0) {
110: if (k+l<t-1) l = l+1;
111: else l = l-1;
112: }
113: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
114: }
115: if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
117: if (eps->reason == EPS_CONVERGED_ITERATING) {
118: if (breakdown) {
119: SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Breakdown in Indefinite Krylov-Schur (beta=%g)",beta);
120: } else {
121: /* Prepare the Rayleigh quotient for restart */
122: DSGetArray(eps->ds,DS_MAT_Q,&Q);
123: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
124: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
125: b = a + ld;
126: r = a + 2*ld;
127: for (i=k;i<k+l;i++) {
128: r[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
129: }
130: b[k+l-1] = r[k+l-1];
131: omega[k+l] = omega[nv];
132: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
133: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
134: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
135: }
136: }
137: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
138: DSGetMat(eps->ds,DS_MAT_Q,&U);
139: BVMultInPlace(eps->V,U,eps->nconv,k+l);
140: MatDestroy(&U);
142: /* Move restart vector and update signature */
143: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
144: BVCopyColumn(eps->V,nv,k+l);
145: DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
146: VecCreateSeq(PETSC_COMM_SELF,k+l,&vomega);
147: VecGetArray(vomega,&aux);
148: for (i=0;i<k+l;i++) aux[i] = omega[i];
149: VecRestoreArray(vomega,&aux);
150: BVSetActiveColumns(eps->V,0,k+l);
151: BVSetSignature(eps->V,vomega);
152: VecDestroy(&vomega);
153: DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
154: }
156: eps->nconv = k;
157: EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
158: }
159: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
160: DSSetState(eps->ds,DS_STATE_RAW);
161: return(0);
162: }