Actual source code: krylovschur.c
slepc-3.5.0 2014-07-29
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 Anal. App. 23(3):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-2014, Universitat Politecnica de Valencia, Spain
27: This file is part of SLEPc.
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 krylovschur.h
48: PetscErrorCode EPSGetArbitraryValues(EPS eps,PetscScalar *rr,PetscScalar *ri)
49: {
51: PetscInt i,newi,ld,n,l;
52: Vec xr=eps->work[0],xi=eps->work[1];
53: PetscScalar re,im,*Zr,*Zi,*X;
56: DSGetLeadingDimension(eps->ds,&ld);
57: DSGetDimensions(eps->ds,&n,NULL,&l,NULL,NULL);
58: for (i=l;i<n;i++) {
59: re = eps->eigr[i];
60: im = eps->eigi[i];
61: STBackTransform(eps->st,1,&re,&im);
62: newi = i;
63: DSVectors(eps->ds,DS_MAT_X,&newi,NULL);
64: DSGetArray(eps->ds,DS_MAT_X,&X);
65: Zr = X+i*ld;
66: if (newi==i+1) Zi = X+newi*ld;
67: else Zi = NULL;
68: EPSComputeRitzVector(eps,Zr,Zi,eps->V,xr,xi);
69: DSRestoreArray(eps->ds,DS_MAT_X,&X);
70: (*eps->arbitrary)(re,im,xr,xi,rr+i,ri+i,eps->arbitraryctx);
71: }
72: return(0);
73: }
77: PetscErrorCode EPSSetUp_KrylovSchur(EPS eps)
78: {
79: PetscErrorCode ierr;
80: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
81: enum { EPS_KS_DEFAULT,EPS_KS_SYMM,EPS_KS_SLICE,EPS_KS_INDEF } variant;
84: /* spectrum slicing requires special treatment of default values */
85: if (eps->which==EPS_ALL) {
86: EPSSetUp_KrylovSchur_Slice(eps);
87: } else {
88: EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
89: if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
90: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
91: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
92: }
94: if (eps->isgeneralized && eps->ishermitian && !eps->ispositive && eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not implemented for indefinite problems");
95: if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
97: if (!eps->extraction) {
98: EPSSetExtraction(eps,EPS_RITZ);
99: } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC)
100: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
102: if (!ctx->keep) ctx->keep = 0.5;
104: EPSAllocateSolution(eps,1);
105: EPS_SetInnerProduct(eps);
106: if (eps->arbitrary) {
107: EPSSetWorkVecs(eps,2);
108: } else if (eps->ishermitian && !eps->ispositive){
109: EPSSetWorkVecs(eps,1);
110: }
112: /* dispatch solve method */
113: if (eps->ishermitian) {
114: if (eps->which==EPS_ALL) {
115: if (eps->isgeneralized && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing not implemented for indefinite problems");
116: else variant = EPS_KS_SLICE;
117: } else if (eps->isgeneralized && !eps->ispositive) {
118: variant = EPS_KS_INDEF;
119: } else {
120: switch (eps->extraction) {
121: case EPS_RITZ: variant = EPS_KS_SYMM; break;
122: case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
123: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
124: }
125: }
126: } else {
127: switch (eps->extraction) {
128: case EPS_RITZ: variant = EPS_KS_DEFAULT; break;
129: case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
130: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
131: }
132: }
133: switch (variant) {
134: case EPS_KS_DEFAULT:
135: eps->ops->solve = EPSSolve_KrylovSchur_Default;
136: DSSetType(eps->ds,DSNHEP);
137: DSAllocate(eps->ds,eps->ncv+1);
138: break;
139: case EPS_KS_SYMM:
140: eps->ops->solve = EPSSolve_KrylovSchur_Symm;
141: DSSetType(eps->ds,DSHEP);
142: DSSetCompact(eps->ds,PETSC_TRUE);
143: DSSetExtraRow(eps->ds,PETSC_TRUE);
144: DSAllocate(eps->ds,eps->ncv+1);
145: break;
146: case EPS_KS_SLICE:
147: eps->ops->solve = EPSSolve_KrylovSchur_Slice;
148: DSSetType(eps->ds,DSHEP);
149: DSSetCompact(eps->ds,PETSC_TRUE);
150: DSAllocate(eps->ds,ctx->ncv+1);
151: break;
152: case EPS_KS_INDEF:
153: eps->ops->solve = EPSSolve_KrylovSchur_Indefinite;
154: DSSetType(eps->ds,DSGHIEP);
155: DSSetCompact(eps->ds,PETSC_TRUE);
156: DSAllocate(eps->ds,eps->ncv+1);
157: break;
158: default: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error");
159: }
160: return(0);
161: }
165: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
166: {
167: PetscErrorCode ierr;
168: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
169: PetscInt i,j,*pj,k,l,nv,ld;
170: Mat U;
171: PetscScalar *S,*Q,*g;
172: PetscReal beta,gamma=1.0;
173: PetscBool breakdown,harmonic;
176: DSGetLeadingDimension(eps->ds,&ld);
177: harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
178: if (harmonic) { PetscMalloc1(ld,&g); }
179: if (eps->arbitrary) pj = &j;
180: else pj = NULL;
182: /* Get the starting Arnoldi vector */
183: EPSGetStartVector(eps,0,NULL);
184: l = 0;
186: /* Restart loop */
187: while (eps->reason == EPS_CONVERGED_ITERATING) {
188: eps->its++;
190: /* Compute an nv-step Arnoldi factorization */
191: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
192: DSGetArray(eps->ds,DS_MAT_A,&S);
193: EPSBasicArnoldi(eps,PETSC_FALSE,S,ld,eps->nconv+l,&nv,&beta,&breakdown);
194: DSRestoreArray(eps->ds,DS_MAT_A,&S);
195: DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
196: if (l==0) {
197: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
198: } else {
199: DSSetState(eps->ds,DS_STATE_RAW);
200: }
201: BVSetActiveColumns(eps->V,eps->nconv,nv);
203: /* Compute translation of Krylov decomposition if harmonic extraction used */
204: if (harmonic) {
205: DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,g,&gamma);
206: }
208: /* Solve projected problem */
209: DSSolve(eps->ds,eps->eigr,eps->eigi);
210: if (eps->arbitrary) {
211: EPSGetArbitraryValues(eps,eps->rr,eps->ri);
212: j=1;
213: }
214: DSSort(eps->ds,eps->eigr,eps->eigi,eps->rr,eps->ri,pj);
216: /* Check convergence */
217: EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,gamma,&k);
218: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
219: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
221: /* Update l */
222: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
223: else {
224: l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
225: #if !defined(PETSC_USE_COMPLEX)
226: DSGetArray(eps->ds,DS_MAT_A,&S);
227: if (S[k+l+(k+l-1)*ld] != 0.0) {
228: if (k+l<nv-1) l = l+1;
229: else l = l-1;
230: }
231: DSRestoreArray(eps->ds,DS_MAT_A,&S);
232: #endif
233: }
235: if (eps->reason == EPS_CONVERGED_ITERATING) {
236: if (breakdown) {
237: /* Start a new Arnoldi factorization */
238: PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
239: if (k<eps->nev) {
240: EPSGetStartVector(eps,k,&breakdown);
241: if (breakdown) {
242: eps->reason = EPS_DIVERGED_BREAKDOWN;
243: PetscInfo(eps,"Unable to generate more start vectors\n");
244: }
245: }
246: } else {
247: /* Undo translation of Krylov decomposition */
248: if (harmonic) {
249: DSSetDimensions(eps->ds,nv,0,k,l);
250: DSTranslateHarmonic(eps->ds,0.0,beta,PETSC_TRUE,g,&gamma);
251: /* gamma u^ = u - U*g~ */
252: BVMultColumn(eps->V,-1.0,1.0,nv,g);
253: BVScaleColumn(eps->V,nv,1.0/gamma);
254: }
255: /* Prepare the Rayleigh quotient for restart */
256: DSGetArray(eps->ds,DS_MAT_A,&S);
257: DSGetArray(eps->ds,DS_MAT_Q,&Q);
258: for (i=k;i<k+l;i++) {
259: S[k+l+i*ld] = Q[nv-1+i*ld]*beta*gamma;
260: }
261: DSRestoreArray(eps->ds,DS_MAT_A,&S);
262: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
263: }
264: }
265: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
266: DSGetMat(eps->ds,DS_MAT_Q,&U);
267: BVMultInPlace(eps->V,U,eps->nconv,k+l);
268: MatDestroy(&U);
270: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
271: BVCopyColumn(eps->V,nv,k+l);
272: }
273: eps->nconv = k;
274: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
275: }
277: if (harmonic) { PetscFree(g); }
278: /* truncate Schur decomposition and change the state to raw so that
279: PSVectors() computes eigenvectors from scratch */
280: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
281: DSSetState(eps->ds,DS_STATE_RAW);
282: return(0);
283: }
287: static PetscErrorCode EPSKrylovSchurSetRestart_KrylovSchur(EPS eps,PetscReal keep)
288: {
289: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
292: if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
293: else {
294: if (keep<0.1 || keep>0.9) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
295: ctx->keep = keep;
296: }
297: return(0);
298: }
302: /*@
303: EPSKrylovSchurSetRestart - Sets the restart parameter for the Krylov-Schur
304: method, in particular the proportion of basis vectors that must be kept
305: after restart.
307: Logically Collective on EPS
309: Input Parameters:
310: + eps - the eigenproblem solver context
311: - keep - the number of vectors to be kept at restart
313: Options Database Key:
314: . -eps_krylovschur_restart - Sets the restart parameter
316: Notes:
317: Allowed values are in the range [0.1,0.9]. The default is 0.5.
319: Level: advanced
321: .seealso: EPSKrylovSchurGetRestart()
322: @*/
323: PetscErrorCode EPSKrylovSchurSetRestart(EPS eps,PetscReal keep)
324: {
330: PetscTryMethod(eps,"EPSKrylovSchurSetRestart_C",(EPS,PetscReal),(eps,keep));
331: return(0);
332: }
336: static PetscErrorCode EPSKrylovSchurGetRestart_KrylovSchur(EPS eps,PetscReal *keep)
337: {
338: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
341: *keep = ctx->keep;
342: return(0);
343: }
347: /*@
348: EPSKrylovSchurGetRestart - Gets the restart parameter used in the
349: Krylov-Schur method.
351: Not Collective
353: Input Parameter:
354: . eps - the eigenproblem solver context
356: Output Parameter:
357: . keep - the restart parameter
359: Level: advanced
361: .seealso: EPSKrylovSchurSetRestart()
362: @*/
363: PetscErrorCode EPSKrylovSchurGetRestart(EPS eps,PetscReal *keep)
364: {
370: PetscTryMethod(eps,"EPSKrylovSchurGetRestart_C",(EPS,PetscReal*),(eps,keep));
371: return(0);
372: }
376: static PetscErrorCode EPSKrylovSchurSetDimensions_KrylovSchur(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
377: {
378: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
381: if (nev<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
382: ctx->nev = nev;
383: if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
384: ctx->ncv = 0;
385: } else {
386: if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
387: ctx->ncv = ncv;
388: }
389: if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
390: ctx->mpd = 0;
391: } else {
392: if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
393: ctx->mpd = mpd;
394: }
395: eps->state = EPS_STATE_INITIAL;
396: return(0);
397: }
401: /*@
402: EPSKrylovSchurSetDimensions - Sets the dimensions used for each subsolve
403: step in case of doing spectrum slicing for a computational interval.
404: The meaning of the parameters is the same as in EPSSetDimensions().
406: Logically Collective on EPS
408: Input Parameters:
409: + eps - the eigenproblem solver context
410: . nev - number of eigenvalues to compute
411: . ncv - the maximum dimension of the subspace to be used by the subsolve
412: - mpd - the maximum dimension allowed for the projected problem
414: Options Database Key:
415: + -eps_krylovschur_nev <nev> - Sets the number of eigenvalues
416: . -eps_krylovschur_ncv <ncv> - Sets the dimension of the subspace
417: - -eps_krylovschur_mpd <mpd> - Sets the maximum projected dimension
419: Level: advanced
421: .seealso: EPSKrylovSchurGetDimensions(), EPSSetDimensions(), EPSSetInterval()
422: @*/
423: PetscErrorCode EPSKrylovSchurSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
424: {
432: PetscTryMethod(eps,"EPSKrylovSchurSetDimensions_C",(EPS,PetscInt,PetscInt,PetscInt),(eps,nev,ncv,mpd));
433: return(0);
434: }
438: static PetscErrorCode EPSKrylovSchurGetDimensions_KrylovSchur(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
439: {
440: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
443: if (nev) *nev = ctx->nev;
444: if (ncv) *ncv = ctx->ncv;
445: if (mpd) *mpd = ctx->mpd;
446: return(0);
447: }
451: /*@
452: EPSKrylovSchurGetDimensions - Gets the dimensions used for each subsolve
453: step in case of doing spectrum slicing for a computational interval.
455: Not Collective
457: Input Parameter:
458: . eps - the eigenproblem solver context
460: Output Parameters:
461: + nev - number of eigenvalues to compute
462: . ncv - the maximum dimension of the subspace to be used by the subsolve
463: - mpd - the maximum dimension allowed for the projected problem
465: Level: advanced
467: .seealso: EPSKrylovSchurSetDimensions()
468: @*/
469: PetscErrorCode EPSKrylovSchurGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
470: {
475: PetscTryMethod(eps,"EPSKrylovSchurGetDimensions_C",(EPS,PetscInt*,PetscInt*,PetscInt*),(eps,nev,ncv,mpd));
476: return(0);
477: }
479: #define SWAP(a,b,t) {t=a;a=b;b=t;}
483: static PetscErrorCode EPSKrylovSchurGetInertias_KrylovSchur(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
484: {
485: PetscErrorCode ierr;
486: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
487: PetscInt i=0,j,tmpi;
488: PetscReal v,tmpr;
489: EPS_shift s;
492: if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
493: if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
494: if (!ctx->sr->s0) { /* EPSSolve not called yet */
495: *n = 2;
496: } else {
497: *n = 1;
498: s = ctx->sr->s0;
499: while (s) {
500: (*n)++;
501: s = s->neighb[1];
502: }
503: }
504: PetscMalloc1(*n,shifts);
505: PetscMalloc1(*n,inertias);
506: if (!ctx->sr->s0) { /* EPSSolve not called yet */
507: (*shifts)[0] = ctx->sr->int0;
508: (*shifts)[1] = ctx->sr->int1;
509: (*inertias)[0] = ctx->sr->inertia0;
510: (*inertias)[1] = ctx->sr->inertia1;
511: } else {
512: s = ctx->sr->s0;
513: while (s) {
514: (*shifts)[i] = s->value;
515: (*inertias)[i++] = s->inertia;
516: s = s->neighb[1];
517: }
518: (*shifts)[i] = ctx->sr->int1;
519: (*inertias)[i] = ctx->sr->inertia1;
520: }
521: /* sort result */
522: for (i=0;i<*n;i++) {
523: v = (*shifts)[i];
524: for (j=i+1;j<*n;j++) {
525: if (v > (*shifts)[j]) {
526: SWAP((*shifts)[i],(*shifts)[j],tmpr);
527: SWAP((*inertias)[i],(*inertias)[j],tmpi);
528: v = (*shifts)[i];
529: }
530: }
531: }
532: /* remove possible duplicate in last position */
533: if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
534: return(0);
535: }
539: /*@C
540: EPSKrylovSchurGetInertias - Gets the values of the shifts and their
541: corresponding inertias in case of doing spectrum slicing for a
542: computational interval.
544: Not Collective
546: Input Parameter:
547: . eps - the eigenproblem solver context
549: Output Parameters:
550: + n - number of shifts, including the endpoints of the interval
551: . shifts - the values of the shifts used internally in the solver
552: - inertias - the values of the inertia in each shift
554: Notes:
555: If called after EPSSolve(), all shifts used internally by the solver are
556: returned (including both endpoints and any intermediate ones). If called
557: before EPSSolve() and after EPSSetUp() then only the information of the
558: endpoints is available.
560: This function is only available for spectrum slicing runs, see EPSSetInterval().
562: The returned arrays should be freed by the user.
564: Level: advanced
565: @*/
566: PetscErrorCode EPSKrylovSchurGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
567: {
573: PetscTryMethod(eps,"EPSKrylovSchurGetInertias_C",(EPS,PetscInt*,PetscReal**,PetscInt**),(eps,n,shifts,inertias));
574: return(0);
575: }
579: PetscErrorCode EPSSetFromOptions_KrylovSchur(EPS eps)
580: {
582: PetscBool flg;
583: PetscReal keep;
584: PetscInt i,j,k;
587: PetscOptionsHead("EPS Krylov-Schur Options");
588: PetscOptionsReal("-eps_krylovschur_restart","Proportion of vectors kept after restart","EPSKrylovSchurSetRestart",0.5,&keep,&flg);
589: if (flg) {
590: EPSKrylovSchurSetRestart(eps,keep);
591: }
592: i = 1;
593: j = k = PETSC_DECIDE;
594: PetscOptionsInt("-eps_krylovschur_nev","Number of eigenvalues to compute in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",40,&i,NULL);
595: PetscOptionsInt("-eps_krylovschur_ncv","Number of basis vectors in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&j,NULL);
596: PetscOptionsInt("-eps_krylovschur_mpd","Maximum dimension of projected problem in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&k,NULL);
597: EPSKrylovSchurSetDimensions(eps,i,j,k);
598: PetscOptionsTail();
599: return(0);
600: }
604: PetscErrorCode EPSView_KrylovSchur(EPS eps,PetscViewer viewer)
605: {
606: PetscErrorCode ierr;
607: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
608: PetscBool isascii;
611: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
612: if (isascii) {
613: PetscViewerASCIIPrintf(viewer," Krylov-Schur: %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
614: if (eps->which==EPS_ALL) {
615: PetscViewerASCIIPrintf(viewer," Krylov-Schur: doing spectrum slicing with nev=%d, ncv=%d, mpd=%d\n",ctx->nev,ctx->ncv,ctx->mpd);
616: }
617: }
618: return(0);
619: }
623: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
624: {
625: PetscErrorCode ierr;
626: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
627: EPS_shift s;
630: /* Reviewing list of shifts to free memory */
631: if (ctx->sr) {
632: s = ctx->sr->s0;
633: if (s) {
634: while (s->neighb[1]) {
635: s = s->neighb[1];
636: PetscFree(s->neighb[0]);
637: }
638: PetscFree(s);
639: }
640: PetscFree(ctx->sr);
641: }
642: return(0);
643: }
647: PetscErrorCode EPSDestroy_KrylovSchur(EPS eps)
648: {
652: PetscFree(eps->data);
653: PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",NULL);
654: PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",NULL);
655: PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",NULL);
656: PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",NULL);
657: PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",NULL);
658: return(0);
659: }
663: PETSC_EXTERN PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
664: {
665: EPS_KRYLOVSCHUR *ctx;
666: PetscErrorCode ierr;
669: PetscNewLog(eps,&ctx);
670: eps->data = (void*)ctx;
671: ctx->nev = 1;
673: eps->ops->setup = EPSSetUp_KrylovSchur;
674: eps->ops->setfromoptions = EPSSetFromOptions_KrylovSchur;
675: eps->ops->destroy = EPSDestroy_KrylovSchur;
676: eps->ops->reset = EPSReset_KrylovSchur;
677: eps->ops->view = EPSView_KrylovSchur;
678: eps->ops->backtransform = EPSBackTransform_Default;
679: eps->ops->computevectors = EPSComputeVectors_Schur;
680: PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",EPSKrylovSchurSetRestart_KrylovSchur);
681: PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",EPSKrylovSchurGetRestart_KrylovSchur);
682: PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",EPSKrylovSchurSetDimensions_KrylovSchur);
683: PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",EPSKrylovSchurGetDimensions_KrylovSchur);
684: PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",EPSKrylovSchurGetInertias_KrylovSchur);
685: return(0);
686: }