Actual source code: krylovschur.c
1: /*
3: SLEPc eigensolver: "krylovschur"
5: Method: Krylov-Schur
7: Algorithm:
9: Single-vector Krylov-Schur method for both symmetric and non-symmetric
10: problems.
12: References:
14: [1] "Krylov-Schur Methods in SLEPc", SLEPc Technical Report STR-7,
15: available at http://www.grycap.upv.es/slepc.
17: [2] G.W. Stewart, "A Krylov-Schur Algorithm for Large Eigenproblems",
18: SIAM J. Matrix Analysis and App., 23(3), pp. 601-614, 2001.
20: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
21: SLEPc - Scalable Library for Eigenvalue Problem Computations
22: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
24: This file is part of SLEPc.
25:
26: SLEPc is free software: you can redistribute it and/or modify it under the
27: terms of version 3 of the GNU Lesser General Public License as published by
28: the Free Software Foundation.
30: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
31: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
32: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
33: more details.
35: You should have received a copy of the GNU Lesser General Public License
36: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: */
40: #include <private/epsimpl.h> /*I "slepceps.h" I*/
41: #include <slepcblaslapack.h>
43: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS);
44: extern PetscErrorCode EPSSolve_KrylovSchur_Harmonic(EPS);
45: extern PetscErrorCode EPSSolve_KrylovSchur_Symm(EPS);
46: extern PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS);
50: PetscErrorCode EPSSetUp_KrylovSchur(EPS eps)
51: {
53: PetscBool issinv;
56: /* spectrum slicing requires special treatment of default values */
57: if (eps->which==EPS_ALL) {
58: if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(((PetscObject)eps)->comm,1,"Must define a computational interval when using EPS_ALL");
59: if (!eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,1,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
60: if (!((PetscObject)(eps->OP))->type_name) { /* default to shift-and-invert */
61: STSetType(eps->OP,STSINVERT);
62: }
63: PetscTypeCompareAny((PetscObject)eps->OP,&issinv,STSINVERT,STCAYLEY,"");
64: if (!issinv) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
65: #if defined(PETSC_USE_REAL_DOUBLE)
66: if (eps->tol==PETSC_DEFAULT) eps->tol = 1e-10; /* use tighter tolerance */
67: #endif
68: if (eps->intb >= PETSC_MAX_REAL) { /* right-open interval */
69: if (eps->inta <= PETSC_MIN_REAL) SETERRQ(((PetscObject)eps)->comm,1,"The defined computational interval should have at least one of their sides bounded");
70: STSetDefaultShift(eps->OP,eps->inta);
71: }
72: else { STSetDefaultShift(eps->OP,eps->intb); }
74: if (eps->nev==1) eps->nev = 40; /* nev not set, use default value */
75: if (eps->nev<10) SETERRQ(((PetscObject)eps)->comm,1,"nev cannot be less than 10 in spectrum slicing runs");
76: eps->ops->backtransform = PETSC_NULL;
77: }
79: /* proceed with the general case */
80: if (eps->ncv) { /* ncv set */
81: if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
82: } else if (eps->mpd) { /* mpd set */
83: eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
84: } else { /* neither set: defaults depend on nev being small or large */
85: if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
86: else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
87: }
88: if (!eps->mpd) eps->mpd = eps->ncv;
89: if (eps->ncv>eps->nev+eps->mpd) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must not be larger than nev+mpd");
90: if (!eps->max_it) {
91: if (eps->which==EPS_ALL) eps->max_it = 100; /* special case for spectrum slicing */
92: else eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
93: }
94: if (!eps->which) { EPSDefaultSetWhich(eps); }
95: if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
96: SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
98: if (!eps->extraction) {
99: EPSSetExtraction(eps,EPS_RITZ);
100: } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC) {
101: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
102: }
104: EPSAllocateSolution(eps);
105: PetscFree(eps->T);
106: if (!eps->ishermitian || eps->extraction==EPS_HARMONIC) {
107: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
108: }
109: EPSDefaultGetWork(eps,1);
111: /* dispatch solve method */
112: if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
113: if (eps->ishermitian) {
114: if (eps->which==EPS_ALL) eps->ops->solve = EPSSolve_KrylovSchur_Slice;
115: else {
116: switch (eps->extraction) {
117: case EPS_RITZ: eps->ops->solve = EPSSolve_KrylovSchur_Symm; break;
118: case EPS_HARMONIC: eps->ops->solve = EPSSolve_KrylovSchur_Harmonic; break;
119: default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
120: }
121: }
122: } else {
123: switch (eps->extraction) {
124: case EPS_RITZ: eps->ops->solve = EPSSolve_KrylovSchur_Default; break;
125: case EPS_HARMONIC: eps->ops->solve = EPSSolve_KrylovSchur_Harmonic; break;
126: default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
127: }
128: }
129: return(0);
130: }
134: /*
135: EPSProjectedKSNonsym - Solves the projected eigenproblem in the Krylov-Schur
136: method (non-symmetric case).
138: On input:
139: l is the number of vectors kept in previous restart (0 means first restart)
140: S is the projected matrix (leading dimension is lds)
142: On output:
143: S has (real) Schur form with diagonal blocks sorted appropriately
144: Q contains the corresponding Schur vectors (order n, leading dimension n)
145: */
146: PetscErrorCode EPSProjectedKSNonsym(EPS eps,PetscInt l,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt n)
147: {
149: PetscInt i;
152: if (l==0) {
153: PetscMemzero(Q,n*n*sizeof(PetscScalar));
154: for (i=0;i<n;i++)
155: Q[i*(n+1)] = 1.0;
156: } else {
157: /* Reduce S to Hessenberg form, S <- Q S Q' */
158: EPSDenseHessenberg(n,eps->nconv,S,lds,Q);
159: }
160: /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
161: EPSDenseSchur(n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);
162: /* Sort the remaining columns of the Schur form */
163: EPSSortDenseSchur(eps,n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);
164: return(0);
165: }
169: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
170: {
172: PetscInt i,k,l,lwork,nv;
173: Vec u=eps->work[0];
174: PetscScalar *S=eps->T,*Q,*work;
175: PetscReal beta;
176: PetscBool breakdown;
179: PetscMemzero(S,eps->ncv*eps->ncv*sizeof(PetscScalar));
180: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&Q);
181: lwork = 7*eps->ncv;
182: PetscMalloc(lwork*sizeof(PetscScalar),&work);
184: /* Get the starting Arnoldi vector */
185: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
186: l = 0;
187:
188: /* Restart loop */
189: while (eps->reason == EPS_CONVERGED_ITERATING) {
190: eps->its++;
192: /* Compute an nv-step Arnoldi factorization */
193: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
194: EPSBasicArnoldi(eps,PETSC_FALSE,S,eps->ncv,eps->V,eps->nconv+l,&nv,u,&beta,&breakdown);
195: VecScale(u,1.0/beta);
197: /* Solve projected problem */
198: EPSProjectedKSNonsym(eps,l,S,eps->ncv,Q,nv);
200: /* Check convergence */
201: EPSKrylovConvergence(eps,PETSC_FALSE,eps->trackall,eps->nconv,nv-eps->nconv,S,eps->ncv,Q,eps->V,nv,beta,1.0,&k,work);
202: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
203: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
204:
205: /* Update l */
206: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
207: else {
208: l = (nv-k)/2;
209: #if !defined(PETSC_USE_COMPLEX)
210: if (S[(k+l-1)*(eps->ncv+1)+1] != 0.0) {
211: if (k+l<nv-1) l = l+1;
212: else l = l-1;
213: }
214: #endif
215: }
217: if (eps->reason == EPS_CONVERGED_ITERATING) {
218: if (breakdown) {
219: /* Start a new Arnoldi factorization */
220: PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);
221: EPSGetStartVector(eps,k,eps->V[k],&breakdown);
222: if (breakdown) {
223: eps->reason = EPS_DIVERGED_BREAKDOWN;
224: PetscInfo(eps,"Unable to generate more start vectors\n");
225: }
226: } else {
227: /* Prepare the Rayleigh quotient for restart */
228: for (i=k;i<k+l;i++) {
229: S[i*eps->ncv+k+l] = Q[(i+1)*nv-1]*beta;
230: }
231: }
232: }
233: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
234: SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,nv,PETSC_FALSE);
236: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
237: VecCopy(u,eps->V[k+l]);
238: }
239: eps->nconv = k;
241: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
242: }
243: PetscFree(Q);
244: PetscFree(work);
245: return(0);
246: }
250: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
251: {
255: PetscFree(eps->T);
256: EPSReset_Default(eps);
257: return(0);
258: }
260: EXTERN_C_BEGIN
263: PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
264: {
266: eps->ops->setup = EPSSetUp_KrylovSchur;
267: eps->ops->reset = EPSReset_KrylovSchur;
268: eps->ops->backtransform = EPSBackTransform_Default;
269: eps->ops->computevectors = EPSComputeVectors_Schur;
270: return(0);
271: }
272: EXTERN_C_END