Actual source code: krylovschur.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  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:       eps->ops->computevectors = EPSComputeVectors_Schur;
137:       DSSetType(eps->ds,DSNHEP);
138:       DSAllocate(eps->ds,eps->ncv+1);
139:       break;
140:     case EPS_KS_SYMM:
141:       eps->ops->solve = EPSSolve_KrylovSchur_Symm;
142:       eps->ops->computevectors = EPSComputeVectors_Hermitian;
143:       DSSetType(eps->ds,DSHEP);
144:       DSSetCompact(eps->ds,PETSC_TRUE);
145:       DSSetExtraRow(eps->ds,PETSC_TRUE);
146:       DSAllocate(eps->ds,eps->ncv+1);
147:       break;
148:     case EPS_KS_SLICE:
149:       eps->ops->solve = EPSSolve_KrylovSchur_Slice;
150:       eps->ops->computevectors = NULL;
151:       DSSetType(eps->ds,DSHEP);
152:       DSSetCompact(eps->ds,PETSC_TRUE);
153:       DSAllocate(eps->ds,ctx->ncv+1);
154:       break;
155:     case EPS_KS_INDEF:
156:       eps->ops->solve = EPSSolve_KrylovSchur_Indefinite;
157:       eps->ops->computevectors = EPSComputeVectors_Indefinite;
158:       DSSetType(eps->ds,DSGHIEP);
159:       DSSetCompact(eps->ds,PETSC_TRUE);
160:       DSAllocate(eps->ds,eps->ncv+1);
161:       break;
162:     default: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error");
163:   }
164:   return(0);
165: }

169: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
170: {
171:   PetscErrorCode  ierr;
172:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
173:   PetscInt        i,j,*pj,k,l,nv,ld;
174:   Mat             U;
175:   PetscScalar     *S,*Q,*g;
176:   PetscReal       beta,gamma=1.0;
177:   PetscBool       breakdown,harmonic;

180:   DSGetLeadingDimension(eps->ds,&ld);
181:   harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
182:   if (harmonic) { PetscMalloc1(ld,&g); }
183:   if (eps->arbitrary) pj = &j;
184:   else pj = NULL;

186:   /* Get the starting Arnoldi vector */
187:   EPSGetStartVector(eps,0,NULL);
188:   l = 0;

190:   /* Restart loop */
191:   while (eps->reason == EPS_CONVERGED_ITERATING) {
192:     eps->its++;

194:     /* Compute an nv-step Arnoldi factorization */
195:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
196:     DSGetArray(eps->ds,DS_MAT_A,&S);
197:     EPSBasicArnoldi(eps,PETSC_FALSE,S,ld,eps->nconv+l,&nv,&beta,&breakdown);
198:     DSRestoreArray(eps->ds,DS_MAT_A,&S);
199:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
200:     if (l==0) {
201:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
202:     } else {
203:       DSSetState(eps->ds,DS_STATE_RAW);
204:     }
205:     BVSetActiveColumns(eps->V,eps->nconv,nv);

207:     /* Compute translation of Krylov decomposition if harmonic extraction used */
208:     if (harmonic) {
209:       DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,g,&gamma);
210:     }

212:     /* Solve projected problem */
213:     DSSolve(eps->ds,eps->eigr,eps->eigi);
214:     if (eps->arbitrary) {
215:       EPSGetArbitraryValues(eps,eps->rr,eps->ri);
216:       j=1;
217:     }
218:     DSSort(eps->ds,eps->eigr,eps->eigi,eps->rr,eps->ri,pj);

220:     /* Check convergence */
221:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,gamma,&k);
222:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
223:     if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;

225:     /* Update l */
226:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
227:     else {
228:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
229: #if !defined(PETSC_USE_COMPLEX)
230:       DSGetArray(eps->ds,DS_MAT_A,&S);
231:       if (S[k+l+(k+l-1)*ld] != 0.0) {
232:         if (k+l<nv-1) l = l+1;
233:         else l = l-1;
234:       }
235:       DSRestoreArray(eps->ds,DS_MAT_A,&S);
236: #endif
237:     }

