Actual source code: krylovschur.c
1: /*
3: SLEPc eigensolver: "krylovschur"
5: Method: Krylov-Schur
7: Algorithm:
9: Single-vector Krylov-Schur method for non-symmetric problems,
10: including harmonic extraction.
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: [3] "Practical Implementation of Harmonic Krylov-Schur", SLEPc Technical
21: Report STR-9, available at http://www.grycap.upv.es/slepc.
23: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24: SLEPc - Scalable Library for Eigenvalue Problem Computations
25: Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain
27: This file is part of SLEPc.
28:
29: SLEPc is free software: you can redistribute it and/or modify it under the
30: terms of version 3 of the GNU Lesser General Public License as published by
31: the Free Software Foundation.
33: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
34: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
35: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
36: more details.
38: You should have received a copy of the GNU Lesser General Public License
39: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
40: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41: */
43: #include <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
44: #include <slepcblaslapack.h>
45: #include krylovschur.h
49: PetscErrorCode EPSGetArbitraryValues(EPS eps,PetscScalar *rr,PetscScalar *ri)
50: {
52: PetscInt i,newi,ld,n,l;
53: Vec xr=eps->work[1],xi=eps->work[2];
54: PetscScalar re,im,*Zr,*Zi,*X;
57: DSGetLeadingDimension(eps->ds,&ld);
58: DSGetDimensions(eps->ds,&n,PETSC_NULL,&l,PETSC_NULL);
59: for (i=l;i<n;i++) {
60: re = eps->eigr[i];
61: im = eps->eigi[i];
62: STBackTransform(eps->OP,1,&re,&im);
63: newi = i;
64: DSVectors(eps->ds,DS_MAT_X,&newi,PETSC_NULL);
65: DSGetArray(eps->ds,DS_MAT_X,&X);
66: Zr = X+i*ld;
67: if (newi==i+1) Zi = X+newi*ld;
68: else Zi = PETSC_NULL;
69: EPSComputeRitzVector(eps,Zr,Zi,eps->V,n,xr,xi);
70: DSRestoreArray(eps->ds,DS_MAT_X,&X);
71: (*eps->arbit_func)(re,im,xr,xi,rr+i,ri+i,eps->arbit_ctx);
72: }
73: return(0);
74: }
78: PetscErrorCode EPSSetUp_KrylovSchur(EPS eps)
79: {
80: PetscErrorCode ierr;
81: PetscBool issinv;
82: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
83: enum { EPS_KS_DEFAULT, EPS_KS_SYMM, EPS_KS_SLICE, EPS_KS_INDEF } variant;
86: /* spectrum slicing requires special treatment of default values */
87: if (eps->which==EPS_ALL) {
88: if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONG,"Must define a computational interval when using EPS_ALL");
89: if (!eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
90: if (eps->arbit_func) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Arbitrary selection of eigenpairs cannot be used with spectrum slicing");
91: if (!((PetscObject)(eps->OP))->type_name) { /* default to shift-and-invert */
92: STSetType(eps->OP,STSINVERT);
93: }
94: PetscObjectTypeCompareAny((PetscObject)eps->OP,&issinv,STSINVERT,STCAYLEY,"");
95: if (!issinv) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
96: #if defined(PETSC_USE_REAL_DOUBLE)
97: if (eps->tol==PETSC_DEFAULT) eps->tol = 1e-10; /* use tighter tolerance */
98: #endif
99: if (eps->intb >= PETSC_MAX_REAL) { /* right-open interval */
100: if (eps->inta <= PETSC_MIN_REAL) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONG,"The defined computational interval should have at least one of their sides bounded");
101: STSetDefaultShift(eps->OP,eps->inta);
102: }
103: else { STSetDefaultShift(eps->OP,eps->intb); }
105: if (eps->nev==1) eps->nev = 40; /* nev not set, use default value */
106: if (eps->nev<10) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONG,"nev cannot be less than 10 in spectrum slicing runs");
107: eps->ops->backtransform = PETSC_NULL;
108: }
110: if (eps->isgeneralized && eps->ishermitian && !eps->ispositive && eps->arbit_func) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not implemented for indefinite problems");
112: /* proceed with the general case */
113: if (eps->ncv) { /* ncv set */
114: if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
115: } else if (eps->mpd) { /* mpd set */
116: eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
117: } else { /* neither set: defaults depend on nev being small or large */
118: if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
119: else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
120: }
121: if (!eps->mpd) eps->mpd = eps->ncv;
122: if (eps->ncv>eps->nev+eps->mpd) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must not be larger than nev+mpd");
123: if (!eps->max_it) {
124: if (eps->which==EPS_ALL) eps->max_it = 100; /* special case for spectrum slicing */
125: else eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
126: }
127: if (!eps->which) { EPSDefaultSetWhich(eps); }
128: if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY)) SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
130: if (!eps->extraction) {
131: EPSSetExtraction(eps,EPS_RITZ);
132: } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC)
133: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
135: if (!ctx->keep) ctx->keep = 0.5;
137: EPSAllocateSolution(eps);
138: if (eps->arbit_func) { EPSDefaultGetWork(eps,3); }
139: else { EPSDefaultGetWork(eps,1); }
141: /* dispatch solve method */
142: if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
143: if (eps->ishermitian) {
144: if (eps->which==EPS_ALL) {
145: if (eps->isgeneralized && !eps->ispositive) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Spectrum slicing not implemented for indefinite problems");
147: else variant = EPS_KS_SLICE;
148: } else if (eps->isgeneralized && !eps->ispositive) {
149: variant = EPS_KS_INDEF;
150: } else {
151: switch (eps->extraction) {
152: case EPS_RITZ: variant = EPS_KS_SYMM; break;
153: case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
154: default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
155: }
156: }
157: } else {
158: switch (eps->extraction) {
159: case EPS_RITZ: variant = EPS_KS_DEFAULT; break;
160: case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
161: default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
162: }
163: }
164: switch (variant) {
165: case EPS_KS_DEFAULT:
166: eps->ops->solve = EPSSolve_KrylovSchur_Default;
167: DSSetType(eps->ds,DSNHEP);
168: break;
169: case EPS_KS_SYMM:
170: eps->ops->solve = EPSSolve_KrylovSchur_Symm;
171: DSSetType(eps->ds,DSHEP);
172: DSSetCompact(eps->ds,PETSC_TRUE);
173: DSSetExtraRow(eps->ds,PETSC_TRUE);
174: break;
175: case EPS_KS_SLICE:
176: eps->ops->solve = EPSSolve_KrylovSchur_Slice;
177: DSSetType(eps->ds,DSHEP);
178: DSSetCompact(eps->ds,PETSC_TRUE);
179: break;
180: case EPS_KS_INDEF:
181: eps->ops->solve = EPSSolve_KrylovSchur_Indefinite;
182: DSSetType(eps->ds,DSGHIEP);
183: DSSetCompact(eps->ds,PETSC_TRUE);
184: break;
185: default: SETERRQ(((PetscObject)eps)->comm,1,"Unexpected error");
186: }
187: DSAllocate(eps->ds,eps->ncv+1);
188: return(0);
189: }
193: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
194: {
195: PetscErrorCode ierr;
196: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
197: PetscInt i,j,*pj,k,l,nv,ld;
198: Vec u=eps->work[0];
199: PetscScalar *S,*Q,*g;
200: PetscReal beta,gamma=1.0;
201: PetscBool breakdown,harmonic;
204: DSGetLeadingDimension(eps->ds,&ld);
205: harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
206: if (harmonic) { PetscMalloc(ld*sizeof(PetscScalar),&g); }
207: if (eps->arbit_func) pj = &j; else pj = PETSC_NULL;
209: /* Get the starting Arnoldi vector */
210: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
211: l = 0;
212:
213: /* Restart loop */
214: while (eps->reason == EPS_CONVERGED_ITERATING) {
215: eps->its++;
217: /* Compute an nv-step Arnoldi factorization */
218: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
219: DSGetArray(eps->ds,DS_MAT_A,&S);
220: EPSBasicArnoldi(eps,PETSC_FALSE,S,ld,eps->V,eps->nconv+l,&nv,u,&beta,&breakdown);
221: VecScale(u,1.0/beta);
222: DSRestoreArray(eps->ds,DS_MAT_A,&S);
223: DSSetDimensions(eps->ds,nv,PETSC_IGNORE,eps->nconv,eps->nconv+l);
224: if (l==0) {
225: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
226: } else {
227: DSSetState(eps->ds,DS_STATE_RAW);
228: }
230: /* Compute translation of Krylov decomposition if harmonic extraction used */
231: if (harmonic) {
232: DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,g,&gamma);
233: }
235: /* Solve projected problem */
236: DSSolve(eps->ds,eps->eigr,eps->eigi);
237: if (eps->arbit_func) { EPSGetArbitraryValues(eps,eps->rr,eps->ri); j=1; }
238: DSSort(eps->ds,eps->eigr,eps->eigi,eps->rr,eps->ri,pj);
240: /* Check convergence */
241: EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,eps->V,nv,beta,gamma,&k);
242: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
243: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
245: /* Update l */
246: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
247: else {
248: l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
249: #if !defined(PETSC_USE_COMPLEX)
250: DSGetArray(eps->ds,DS_MAT_A,&S);
251: if (S[k+l+(k+l-1)*ld] != 0.0) {
252: if (k+l<nv-1) l = l+1;
253: else l = l-1;
254: }
255: DSRestoreArray(eps->ds,DS_MAT_A,&S);
256: #endif
257: }
259: if (eps->reason == EPS_CONVERGED_ITERATING) {
260: if (breakdown) {
261: /* Start a new Arnoldi factorization */
262: PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);
263: EPSGetStartVector(eps,k,eps->V[k],&breakdown);
264: if (breakdown) {
265: eps->reason = EPS_DIVERGED_BREAKDOWN;
266: PetscInfo(eps,"Unable to generate more start vectors\n");
267: }
268: } else {
269: /* Undo translation of Krylov decomposition */
270: if (harmonic) {
271: DSSetDimensions(eps->ds,nv,PETSC_IGNORE,k,l);
272: DSTranslateHarmonic(eps->ds,0.0,beta,PETSC_TRUE,g,&gamma);
273: /* gamma u^ = u - U*g~ */
274: SlepcVecMAXPBY(u,1.0,-1.0,nv,g,eps->V);
275: VecScale(u,1.0/gamma);
276: }
277: /* Prepare the Rayleigh quotient for restart */
278: DSGetArray(eps->ds,DS_MAT_A,&S);
279: DSGetArray(eps->ds,DS_MAT_Q,&Q);
280: for (i=k;i<k+l;i++) {
281: S[k+l+i*ld] = Q[nv-1+i*ld]*beta*gamma;
282: }
283: DSRestoreArray(eps->ds,DS_MAT_A,&S);
284: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
285: }
286: }
287: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
288: DSGetArray(eps->ds,DS_MAT_Q,&Q);
289: SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,ld,PETSC_FALSE);
290: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
292: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
293: VecCopy(u,eps->V[k+l]);
294: }
295: eps->nconv = k;
296: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
297: }
299: if (harmonic) { PetscFree(g); }
300: /* truncate Schur decomposition and change the state to raw so that
301: PSVectors() computes eigenvectors from scratch */
302: DSSetDimensions(eps->ds,eps->nconv,PETSC_IGNORE,0,0);
303: DSSetState(eps->ds,DS_STATE_RAW);
304: return(0);
305: }
307: EXTERN_C_BEGIN
310: PetscErrorCode EPSKrylovSchurSetRestart_KrylovSchur(EPS eps,PetscReal keep)
311: {
312: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
315: if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
316: else {
317: if (keep<0.1 || keep>0.9) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
318: ctx->keep = keep;
319: }
320: return(0);
321: }
322: EXTERN_C_END
326: /*@
327: EPSKrylovSchurSetRestart - Sets the restart parameter for the Krylov-Schur
328: method, in particular the proportion of basis vectors that must be kept
329: after restart.
331: Logically Collective on EPS
333: Input Parameters:
334: + eps - the eigenproblem solver context
335: - keep - the number of vectors to be kept at restart
337: Options Database Key:
338: . -eps_krylovschur_restart - Sets the restart parameter
340: Notes:
341: Allowed values are in the range [0.1,0.9]. The default is 0.5.
342:
343: Level: advanced
345: .seealso: EPSKrylovSchurGetRestart()
346: @*/
347: PetscErrorCode EPSKrylovSchurSetRestart(EPS eps,PetscReal keep)
348: {
354: PetscTryMethod(eps,"EPSKrylovSchurSetRestart_C",(EPS,PetscReal),(eps,keep));
355: return(0);
356: }
358: EXTERN_C_BEGIN
361: PetscErrorCode EPSKrylovSchurGetRestart_KrylovSchur(EPS eps,PetscReal *keep)
362: {
363: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
366: *keep = ctx->keep;
367: return(0);
368: }
369: EXTERN_C_END
373: /*@
374: EPSKrylovSchurGetRestart - Gets the restart parameter used in the
375: Krylov-Schur method.
377: Not Collective
379: Input Parameter:
380: . eps - the eigenproblem solver context
382: Output Parameter:
383: . keep - the restart parameter
385: Level: advanced
387: .seealso: EPSKrylovSchurSetRestart()
388: @*/
389: PetscErrorCode EPSKrylovSchurGetRestart(EPS eps,PetscReal *keep)
390: {
396: PetscTryMethod(eps,"EPSKrylovSchurGetRestart_C",(EPS,PetscReal*),(eps,keep));
397: return(0);
398: }
402: PetscErrorCode EPSSetFromOptions_KrylovSchur(EPS eps)
403: {
405: PetscBool flg;
406: PetscReal keep;
409: PetscOptionsHead("EPS Krylov-Schur Options");
410: PetscOptionsReal("-eps_krylovschur_restart","Proportion of vectors kept after restart","EPSKrylovSchurSetRestart",0.5,&keep,&flg);
411: if (flg) { EPSKrylovSchurSetRestart(eps,keep); }
412: PetscOptionsTail();
413: return(0);
414: }
418: PetscErrorCode EPSView_KrylovSchur(EPS eps,PetscViewer viewer)
419: {
420: PetscErrorCode ierr;
421: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
422: PetscBool isascii;
425: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
426: if (!isascii) SETERRQ1(((PetscObject)eps)->comm,1,"Viewer type %s not supported for EPS Krylov-Schur",((PetscObject)viewer)->type_name);
427: PetscViewerASCIIPrintf(viewer," Krylov-Schur: %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
428: return(0);
429: }
433: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
434: {
435: PetscErrorCode ierr;
436: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
439: ctx->keep = 0.0;
440: EPSReset_Default(eps);
441: return(0);
442: }
446: PetscErrorCode EPSDestroy_KrylovSchur(EPS eps)
447: {
451: PetscFree(eps->data);
452: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSKrylovSchurSetRestart_C","",PETSC_NULL);
453: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSKrylovSchurGetRestart_C","",PETSC_NULL);
454: return(0);
455: }
457: EXTERN_C_BEGIN
460: PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
461: {
465: PetscNewLog(eps,EPS_KRYLOVSCHUR,&eps->data);
466: eps->ops->setup = EPSSetUp_KrylovSchur;
467: eps->ops->setfromoptions = EPSSetFromOptions_KrylovSchur;
468: eps->ops->destroy = EPSDestroy_KrylovSchur;
469: eps->ops->reset = EPSReset_KrylovSchur;
470: eps->ops->view = EPSView_KrylovSchur;
471: eps->ops->backtransform = EPSBackTransform_Default;
472: eps->ops->computevectors = EPSComputeVectors_Schur;
473: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSKrylovSchurSetRestart_C","EPSKrylovSchurSetRestart_KrylovSchur",EPSKrylovSchurSetRestart_KrylovSchur);
474: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSKrylovSchurGetRestart_C","EPSKrylovSchurGetRestart_KrylovSchur",EPSKrylovSchurGetRestart_KrylovSchur);
475: return(0);
476: }
477: EXTERN_C_END