239:     if (eps->reason == EPS_CONVERGED_ITERATING) {
240:       if (breakdown) {
241:         /* Start a new Arnoldi factorization */
242:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
243:         if (k<eps->nev) {
244:           EPSGetStartVector(eps,k,&breakdown);
245:           if (breakdown) {
246:             eps->reason = EPS_DIVERGED_BREAKDOWN;
247:             PetscInfo(eps,"Unable to generate more start vectors\n");
248:           }
249:         }
250:       } else {
251:         /* Undo translation of Krylov decomposition */
252:         if (harmonic) {
253:           DSSetDimensions(eps->ds,nv,0,k,l);
254:           DSTranslateHarmonic(eps->ds,0.0,beta,PETSC_TRUE,g,&gamma);
255:           /* gamma u^ = u - U*g~ */
256:           BVMultColumn(eps->V,-1.0,1.0,nv,g);
257:           BVScaleColumn(eps->V,nv,1.0/gamma);
258:         }
259:         /* Prepare the Rayleigh quotient for restart */
260:         DSGetArray(eps->ds,DS_MAT_A,&S);
261:         DSGetArray(eps->ds,DS_MAT_Q,&Q);
262:         for (i=k;i<k+l;i++) {
263:           S[k+l+i*ld] = Q[nv-1+i*ld]*beta*gamma;
264:         }
265:         DSRestoreArray(eps->ds,DS_MAT_A,&S);
266:         DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
267:       }
268:     }
269:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
270:     DSGetMat(eps->ds,DS_MAT_Q,&U);
271:     BVMultInPlace(eps->V,U,eps->nconv,k+l);
272:     MatDestroy(&U);

274:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
275:       BVCopyColumn(eps->V,nv,k+l);
276:     }
277:     eps->nconv = k;
278:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
279:   }

281:   if (harmonic) { PetscFree(g); }
282:   /* truncate Schur decomposition and change the state to raw so that
283:      PSVectors() computes eigenvectors from scratch */
284:   DSSetDimensions(eps->ds,eps->nconv,0,0,0);
285:   DSSetState(eps->ds,DS_STATE_RAW);
286:   return(0);
287: }

291: static PetscErrorCode EPSKrylovSchurSetRestart_KrylovSchur(EPS eps,PetscReal keep)
292: {
293:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

296:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
297:   else {
298:     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]");
299:     ctx->keep = keep;
300:   }
301:   return(0);
302: }

306: /*@
307:    EPSKrylovSchurSetRestart - Sets the restart parameter for the Krylov-Schur
308:    method, in particular the proportion of basis vectors that must be kept
309:    after restart.

311:    Logically Collective on EPS

313:    Input Parameters:
314: +  eps - the eigenproblem solver context
315: -  keep - the number of vectors to be kept at restart

317:    Options Database Key:
318: .  -eps_krylovschur_restart - Sets the restart parameter

320:    Notes:
321:    Allowed values are in the range [0.1,0.9]. The default is 0.5.

323:    Level: advanced

325: .seealso: EPSKrylovSchurGetRestart()
326: @*/
327: PetscErrorCode EPSKrylovSchurSetRestart(EPS eps,PetscReal keep)
328: {

334:   PetscTryMethod(eps,"EPSKrylovSchurSetRestart_C",(EPS,PetscReal),(eps,keep));
335:   return(0);
336: }

340: static PetscErrorCode EPSKrylovSchurGetRestart_KrylovSchur(EPS eps,PetscReal *keep)
341: {
342:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

345:   *keep = ctx->keep;
346:   return(0);
347: }

351: /*@
352:    EPSKrylovSchurGetRestart - Gets the restart parameter used in the
353:    Krylov-Schur method.

355:    Not Collective

357:    Input Parameter:
358: .  eps - the eigenproblem solver context

360:    Output Parameter:
361: .  keep - the restart parameter

363:    Level: advanced

365: .seealso: EPSKrylovSchurSetRestart()
366: @*/
367: PetscErrorCode EPSKrylovSchurGetRestart(EPS eps,PetscReal *keep)
368: {

374:   PetscTryMethod(eps,"EPSKrylovSchurGetRestart_C",(EPS,PetscReal*),(eps,keep));
375:   return(0);
376: }

380: static PetscErrorCode EPSKrylovSchurSetDimensions_KrylovSchur(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
381: {
382:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

385:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
386:   ctx->nev = nev;
387:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
388:     ctx->ncv = 0;
389:   } else {
390:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
391:     ctx->ncv = ncv;
392:   }
393:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
394:     ctx->mpd = 0;
395:   } else {
396:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
397:     ctx->mpd = mpd;
398:   }
399:   eps->state = EPS_STATE_INITIAL;
400:   return(0);
401: }

405: /*@
406:    EPSKrylovSchurSetDimensions - Sets the dimensions used for each subsolve
407:    step in case of doing spectrum slicing for a computational interval.
408:    The meaning of the parameters is the same as in EPSSetDimensions().

410:    Logically Collective on EPS

412:    Input Parameters:
413: +  eps - the eigenproblem solver context
414: .  nev - number of eigenvalues to compute
415: .  ncv - the maximum dimension of the subspace to be used by the subsolve
416: -  mpd - the maximum dimension allowed for the projected problem

418:    Options Database Key:
419: +  -eps_krylovschur_nev <nev> - Sets the number of eigenvalues
420: .  -eps_krylovschur_ncv <ncv> - Sets the dimension of the subspace
421: -  -eps_krylovschur_mpd <mpd> - Sets the maximum projected dimension

423:    Level: advanced

425: .seealso: EPSKrylovSchurGetDimensions(), EPSSetDimensions(), EPSSetInterval()
426: @*/
427: PetscErrorCode EPSKrylovSchurSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
428: {

436:   PetscTryMethod(eps,"EPSKrylovSchurSetDimensions_C",(EPS,PetscInt,PetscInt,PetscInt),(eps,nev,ncv,mpd));
437:   return(0);
438: }

442: static PetscErrorCode EPSKrylovSchurGetDimensions_KrylovSchur(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
443: {
444:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

447:   if (nev) *nev = ctx->nev;
448:   if (ncv) *ncv = ctx->ncv;
449:   if (mpd) *mpd = ctx->mpd;
450:   return(0);
451: }

455: /*@
456:    EPSKrylovSchurGetDimensions - Gets the dimensions used for each subsolve
457:    step in case of doing spectrum slicing for a computational interval.

459:    Not Collective

461:    Input Parameter:
462: .  eps - the eigenproblem solver context

464:    Output Parameters:
465: +  nev - number of eigenvalues to compute
466: .  ncv - the maximum dimension of the subspace to be used by the subsolve
467: -  mpd - the maximum dimension allowed for the projected problem

469:    Level: advanced

471: .seealso: EPSKrylovSchurSetDimensions()
472: @*/
473: PetscErrorCode EPSKrylovSchurGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
474: {

479:   PetscTryMethod(eps,"EPSKrylovSchurGetDimensions_C",(EPS,PetscInt*,PetscInt*,PetscInt*),(eps,nev,ncv,mpd));
480:   return(0);
481: }

483: #define SWAP(a,b,t) {t=a;a=b;b=t;}

487: static PetscErrorCode EPSKrylovSchurGetInertias_KrylovSchur(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
488: {
489:   PetscErrorCode  ierr;
490:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
491:   PetscInt        i=0,j,tmpi;
492:   PetscReal       v,tmpr;
493:   EPS_shift       s;

496:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
497:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
498:   if (!ctx->sr->s0) {  /* EPSSolve not called yet */
499:     *n = 2;
500:   } else {
501:     *n = 1;
502:     s = ctx->sr->s0;
503:     while (s) {
504:       (*n)++;
505:       s = s->neighb[1];
506:     }
507:   }
508:   PetscMalloc1(*n,shifts);
509:   PetscMalloc1(*n,inertias);
510:   if (!ctx->sr->s0) {  /* EPSSolve not called yet */
511:     (*shifts)[0]   = ctx->sr->int0;
512:     (*shifts)[1]   = ctx->sr->int1;
513:     (*inertias)[0] = ctx->sr->inertia0;
514:     (*inertias)[1] = ctx->sr->inertia1;
515:   } else {
516:     s = ctx->sr->s0;
517:     while (s) {
518:       (*shifts)[i]     = s->value;
519:       (*inertias)[i++] = s->inertia;
520:       s = s->neighb[1];
521:     }
522:     (*shifts)[i]   = ctx->sr->int1;
523:     (*inertias)[i] = ctx->sr->inertia1;
524:   }
525:   /* sort result */
526:   for (i=0;i<*n;i++) {
527:     v = (*shifts)[i];
528:     for (j=i+1;j<*n;j++) {
529:       if (v > (*shifts)[j]) {
530:         SWAP((*shifts)[i],(*shifts)[j],tmpr);
531:         SWAP((*inertias)[i],(*inertias)[j],tmpi);
532:         v = (*shifts)[i];
533:       }
534:     }
535:   }
536:   /* remove possible duplicate in last position */
537:   if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
538:   return(0);
539: }

543: /*@C
544:    EPSKrylovSchurGetInertias - Gets the values of the shifts and their
545:    corresponding inertias in case of doing spectrum slicing for a
546:    computational interval.

548:    Not Collective

550:    Input Parameter:
551: .  eps - the eigenproblem solver context

553:    Output Parameters:
554: +  n        - number of shifts, including the endpoints of the interval
555: .  shifts   - the values of the shifts used internally in the solver
556: -  inertias - the values of the inertia in each shift

558:    Notes:
559:    If called after EPSSolve(), all shifts used internally by the solver are
560:    returned (including both endpoints and any intermediate ones). If called
561:    before EPSSolve() and after EPSSetUp() then only the information of the
562:    endpoints is available.

564:    This function is only available for spectrum slicing runs, see EPSSetInterval().

566:    The returned arrays should be freed by the user.

568:    Level: advanced
569: @*/
570: PetscErrorCode EPSKrylovSchurGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
571: {

577:   PetscTryMethod(eps,"EPSKrylovSchurGetInertias_C",(EPS,PetscInt*,PetscReal**,PetscInt**),(eps,n,shifts,inertias));
578:   return(0);
579: }

583: PetscErrorCode EPSSetFromOptions_KrylovSchur(EPS eps)
584: {
586:   PetscBool      flg;
587:   PetscReal      keep;
588:   PetscInt       i,j,k;

591:   PetscOptionsHead("EPS Krylov-Schur Options");
592:   PetscOptionsReal("-eps_krylovschur_restart","Proportion of vectors kept after restart","EPSKrylovSchurSetRestart",0.5,&keep,&flg);
593:   if (flg) {
594:     EPSKrylovSchurSetRestart(eps,keep);
595:   }
596:   i = 1;
597:   j = k = PETSC_DECIDE;
598:   PetscOptionsInt("-eps_krylovschur_nev","Number of eigenvalues to compute in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",40,&i,NULL);
599:   PetscOptionsInt("-eps_krylovschur_ncv","Number of basis vectors in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&j,NULL);
600:   PetscOptionsInt("-eps_krylovschur_mpd","Maximum dimension of projected problem in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&k,NULL);
601:   EPSKrylovSchurSetDimensions(eps,i,j,k);
602:   PetscOptionsTail();
603:   return(0);
604: }

608: PetscErrorCode EPSView_KrylovSchur(EPS eps,PetscViewer viewer)
609: {
610:   PetscErrorCode  ierr;
611:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
612:   PetscBool       isascii;

615:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
616:   if (isascii) {
617:     PetscViewerASCIIPrintf(viewer,"  Krylov-Schur: %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
618:     if (eps->which==EPS_ALL) {
619:       PetscViewerASCIIPrintf(viewer,"  Krylov-Schur: doing spectrum slicing with nev=%D, ncv=%D, mpd=%D\n",ctx->nev,ctx->ncv,ctx->mpd);
620:     }
621:   }
622:   return(0);
623: }

627: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
628: {
629:   PetscErrorCode  ierr;
630:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
631:   EPS_shift       s;

634:   /* Reviewing list of shifts to free memory */
635:   if (ctx->sr) {
636:     s = ctx->sr->s0;
637:     if (s) {
638:       while (s->neighb[1]) {
639:         s = s->neighb[1];
640:         PetscFree(s->neighb[0]);
641:       }
642:       PetscFree(s);
643:     }
644:     PetscFree(ctx->sr);
645:   }
646:   return(0);
647: }

651: PetscErrorCode EPSDestroy_KrylovSchur(EPS eps)
652: {

656:   PetscFree(eps->data);
657:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",NULL);
658:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",NULL);
659:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",NULL);
660:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",NULL);
661:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",NULL);
662:   return(0);
663: }

667: PETSC_EXTERN PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
668: {
669:   EPS_KRYLOVSCHUR *ctx;
670:   PetscErrorCode  ierr;

673:   PetscNewLog(eps,&ctx);
674:   eps->data = (void*)ctx;
675:   ctx->nev = 1;

677:   eps->ops->setup          = EPSSetUp_KrylovSchur;
678:   eps->ops->setfromoptions = EPSSetFromOptions_KrylovSchur;
679:   eps->ops->destroy        = EPSDestroy_KrylovSchur;
680:   eps->ops->reset          = EPSReset_KrylovSchur;
681:   eps->ops->view           = EPSView_KrylovSchur;
682:   eps->ops->backtransform  = EPSBackTransform_Default;
683:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",EPSKrylovSchurSetRestart_KrylovSchur);
684:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",EPSKrylovSchurGetRestart_KrylovSchur);
685:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",EPSKrylovSchurSetDimensions_KrylovSchur);
686:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",EPSKrylovSchurGetDimensions_KrylovSchur);
687:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",EPSKrylovSchurGetInertias_KrylovSchur);
688:   return(0);
689: }