Actual source code: krylovschur.c

slepc-3.15.2 2021-09-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc eigensolver: "krylovschur"

 13:    Method: Krylov-Schur

 15:    Algorithm:

 17:        Single-vector Krylov-Schur method for non-symmetric problems,
 18:        including harmonic extraction.

 20:    References:

 22:        [1] "Krylov-Schur Methods in SLEPc", SLEPc Technical Report STR-7,
 23:            available at https://slepc.upv.es.

 25:        [2] G.W. Stewart, "A Krylov-Schur Algorithm for Large Eigenproblems",
 26:            SIAM J. Matrix Anal. App. 23(3):601-614, 2001.

 28:        [3] "Practical Implementation of Harmonic Krylov-Schur", SLEPc Technical
 29:             Report STR-9, available at https://slepc.upv.es.
 30: */

 32: #include <slepc/private/epsimpl.h>
 33: #include "krylovschur.h"

 35: PetscErrorCode EPSGetArbitraryValues(EPS eps,PetscScalar *rr,PetscScalar *ri)
 36: {
 38:   PetscInt       i,newi,ld,n,l;
 39:   Vec            xr=eps->work[0],xi=eps->work[1];
 40:   PetscScalar    re,im,*Zr,*Zi,*X;

 43:   DSGetLeadingDimension(eps->ds,&ld);
 44:   DSGetDimensions(eps->ds,&n,NULL,&l,NULL,NULL);
 45:   for (i=l;i<n;i++) {
 46:     re = eps->eigr[i];
 47:     im = eps->eigi[i];
 48:     STBackTransform(eps->st,1,&re,&im);
 49:     newi = i;
 50:     DSVectors(eps->ds,DS_MAT_X,&newi,NULL);
 51:     DSGetArray(eps->ds,DS_MAT_X,&X);
 52:     Zr = X+i*ld;
 53:     if (newi==i+1) Zi = X+newi*ld;
 54:     else Zi = NULL;
 55:     EPSComputeRitzVector(eps,Zr,Zi,eps->V,xr,xi);
 56:     DSRestoreArray(eps->ds,DS_MAT_X,&X);
 57:     (*eps->arbitrary)(re,im,xr,xi,rr+i,ri+i,eps->arbitraryctx);
 58:   }
 59:   return(0);
 60: }

 62: static PetscErrorCode EstimateRange(Mat A,PetscReal *left,PetscReal *right)
 63: {
 65:   PetscInt       nconv;
 66:   PetscScalar    eig0;
 67:   PetscReal      tol=1e-3,errest=tol;
 68:   EPS            eps;

 71:   *left = 0.0; *right = 0.0;
 72:   EPSCreate(PetscObjectComm((PetscObject)A),&eps);
 73:   EPSSetOptionsPrefix(eps,"eps_filter_");
 74:   EPSSetOperators(eps,A,NULL);
 75:   EPSSetProblemType(eps,EPS_HEP);
 76:   EPSSetTolerances(eps,tol,50);
 77:   EPSSetConvergenceTest(eps,EPS_CONV_ABS);
 78:   EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
 79:   EPSSolve(eps);
 80:   EPSGetConverged(eps,&nconv);
 81:   if (nconv>0) {
 82:     EPSGetEigenvalue(eps,0,&eig0,NULL);
 83:     EPSGetErrorEstimate(eps,0,&errest);
 84:   } else eig0 = eps->eigr[0];
 85:   *left = PetscRealPart(eig0)-errest;
 86:   EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
 87:   EPSSolve(eps);
 88:   EPSGetConverged(eps,&nconv);
 89:   if (nconv>0) {
 90:     EPSGetEigenvalue(eps,0,&eig0,NULL);
 91:     EPSGetErrorEstimate(eps,0,&errest);
 92:   } else eig0 = eps->eigr[0];
 93:   *right = PetscRealPart(eig0)+errest;
 94:   EPSDestroy(&eps);
 95:   return(0);
 96: }

 98: static PetscErrorCode EPSSetUp_KrylovSchur_Filter(EPS eps)
 99: {
100:   PetscErrorCode  ierr;
101:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
102:   PetscBool       estimaterange=PETSC_TRUE;
103:   PetscReal       rleft,rright;
104:   Mat             A;

107:   EPSCheckHermitianCondition(eps,PETSC_TRUE," with polynomial filter");
108:   EPSCheckStandardCondition(eps,PETSC_TRUE," with polynomial filter");
109:   if (eps->intb >= PETSC_MAX_REAL && eps->inta <= PETSC_MIN_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The defined computational interval should have at least one of their sides bounded");
110:   EPSCheckUnsupportedCondition(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_EXTRACTION,PETSC_TRUE," with polynomial filter");
111:   if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL*1e-2;  /* use tighter tolerance */
112:   STFilterSetInterval(eps->st,eps->inta,eps->intb);
113:   if (!ctx->estimatedrange) {
114:     STFilterGetRange(eps->st,&rleft,&rright);
115:     estimaterange = (!rleft && !rright)? PETSC_TRUE: PETSC_FALSE;
116:   }
117:   if (estimaterange) { /* user did not set a range */
118:     STGetMatrix(eps->st,0,&A);
119:     EstimateRange(A,&rleft,&rright);
120:     PetscInfo2(eps,"Setting eigenvalue range to [%g,%g]\n",(double)rleft,(double)rright);
121:     STFilterSetRange(eps->st,rleft,rright);
122:     ctx->estimatedrange = PETSC_TRUE;
123:   }
124:   if (eps->ncv==PETSC_DEFAULT && eps->nev==1) eps->nev = 40;  /* user did not provide nev estimation */
125:   EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
126:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
127:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
128:   return(0);
129: }

131: PetscErrorCode EPSSetUp_KrylovSchur(EPS eps)
132: {
133:   PetscErrorCode    ierr;
134:   PetscReal         eta;
135:   PetscBool         isfilt=PETSC_FALSE;
136:   BVOrthogType      otype;
137:   BVOrthogBlockType obtype;
138:   EPS_KRYLOVSCHUR   *ctx = (EPS_KRYLOVSCHUR*)eps->data;
139:   enum { EPS_KS_DEFAULT,EPS_KS_SYMM,EPS_KS_SLICE,EPS_KS_FILTER,EPS_KS_INDEF,EPS_KS_TWOSIDED } variant;

142:   if (eps->which==EPS_ALL) {  /* default values in case of spectrum slicing or polynomial filter  */
143:     PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
144:     if (isfilt) {
145:       EPSSetUp_KrylovSchur_Filter(eps);
146:     } else {
147:       EPSSetUp_KrylovSchur_Slice(eps);
148:     }
149:   } else {
150:     EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
151:     if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
152:     if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
153:     if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
154:   }
155:   if (!ctx->lock && eps->mpd<eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");

157:   EPSCheckDefiniteCondition(eps,eps->arbitrary," with arbitrary selection of eigenpairs");

159:   if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");

161:   if (!ctx->keep) ctx->keep = 0.5;

163:   EPSAllocateSolution(eps,1);
164:   EPS_SetInnerProduct(eps);
165:   if (eps->arbitrary) {
166:     EPSSetWorkVecs(eps,2);
167:   } else if (eps->ishermitian && !eps->ispositive){
168:     EPSSetWorkVecs(eps,1);
169:   }

171:   /* dispatch solve method */
172:   if (eps->ishermitian) {
173:     if (eps->which==EPS_ALL) {
174:       if (eps->isgeneralized && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing not implemented for indefinite problems");
175:       else variant = isfilt? EPS_KS_FILTER: EPS_KS_SLICE;
176:     } else if (eps->isgeneralized && !eps->ispositive) {
177:       variant = EPS_KS_INDEF;
178:     } else {
179:       switch (eps->extraction) {
180:         case EPS_RITZ:     variant = EPS_KS_SYMM; break;
181:         case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
182:         default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
183:       }
184:     }
185:   } else if (eps->twosided) {
186:     variant = EPS_KS_TWOSIDED;
187:   } else {
188:     switch (eps->extraction) {
189:       case EPS_RITZ:     variant = EPS_KS_DEFAULT; break;
190:       case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
191:       default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
192:     }
193:   }
194:   switch (variant) {
195:     case EPS_KS_DEFAULT:
196:       eps->ops->solve = EPSSolve_KrylovSchur_Default;
197:       eps->ops->computevectors = EPSComputeVectors_Schur;
198:       DSSetType(eps->ds,DSNHEP);
199:       DSSetExtraRow(eps->ds,PETSC_TRUE);
200:       DSAllocate(eps->ds,eps->ncv+1);
201:       break;
202:     case EPS_KS_SYMM:
203:     case EPS_KS_FILTER:
204:       eps->ops->solve = EPSSolve_KrylovSchur_Default;
205:       eps->ops->computevectors = EPSComputeVectors_Hermitian;
206:       DSSetType(eps->ds,DSHEP);
207:       DSSetCompact(eps->ds,PETSC_TRUE);
208:       DSSetExtraRow(eps->ds,PETSC_TRUE);
209:       DSAllocate(eps->ds,eps->ncv+1);
210:       break;
211:     case EPS_KS_SLICE:
212:       eps->ops->solve = EPSSolve_KrylovSchur_Slice;
213:       eps->ops->computevectors = EPSComputeVectors_Slice;
214:       break;
215:     case EPS_KS_INDEF:
216:       eps->ops->solve = EPSSolve_KrylovSchur_Indefinite;
217:       eps->ops->computevectors = EPSComputeVectors_Indefinite;
218:       DSSetType(eps->ds,DSGHIEP);
219:       DSSetCompact(eps->ds,PETSC_TRUE);
220:       DSSetExtraRow(eps->ds,PETSC_TRUE);
221:       DSAllocate(eps->ds,eps->ncv+1);
222:       /* force reorthogonalization for pseudo-Lanczos */
223:       BVGetOrthogonalization(eps->V,&otype,NULL,&eta,&obtype);
224:       BVSetOrthogonalization(eps->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
225:       break;
226:     case EPS_KS_TWOSIDED:
227:       eps->ops->solve = EPSSolve_KrylovSchur_TwoSided;
228:       eps->ops->computevectors = EPSComputeVectors_Schur;
229:       DSSetType(eps->ds,DSNHEPTS);
230:       DSAllocate(eps->ds,eps->ncv+1);
231:       DSSetExtraRow(eps->ds,PETSC_TRUE);
232:       break;
233:     default: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error");
234:   }
235:   return(0);
236: }

238: PetscErrorCode EPSSetUpSort_KrylovSchur(EPS eps)
239: {
240:   PetscErrorCode  ierr;
241:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
242:   SlepcSC         sc;
243:   PetscBool       isfilt;

246:   EPSSetUpSort_Default(eps);
247:   if (eps->which==EPS_ALL) {
248:     PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
249:     if (isfilt) {
250:       DSGetSlepcSC(eps->ds,&sc);
251:       sc->rg            = NULL;
252:       sc->comparison    = SlepcCompareLargestReal;
253:       sc->comparisonctx = NULL;
254:       sc->map           = NULL;
255:       sc->mapobj        = NULL;
256:     } else {
257:       if (!ctx->global && ctx->sr->numEigs>0) {
258:         DSGetSlepcSC(eps->ds,&sc);
259:         sc->rg            = NULL;
260:         sc->comparison    = SlepcCompareLargestMagnitude;
261:         sc->comparisonctx = NULL;
262:         sc->map           = NULL;
263:         sc->mapobj        = NULL;
264:       }
265:     }
266:   }
267:   return(0);
268: }

270: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
271: {
272:   PetscErrorCode  ierr;
273:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
274:   PetscInt        j,*pj,k,l,nv,ld,nconv;
275:   Mat             U,Op;
276:   PetscScalar     *S,*g;
277:   PetscReal       beta,gamma=1.0,*a,*b;
278:   PetscBool       breakdown,harmonic,hermitian;

281:   DSGetLeadingDimension(eps->ds,&ld);
282:   harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
283:   hermitian = (eps->ishermitian && !harmonic)?PETSC_TRUE:PETSC_FALSE;
284:   if (harmonic) { PetscMalloc1(ld,&g); }
285:   if (eps->arbitrary) pj = &j;
286:   else pj = NULL;

288:   /* Get the starting Arnoldi vector */
289:   EPSGetStartVector(eps,0,NULL);
290:   l = 0;

292:   /* Restart loop */
293:   while (eps->reason == EPS_CONVERGED_ITERATING) {
294:     eps->its++;

296:     /* Compute an nv-step Arnoldi factorization */
297:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
298:     STGetOperator(eps->st,&Op);
299:     if (hermitian) {
300:       DSGetArrayReal(eps->ds,DS_MAT_T,&a);
301:       b = a + ld;
302:       BVMatLanczos(eps->V,Op,a,b,eps->nconv+l,&nv,&breakdown);
303:       beta = b[nv-1];
304:       DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
305:     } else {
306:       DSGetArray(eps->ds,DS_MAT_A,&S);
307:       BVMatArnoldi(eps->V,Op,S,ld,eps->nconv+l,&nv,&beta,&breakdown);
308:       DSRestoreArray(eps->ds,DS_MAT_A,&S);
309:     }
310:     STRestoreOperator(eps->st,&Op);
311:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
312:     if (l==0) {
313:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
314:     } else {
315:       DSSetState(eps->ds,DS_STATE_RAW);
316:     }
317:     BVSetActiveColumns(eps->V,eps->nconv,nv);

319:     /* Compute translation of Krylov decomposition if harmonic extraction used */
320:     if (harmonic) {
321:       DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,g,&gamma);
322:     }

324:     /* Solve projected problem */
325:     DSSolve(eps->ds,eps->eigr,eps->eigi);
326:     if (eps->arbitrary) {
327:       EPSGetArbitraryValues(eps,eps->rr,eps->ri);
328:       j=1;
329:     }
330:     DSSort(eps->ds,eps->eigr,eps->eigi,eps->rr,eps->ri,pj);
331:     DSUpdateExtraRow(eps->ds);
332:     DSSynchronize(eps->ds,eps->eigr,eps->eigi);

334:     /* Check convergence */
335:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,0.0,gamma,&k);
336:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
337:     nconv = k;

339:     /* Update l */
340:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
341:     else {
342:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
343:       if (!hermitian) { DSGetTruncateSize(eps->ds,k,nv,&l); }
344:     }
345:     if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
346:     if (l) { PetscInfo1(eps,"Preparing to restart keeping l=%D vectors\n",l); }

348:     if (eps->reason == EPS_CONVERGED_ITERATING) {
349:       if (breakdown || k==nv) {
350:         /* Start a new Arnoldi factorization */
351:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
352:         if (k<eps->nev) {
353:           EPSGetStartVector(eps,k,&breakdown);
354:           if (breakdown) {
355:             eps->reason = EPS_DIVERGED_BREAKDOWN;
356:             PetscInfo(eps,"Unable to generate more start vectors\n");
357:           }
358:         }
359:       } else {
360:         /* Undo translation of Krylov decomposition */
361:         if (harmonic) {
362:           DSSetDimensions(eps->ds,nv,0,k,l);
363:           DSTranslateHarmonic(eps->ds,0.0,beta,PETSC_TRUE,g,&gamma);
364:           /* gamma u^ = u - U*g~ */
365:           BVSetActiveColumns(eps->V,0,nv);
366:           BVMultColumn(eps->V,-1.0,1.0,nv,g);
367:           BVScaleColumn(eps->V,nv,1.0/gamma);
368:           BVSetActiveColumns(eps->V,eps->nconv,nv);
369:           DSSetDimensions(eps->ds,nv,0,k,nv);
370:         }
371:         /* Prepare the Rayleigh quotient for restart */
372:         DSTruncate(eps->ds,k+l,PETSC_FALSE);
373:       }
374:     }
375:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
376:     DSGetMat(eps->ds,DS_MAT_Q,&U);
377:     BVMultInPlace(eps->V,U,eps->nconv,k+l);
378:     MatDestroy(&U);

380:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
381:       BVCopyColumn(eps->V,nv,k+l);
382:     }
383:     eps->nconv = k;
384:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
385:   }

387:   if (harmonic) { PetscFree(g); }
388:   DSTruncate(eps->ds,eps->nconv,PETSC_TRUE);
389:   return(0);
390: }

392: static PetscErrorCode EPSKrylovSchurSetRestart_KrylovSchur(EPS eps,PetscReal keep)
393: {
394:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

397:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
398:   else {
399:     if (keep<0.1 || keep>0.9) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument %g must be in the range [0.1,0.9]",keep);
400:     ctx->keep = keep;
401:   }
402:   return(0);
403: }

405: /*@
406:    EPSKrylovSchurSetRestart - Sets the restart parameter for the Krylov-Schur
407:    method, in particular the proportion of basis vectors that must be kept
408:    after restart.

410:    Logically Collective on eps

412:    Input Parameters:
413: +  eps - the eigenproblem solver context
414: -  keep - the number of vectors to be kept at restart

416:    Options Database Key:
417: .  -eps_krylovschur_restart - Sets the restart parameter

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

422:    Level: advanced

424: .seealso: EPSKrylovSchurGetRestart()
425: @*/
426: PetscErrorCode EPSKrylovSchurSetRestart(EPS eps,PetscReal keep)
427: {

433:   PetscTryMethod(eps,"EPSKrylovSchurSetRestart_C",(EPS,PetscReal),(eps,keep));
434:   return(0);
435: }

437: static PetscErrorCode EPSKrylovSchurGetRestart_KrylovSchur(EPS eps,PetscReal *keep)
438: {
439:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

442:   *keep = ctx->keep;
443:   return(0);
444: }

446: /*@
447:    EPSKrylovSchurGetRestart - Gets the restart parameter used in the
448:    Krylov-Schur method.

450:    Not Collective

452:    Input Parameter:
453: .  eps - the eigenproblem solver context

455:    Output Parameter:
456: .  keep - the restart parameter

458:    Level: advanced

460: .seealso: EPSKrylovSchurSetRestart()
461: @*/
462: PetscErrorCode EPSKrylovSchurGetRestart(EPS eps,PetscReal *keep)
463: {

469:   PetscUseMethod(eps,"EPSKrylovSchurGetRestart_C",(EPS,PetscReal*),(eps,keep));
470:   return(0);
471: }

473: static PetscErrorCode EPSKrylovSchurSetLocking_KrylovSchur(EPS eps,PetscBool lock)
474: {
475:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

478:   ctx->lock = lock;
479:   return(0);
480: }

482: /*@
483:    EPSKrylovSchurSetLocking - Choose between locking and non-locking variants of
484:    the Krylov-Schur method.

486:    Logically Collective on eps

488:    Input Parameters:
489: +  eps  - the eigenproblem solver context
490: -  lock - true if the locking variant must be selected

492:    Options Database Key:
493: .  -eps_krylovschur_locking - Sets the locking flag

495:    Notes:
496:    The default is to lock converged eigenpairs when the method restarts.
497:    This behaviour can be changed so that all directions are kept in the
498:    working subspace even if already converged to working accuracy (the
499:    non-locking variant).

501:    Level: advanced

503: .seealso: EPSKrylovSchurGetLocking()
504: @*/
505: PetscErrorCode EPSKrylovSchurSetLocking(EPS eps,PetscBool lock)
506: {

512:   PetscTryMethod(eps,"EPSKrylovSchurSetLocking_C",(EPS,PetscBool),(eps,lock));
513:   return(0);
514: }

516: static PetscErrorCode EPSKrylovSchurGetLocking_KrylovSchur(EPS eps,PetscBool *lock)
517: {
518:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

521:   *lock = ctx->lock;
522:   return(0);
523: }

525: /*@
526:    EPSKrylovSchurGetLocking - Gets the locking flag used in the Krylov-Schur
527:    method.

529:    Not Collective

531:    Input Parameter:
532: .  eps - the eigenproblem solver context

534:    Output Parameter:
535: .  lock - the locking flag

537:    Level: advanced

539: .seealso: EPSKrylovSchurSetLocking()
540: @*/
541: PetscErrorCode EPSKrylovSchurGetLocking(EPS eps,PetscBool *lock)
542: {

548:   PetscUseMethod(eps,"EPSKrylovSchurGetLocking_C",(EPS,PetscBool*),(eps,lock));
549:   return(0);
550: }

552: static PetscErrorCode EPSKrylovSchurSetPartitions_KrylovSchur(EPS eps,PetscInt npart)
553: {
554:   PetscErrorCode  ierr;
555:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
556:   PetscMPIInt     size;

559:   if (ctx->npart!=npart) {
560:     if (ctx->commset) { PetscSubcommDestroy(&ctx->subc); }
561:     EPSDestroy(&ctx->eps);
562:   }
563:   if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
564:     ctx->npart = 1;
565:   } else {
566:     MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size);CHKERRMPI(ierr);
567:     if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
568:     ctx->npart = npart;
569:   }
570:   eps->state = EPS_STATE_INITIAL;
571:   return(0);
572: }

574: /*@
575:    EPSKrylovSchurSetPartitions - Sets the number of partitions for the
576:    case of doing spectrum slicing for a computational interval with the
577:    communicator split in several sub-communicators.

579:    Logically Collective on eps

581:    Input Parameters:
582: +  eps   - the eigenproblem solver context
583: -  npart - number of partitions

585:    Options Database Key:
586: .  -eps_krylovschur_partitions <npart> - Sets the number of partitions

588:    Notes:
589:    By default, npart=1 so all processes in the communicator participate in
590:    the processing of the whole interval. If npart>1 then the interval is
591:    divided into npart subintervals, each of them being processed by a
592:    subset of processes.

594:    The interval is split proportionally unless the separation points are
595:    specified with EPSKrylovSchurSetSubintervals().

597:    Level: advanced

599: .seealso: EPSKrylovSchurSetSubintervals(), EPSSetInterval()
600: @*/
601: PetscErrorCode EPSKrylovSchurSetPartitions(EPS eps,PetscInt npart)
602: {

608:   PetscTryMethod(eps,"EPSKrylovSchurSetPartitions_C",(EPS,PetscInt),(eps,npart));
609:   return(0);
610: }

612: static PetscErrorCode EPSKrylovSchurGetPartitions_KrylovSchur(EPS eps,PetscInt *npart)
613: {
614:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

617:   *npart  = ctx->npart;
618:   return(0);
619: }

621: /*@
622:    EPSKrylovSchurGetPartitions - Gets the number of partitions of the
623:    communicator in case of spectrum slicing.

625:    Not Collective

627:    Input Parameter:
628: .  eps - the eigenproblem solver context

630:    Output Parameter:
631: .  npart - number of partitions

633:    Level: advanced

635: .seealso: EPSKrylovSchurSetPartitions()
636: @*/
637: PetscErrorCode EPSKrylovSchurGetPartitions(EPS eps,PetscInt *npart)
638: {

644:   PetscUseMethod(eps,"EPSKrylovSchurGetPartitions_C",(EPS,PetscInt*),(eps,npart));
645:   return(0);
646: }

648: static PetscErrorCode EPSKrylovSchurSetDetectZeros_KrylovSchur(EPS eps,PetscBool detect)
649: {
650:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

653:   ctx->detect = detect;
654:   eps->state  = EPS_STATE_INITIAL;
655:   return(0);
656: }

658: /*@
659:    EPSKrylovSchurSetDetectZeros - Sets a flag to enforce detection of
660:    zeros during the factorizations throughout the spectrum slicing computation.

662:    Logically Collective on eps

664:    Input Parameters:
665: +  eps    - the eigenproblem solver context
666: -  detect - check for zeros

668:    Options Database Key:
669: .  -eps_krylovschur_detect_zeros - Check for zeros; this takes an optional
670:    bool value (0/1/no/yes/true/false)

672:    Notes:
673:    A zero in the factorization indicates that a shift coincides with an eigenvalue.

675:    This flag is turned off by default, and may be necessary in some cases,
676:    especially when several partitions are being used. This feature currently
677:    requires an external package for factorizations with support for zero
678:    detection, e.g. MUMPS.

680:    Level: advanced

682: .seealso: EPSKrylovSchurSetPartitions(), EPSSetInterval()
683: @*/
684: PetscErrorCode EPSKrylovSchurSetDetectZeros(EPS eps,PetscBool detect)
685: {

691:   PetscTryMethod(eps,"EPSKrylovSchurSetDetectZeros_C",(EPS,PetscBool),(eps,detect));
692:   return(0);
693: }

695: static PetscErrorCode EPSKrylovSchurGetDetectZeros_KrylovSchur(EPS eps,PetscBool *detect)
696: {
697:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

700:   *detect = ctx->detect;
701:   return(0);
702: }

704: /*@
705:    EPSKrylovSchurGetDetectZeros - Gets the flag that enforces zero detection
706:    in spectrum slicing.

708:    Not Collective

710:    Input Parameter:
711: .  eps - the eigenproblem solver context

713:    Output Parameter:
714: .  detect - whether zeros detection is enforced during factorizations

716:    Level: advanced

718: .seealso: EPSKrylovSchurSetDetectZeros()
719: @*/
720: PetscErrorCode EPSKrylovSchurGetDetectZeros(EPS eps,PetscBool *detect)
721: {

727:   PetscUseMethod(eps,"EPSKrylovSchurGetDetectZeros_C",(EPS,PetscBool*),(eps,detect));
728:   return(0);
729: }

731: static PetscErrorCode EPSKrylovSchurSetDimensions_KrylovSchur(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
732: {
733:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

736:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
737:   ctx->nev = nev;
738:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
739:     ctx->ncv = PETSC_DEFAULT;
740:   } else {
741:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
742:     ctx->ncv = ncv;
743:   }
744:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
745:     ctx->mpd = PETSC_DEFAULT;
746:   } else {
747:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
748:     ctx->mpd = mpd;
749:   }
750:   eps->state = EPS_STATE_INITIAL;
751:   return(0);
752: }

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

759:    Logically Collective on eps

761:    Input Parameters:
762: +  eps - the eigenproblem solver context
763: .  nev - number of eigenvalues to compute
764: .  ncv - the maximum dimension of the subspace to be used by the subsolve
765: -  mpd - the maximum dimension allowed for the projected problem

767:    Options Database Key:
768: +  -eps_krylovschur_nev <nev> - Sets the number of eigenvalues
769: .  -eps_krylovschur_ncv <ncv> - Sets the dimension of the subspace
770: -  -eps_krylovschur_mpd <mpd> - Sets the maximum projected dimension

772:    Level: advanced

774: .seealso: EPSKrylovSchurGetDimensions(), EPSSetDimensions(), EPSSetInterval()
775: @*/
776: PetscErrorCode EPSKrylovSchurSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
777: {

785:   PetscTryMethod(eps,"EPSKrylovSchurSetDimensions_C",(EPS,PetscInt,PetscInt,PetscInt),(eps,nev,ncv,mpd));
786:   return(0);
787: }

789: static PetscErrorCode EPSKrylovSchurGetDimensions_KrylovSchur(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
790: {
791:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

794:   if (nev) *nev = ctx->nev;
795:   if (ncv) *ncv = ctx->ncv;
796:   if (mpd) *mpd = ctx->mpd;
797:   return(0);
798: }

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

804:    Not Collective

806:    Input Parameter:
807: .  eps - the eigenproblem solver context

809:    Output Parameters:
810: +  nev - number of eigenvalues to compute
811: .  ncv - the maximum dimension of the subspace to be used by the subsolve
812: -  mpd - the maximum dimension allowed for the projected problem

814:    Level: advanced

816: .seealso: EPSKrylovSchurSetDimensions()
817: @*/
818: PetscErrorCode EPSKrylovSchurGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
819: {

824:   PetscUseMethod(eps,"EPSKrylovSchurGetDimensions_C",(EPS,PetscInt*,PetscInt*,PetscInt*),(eps,nev,ncv,mpd));
825:   return(0);
826: }

828: static PetscErrorCode EPSKrylovSchurSetSubintervals_KrylovSchur(EPS eps,PetscReal* subint)
829: {
830:   PetscErrorCode  ierr;
831:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
832:   PetscInt        i;

835:   if (subint[0]!=eps->inta || subint[ctx->npart]!=eps->intb) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"First and last values must match the endpoints of EPSSetInterval()");
836:   for (i=0;i<ctx->npart;i++) if (subint[i]>subint[i+1]) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Array must contain values in ascending order");
837:   if (ctx->subintervals) { PetscFree(ctx->subintervals); }
838:   PetscMalloc1(ctx->npart+1,&ctx->subintervals);
839:   for (i=0;i<ctx->npart+1;i++) ctx->subintervals[i] = subint[i];
840:   ctx->subintset = PETSC_TRUE;
841:   eps->state = EPS_STATE_INITIAL;
842:   return(0);
843: }

845: /*@C
846:    EPSKrylovSchurSetSubintervals - Sets the points that delimit the
847:    subintervals to be used in spectrum slicing with several partitions.

849:    Logically Collective on eps

851:    Input Parameters:
852: +  eps    - the eigenproblem solver context
853: -  subint - array of real values specifying subintervals

855:    Notes:
856:    This function must be called after EPSKrylovSchurSetPartitions(). For npart
857:    partitions, the argument subint must contain npart+1 real values sorted in
858:    ascending order: subint_0, subint_1, ..., subint_npart, where the first
859:    and last values must coincide with the interval endpoints set with
860:    EPSSetInterval().

862:    The subintervals are then defined by two consecutive points: [subint_0,subint_1],
863:    [subint_1,subint_2], and so on.

865:    Level: advanced

867: .seealso: EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubintervals(), EPSSetInterval()
868: @*/
869: PetscErrorCode EPSKrylovSchurSetSubintervals(EPS eps,PetscReal *subint)
870: {

875:   PetscTryMethod(eps,"EPSKrylovSchurSetSubintervals_C",(EPS,PetscReal*),(eps,subint));
876:   return(0);
877: }

879: static PetscErrorCode EPSKrylovSchurGetSubintervals_KrylovSchur(EPS eps,PetscReal **subint)
880: {
881:   PetscErrorCode  ierr;
882:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
883:   PetscInt        i;

886:   if (!ctx->subintset) {
887:     if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
888:     if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
889:   }
890:   PetscMalloc1(ctx->npart+1,subint);
891:   for (i=0;i<=ctx->npart;i++) (*subint)[i] = ctx->subintervals[i];
892:   return(0);
893: }

895: /*@C
896:    EPSKrylovSchurGetSubintervals - Returns the points that delimit the
897:    subintervals used in spectrum slicing with several partitions.

899:    Logically Collective on eps

901:    Input Parameter:
902: .  eps    - the eigenproblem solver context

904:    Output Parameter:
905: .  subint - array of real values specifying subintervals

907:    Notes:
908:    If the user passed values with EPSKrylovSchurSetSubintervals(), then the
909:    same values are returned. Otherwise, the values computed internally are
910:    obtained.

912:    This function is only available for spectrum slicing runs.

914:    The returned array has length npart+1 (see EPSKrylovSchurGetPartitions())
915:    and should be freed by the user.

917:    Fortran Notes:
918:    The calling sequence from Fortran is
919: .vb
920:    EPSKrylovSchurGetSubintervals(eps,subint,ierr)
921:    double precision subint(npart+1) output
922: .ve

924:    Level: advanced

926: .seealso: EPSKrylovSchurSetSubintervals(), EPSKrylovSchurGetPartitions(), EPSSetInterval()
927: @*/
928: PetscErrorCode EPSKrylovSchurGetSubintervals(EPS eps,PetscReal **subint)
929: {

935:   PetscUseMethod(eps,"EPSKrylovSchurGetSubintervals_C",(EPS,PetscReal**),(eps,subint));
936:   return(0);
937: }

939: static PetscErrorCode EPSKrylovSchurGetInertias_KrylovSchur(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
940: {
941:   PetscErrorCode  ierr;
942:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
943:   PetscInt        i,numsh;
944:   EPS_SR          sr = ctx->sr;

947:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
948:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
949:   switch (eps->state) {
950:   case EPS_STATE_INITIAL:
951:     break;
952:   case EPS_STATE_SETUP:
953:     numsh = ctx->npart+1;
954:     if (n) *n = numsh;
955:     if (shifts) {
956:       PetscMalloc1(numsh,shifts);
957:       (*shifts)[0] = eps->inta;
958:       if (ctx->npart==1) (*shifts)[1] = eps->intb;
959:       else for (i=1;i<numsh;i++) (*shifts)[i] = ctx->subintervals[i];
960:     }
961:     if (inertias) {
962:       PetscMalloc1(numsh,inertias);
963:       (*inertias)[0] = (sr->dir==1)?sr->inertia0:sr->inertia1;
964:       if (ctx->npart==1) (*inertias)[1] = (sr->dir==1)?sr->inertia1:sr->inertia0;
965:       else for (i=1;i<numsh;i++) (*inertias)[i] = (*inertias)[i-1]+ctx->nconv_loc[i-1];
966:     }
967:     break;
968:   case EPS_STATE_SOLVED:
969:   case EPS_STATE_EIGENVECTORS:
970:     numsh = ctx->nshifts;
971:     if (n) *n = numsh;
972:     if (shifts) {
973:       PetscMalloc1(numsh,shifts);
974:       for (i=0;i<numsh;i++) (*shifts)[i] = ctx->shifts[i];
975:     }
976:     if (inertias) {
977:       PetscMalloc1(numsh,inertias);
978:       for (i=0;i<numsh;i++) (*inertias)[i] = ctx->inertias[i];
979:     }
980:     break;
981:   }
982:   return(0);
983: }

985: /*@C
986:    EPSKrylovSchurGetInertias - Gets the values of the shifts and their
987:    corresponding inertias in case of doing spectrum slicing for a
988:    computational interval.

990:    Not Collective

992:    Input Parameter:
993: .  eps - the eigenproblem solver context

995:    Output Parameters:
996: +  n        - number of shifts, including the endpoints of the interval
997: .  shifts   - the values of the shifts used internally in the solver
998: -  inertias - the values of the inertia in each shift

1000:    Notes:
1001:    If called after EPSSolve(), all shifts used internally by the solver are
1002:    returned (including both endpoints and any intermediate ones). If called
1003:    before EPSSolve() and after EPSSetUp() then only the information of the
1004:    endpoints of subintervals is available.

1006:    This function is only available for spectrum slicing runs.

1008:    The returned arrays should be freed by the user. Can pass NULL in any of
1009:    the two arrays if not required.

1011:    Fortran Notes:
1012:    The calling sequence from Fortran is
1013: .vb
1014:    EPSKrylovSchurGetInertias(eps,n,shifts,inertias,ierr)
1015:    integer n
1016:    double precision shifts(*)
1017:    integer inertias(*)
1018: .ve
1019:    The arrays should be at least of length n. The value of n can be determined
1020:    by an initial call
1021: .vb
1022:    EPSKrylovSchurGetInertias(eps,n,PETSC_NULL_REAL,PETSC_NULL_INTEGER,ierr)
1023: .ve

1025:    Level: advanced

1027: .seealso: EPSSetInterval(), EPSKrylovSchurSetSubintervals()
1028: @*/
1029: PetscErrorCode EPSKrylovSchurGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
1030: {

1036:   PetscUseMethod(eps,"EPSKrylovSchurGetInertias_C",(EPS,PetscInt*,PetscReal**,PetscInt**),(eps,n,shifts,inertias));
1037:   return(0);
1038: }

1040: static PetscErrorCode EPSKrylovSchurGetSubcommInfo_KrylovSchur(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
1041: {
1042:   PetscErrorCode  ierr;
1043:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1044:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

1047:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1048:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1049:   if (k) *k = (ctx->npart==1)? 0: ctx->subc->color;
1050:   if (n) *n = sr->numEigs;
1051:   if (v) {
1052:     BVCreateVec(sr->V,v);
1053:   }
1054:   return(0);
1055: }

1057: /*@C
1058:    EPSKrylovSchurGetSubcommInfo - Gets information related to the case of
1059:    doing spectrum slicing for a computational interval with multiple
1060:    communicators.

1062:    Collective on the subcommunicator (if v is given)

1064:    Input Parameter:
1065: .  eps - the eigenproblem solver context

1067:    Output Parameters:
1068: +  k - index of the subinterval for the calling process
1069: .  n - number of eigenvalues found in the k-th subinterval
1070: -  v - a vector owned by processes in the subcommunicator with dimensions
1071:        compatible for locally computed eigenvectors (or NULL)

1073:    Notes:
1074:    This function is only available for spectrum slicing runs.

1076:    The returned Vec should be destroyed by the user.

1078:    Level: advanced

1080: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommPairs()
1081: @*/
1082: PetscErrorCode EPSKrylovSchurGetSubcommInfo(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
1083: {

1088:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommInfo_C",(EPS,PetscInt*,PetscInt*,Vec*),(eps,k,n,v));
1089:   return(0);
1090: }

1092: static PetscErrorCode EPSKrylovSchurGetSubcommPairs_KrylovSchur(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1093: {
1094:   PetscErrorCode  ierr;
1095:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1096:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

1099:   EPSCheckSolved(eps,1);
1100:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1101:   if (i<0 || i>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
1102:   if (eig) *eig = sr->eigr[sr->perm[i]];
1103:   if (v) { BVCopyVec(sr->V,sr->perm[i],v); }
1104:   return(0);
1105: }

1107: /*@
1108:    EPSKrylovSchurGetSubcommPairs - Gets the i-th eigenpair stored
1109:    internally in the subcommunicator to which the calling process belongs.

1111:    Collective on the subcommunicator (if v is given)

1113:    Input Parameter:
1114: +  eps - the eigenproblem solver context
1115: -  i   - index of the solution

1117:    Output Parameters:
1118: +  eig - the eigenvalue
1119: -  v   - the eigenvector

1121:    Notes:
1122:    It is allowed to pass NULL for v if the eigenvector is not required.
1123:    Otherwise, the caller must provide a valid Vec objects, i.e.,
1124:    it must be created by the calling program with EPSKrylovSchurGetSubcommInfo().

1126:    The index i should be a value between 0 and n-1, where n is the number of
1127:    vectors in the local subinterval, see EPSKrylovSchurGetSubcommInfo().

1129:    Level: advanced

1131: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo(), EPSKrylovSchurGetSubcommMats()
1132: @*/
1133: PetscErrorCode EPSKrylovSchurGetSubcommPairs(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1134: {

1140:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommPairs_C",(EPS,PetscInt,PetscScalar*,Vec),(eps,i,eig,v));
1141:   return(0);
1142: }

1144: static PetscErrorCode EPSKrylovSchurGetSubcommMats_KrylovSchur(EPS eps,Mat *A,Mat *B)
1145: {
1146:   PetscErrorCode  ierr;
1147:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

1150:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1151:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1152:   EPSGetOperators(ctx->eps,A,B);
1153:   return(0);
1154: }

1156: /*@C
1157:    EPSKrylovSchurGetSubcommMats - Gets the eigenproblem matrices stored
1158:    internally in the subcommunicator to which the calling process belongs.

1160:    Collective on the subcommunicator

1162:    Input Parameter:
1163: .  eps - the eigenproblem solver context

1165:    Output Parameters:
1166: +  A  - the matrix associated with the eigensystem
1167: -  B  - the second matrix in the case of generalized eigenproblems

1169:    Notes:
1170:    This is the analog of EPSGetOperators(), but returns the matrices distributed
1171:    differently (in the subcommunicator rather than in the parent communicator).

1173:    These matrices should not be modified by the user.

1175:    Level: advanced

1177: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo()
1178: @*/
1179: PetscErrorCode EPSKrylovSchurGetSubcommMats(EPS eps,Mat *A,Mat *B)
1180: {

1185:   PetscTryMethod(eps,"EPSKrylovSchurGetSubcommMats_C",(EPS,Mat*,Mat*),(eps,A,B));
1186:   return(0);
1187: }

1189: static PetscErrorCode EPSKrylovSchurUpdateSubcommMats_KrylovSchur(EPS eps,PetscScalar a,PetscScalar ap,Mat Au,PetscScalar b,PetscScalar bp, Mat Bu,MatStructure str,PetscBool globalup)
1190: {
1191:   PetscErrorCode  ierr;
1192:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data,*subctx;
1193:   Mat             A,B=NULL,Ag,Bg=NULL;
1194:   PetscBool       reuse=PETSC_TRUE;

1197:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1198:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1199:   EPSGetOperators(eps,&Ag,&Bg);
1200:   EPSGetOperators(ctx->eps,&A,&B);

1202:   MatScale(A,a);
1203:   if (Au) {
1204:     MatAXPY(A,ap,Au,str);
1205:   }
1206:   if (B) MatScale(B,b);
1207:   if (Bu) {
1208:     MatAXPY(B,bp,Bu,str);
1209:   }
1210:   EPSSetOperators(ctx->eps,A,B);

1212:   /* Update stored matrix state */
1213:   subctx = (EPS_KRYLOVSCHUR*)ctx->eps->data;
1214:   PetscObjectStateGet((PetscObject)A,&subctx->Astate);
1215:   if (B) { PetscObjectStateGet((PetscObject)B,&subctx->Bstate); }

1217:   /* Update matrices in the parent communicator if requested by user */
1218:   if (globalup) {
1219:     if (ctx->npart>1) {
1220:       if (!ctx->isrow) {
1221:         MatGetOwnershipIS(Ag,&ctx->isrow,&ctx->iscol);
1222:         reuse = PETSC_FALSE;
1223:       }
1224:       if (str==DIFFERENT_NONZERO_PATTERN || str==UNKNOWN_NONZERO_PATTERN) reuse = PETSC_FALSE;
1225:       if (ctx->submata && !reuse) {
1226:         MatDestroyMatrices(1,&ctx->submata);
1227:       }
1228:       MatCreateSubMatrices(A,1,&ctx->isrow,&ctx->iscol,(reuse)?MAT_REUSE_MATRIX:MAT_INITIAL_MATRIX,&ctx->submata);
1229:       MatCreateMPIMatConcatenateSeqMat(((PetscObject)Ag)->comm,ctx->submata[0],PETSC_DECIDE,MAT_REUSE_MATRIX,&Ag);
1230:       if (B) {
1231:         if (ctx->submatb && !reuse) {
1232:           MatDestroyMatrices(1,&ctx->submatb);
1233:         }
1234:         MatCreateSubMatrices(B,1,&ctx->isrow,&ctx->iscol,(reuse)?MAT_REUSE_MATRIX:MAT_INITIAL_MATRIX,&ctx->submatb);
1235:         MatCreateMPIMatConcatenateSeqMat(((PetscObject)Bg)->comm,ctx->submatb[0],PETSC_DECIDE,MAT_REUSE_MATRIX,&Bg);
1236:       }
1237:     }
1238:     PetscObjectStateGet((PetscObject)Ag,&ctx->Astate);
1239:     if (Bg) { PetscObjectStateGet((PetscObject)Bg,&ctx->Bstate); }
1240:   }
1241:   EPSSetOperators(eps,Ag,Bg);
1242:   return(0);
1243: }

1245: /*@
1246:    EPSKrylovSchurUpdateSubcommMats - Update the eigenproblem matrices stored
1247:    internally in the subcommunicator to which the calling process belongs.

1249:    Collective on eps

1251:    Input Parameters:
1252: +  eps - the eigenproblem solver context
1253: .  s   - scalar that multiplies the existing A matrix
1254: .  a   - scalar used in the axpy operation on A
1255: .  Au  - matrix used in the axpy operation on A
1256: .  t   - scalar that multiplies the existing B matrix
1257: .  b   - scalar used in the axpy operation on B
1258: .  Bu  - matrix used in the axpy operation on B
1259: .  str - structure flag
1260: -  globalup - flag indicating if global matrices must be updated

1262:    Notes:
1263:    This function modifies the eigenproblem matrices at the subcommunicator level,
1264:    and optionally updates the global matrices in the parent communicator. The updates
1265:    are expressed as A <-- s*A + a*Au,  B <-- t*B + b*Bu.

1267:    It is possible to update one of the matrices, or both.

1269:    The matrices Au and Bu must be equal in all subcommunicators.

1271:    The str flag is passed to the MatAXPY() operations to perform the updates.

1273:    If globalup is true, communication is carried out to reconstruct the updated
1274:    matrices in the parent communicator. The user must be warned that if global
1275:    matrices are not in sync with subcommunicator matrices, the errors computed
1276:    by EPSComputeError() will be wrong even if the computed solution is correct
1277:    (the synchronization may be done only once at the end).

1279:    Level: advanced

1281: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommMats()
1282: @*/
1283: PetscErrorCode EPSKrylovSchurUpdateSubcommMats(EPS eps,PetscScalar s,PetscScalar a,Mat Au,PetscScalar t,PetscScalar b,Mat Bu,MatStructure str,PetscBool globalup)
1284: {

1297:   PetscTryMethod(eps,"EPSKrylovSchurUpdateSubcommMats_C",(EPS,PetscScalar,PetscScalar,Mat,PetscScalar,PetscScalar,Mat,MatStructure,PetscBool),(eps,s,a,Au,t,b,Bu,str,globalup));
1298:   return(0);
1299: }

1301: PetscErrorCode EPSKrylovSchurGetChildEPS(EPS eps,EPS *child)
1302: {
1303:   PetscErrorCode   ierr;
1304:   EPS_KRYLOVSCHUR  *ctx=(EPS_KRYLOVSCHUR*)eps->data,*ctx_local;
1305:   Mat              A,B=NULL,Ar=NULL,Br=NULL;
1306:   PetscMPIInt      rank;
1307:   PetscObjectState Astate,Bstate=0;
1308:   PetscObjectId    Aid,Bid=0;
1309:   STType           sttype;
1310:   PetscInt         nmat;
1311:   const char       *prefix;

1314:   EPSGetOperators(eps,&A,&B);
1315:   if (ctx->npart==1) {
1316:     if (!ctx->eps) {EPSCreate(((PetscObject)eps)->comm,&ctx->eps);}
1317:     EPSGetOptionsPrefix(eps,&prefix);
1318:     EPSSetOptionsPrefix(ctx->eps,prefix);
1319:     EPSSetOperators(ctx->eps,A,B);
1320:   } else {
1321:     PetscObjectStateGet((PetscObject)A,&Astate);
1322:     PetscObjectGetId((PetscObject)A,&Aid);
1323:     if (B) {
1324:       PetscObjectStateGet((PetscObject)B,&Bstate);
1325:       PetscObjectGetId((PetscObject)B,&Bid);
1326:     }
1327:     if (!ctx->subc) {
1328:       /* Create context for subcommunicators */
1329:       PetscSubcommCreate(PetscObjectComm((PetscObject)eps),&ctx->subc);
1330:       PetscSubcommSetNumber(ctx->subc,ctx->npart);
1331:       PetscSubcommSetType(ctx->subc,PETSC_SUBCOMM_CONTIGUOUS);
1332:       PetscLogObjectMemory((PetscObject)eps,sizeof(PetscSubcomm));

1334:       /* Duplicate matrices */
1335:       MatCreateRedundantMatrix(A,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Ar);
1336:       ctx->Astate = Astate;
1337:       ctx->Aid = Aid;
1338:       MatPropagateSymmetryOptions(A,Ar);
1339:       if (B) {
1340:         MatCreateRedundantMatrix(B,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Br);
1341:         ctx->Bstate = Bstate;
1342:         ctx->Bid = Bid;
1343:         MatPropagateSymmetryOptions(B,Br);
1344:       }
1345:     } else {
1346:       if (ctx->Astate != Astate || (B && ctx->Bstate != Bstate) || ctx->Aid != Aid || (B && ctx->Bid != Bid)) {
1347:         STGetNumMatrices(ctx->eps->st,&nmat);
1348:         if (nmat) {EPSGetOperators(ctx->eps,&Ar,&Br);}
1349:         MatCreateRedundantMatrix(A,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Ar);
1350:         ctx->Astate = Astate;
1351:         ctx->Aid = Aid;
1352:         MatPropagateSymmetryOptions(A,Ar);
1353:         if (B) {
1354:           MatCreateRedundantMatrix(B,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Br);
1355:           ctx->Bstate = Bstate;
1356:           ctx->Bid = Bid;
1357:           MatPropagateSymmetryOptions(B,Br);
1358:         }
1359:         EPSSetOperators(ctx->eps,Ar,Br);
1360:         MatDestroy(&Ar);
1361:         MatDestroy(&Br);
1362:       }
1363:     }

1365:     /* Create auxiliary EPS */
1366:     if (!ctx->eps) {
1367:       EPSCreate(PetscSubcommChild(ctx->subc),&ctx->eps);
1368:       EPSGetOptionsPrefix(eps,&prefix);
1369:       EPSSetOptionsPrefix(ctx->eps,prefix);
1370:       EPSSetOperators(ctx->eps,Ar,Br);
1371:       MatDestroy(&Ar);
1372:       MatDestroy(&Br);
1373:     }
1374:     /* Create subcommunicator grouping processes with same rank */
1375:     if (ctx->commset) { MPI_Comm_free(&ctx->commrank);CHKERRMPI(ierr); }
1376:     MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);CHKERRMPI(ierr);
1377:     MPI_Comm_split(((PetscObject)eps)->comm,rank,ctx->subc->color,&ctx->commrank);CHKERRMPI(ierr);
1378:     ctx->commset = PETSC_TRUE;
1379:   }
1380:   EPSSetType(ctx->eps,((PetscObject)eps)->type_name);
1381:   STGetType(eps->st,&sttype);
1382:   STSetType(ctx->eps->st,sttype);

1384:   ctx_local = (EPS_KRYLOVSCHUR*)ctx->eps->data;
1385:   ctx_local->npart = ctx->npart;
1386:   ctx_local->global = PETSC_FALSE;
1387:   ctx_local->eps = eps;
1388:   ctx_local->subc = ctx->subc;
1389:   ctx_local->commrank = ctx->commrank;
1390:   *child = ctx->eps;
1391:   return(0);
1392: }

1394: static PetscErrorCode EPSKrylovSchurGetKSP_KrylovSchur(EPS eps,KSP *ksp)
1395: {
1396:   PetscErrorCode  ierr;
1397:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1398:   ST              st;
1399:   PetscBool       isfilt;

1402:   PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
1403:   if (eps->which!=EPS_ALL || isfilt) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations with spectrum slicing");
1404:   EPSKrylovSchurGetChildEPS(eps,&ctx->eps);
1405:   EPSGetST(ctx->eps,&st);
1406:   STGetOperator(st,NULL);
1407:   STGetKSP(st,ksp);
1408:   return(0);
1409: }

1411: /*@
1412:    EPSKrylovSchurGetKSP - Retrieve the linear solver object associated with the
1413:    internal EPS object in case of doing spectrum slicing for a computational interval.

1415:    Collective on eps

1417:    Input Parameter:
1418: .  eps - the eigenproblem solver context

1420:    Output Parameters:
1421: -  ksp - the internal KSP object

1423:    Notes:
1424:    When invoked to compute all eigenvalues in an interval with spectrum
1425:    slicing, EPSKRYLOVSCHUR creates another EPS object internally that is
1426:    used to compute eigenvalues by chunks near selected shifts. This function
1427:    allows access to the KSP object associated to this internal EPS object.

1429:    This function is only available for spectrum slicing runs. In case of
1430:    having more than one partition, the returned KSP will be different
1431:    in MPI processes belonging to different partitions. Hence, if required,
1432:    EPSKrylovSchurSetPartitions() must be called BEFORE this function.

1434:    Level: advanced

1436: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions()
1437: @*/
1438: PetscErrorCode EPSKrylovSchurGetKSP(EPS eps,KSP *ksp)
1439: {

1444:   PetscUseMethod(eps,"EPSKrylovSchurGetKSP_C",(EPS,KSP*),(eps,ksp));
1445:   return(0);
1446: }

1448: PetscErrorCode EPSSetFromOptions_KrylovSchur(PetscOptionItems *PetscOptionsObject,EPS eps)
1449: {
1450:   PetscErrorCode  ierr;
1451:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1452:   PetscBool       flg,lock,b,f1,f2,f3,isfilt;
1453:   PetscReal       keep;
1454:   PetscInt        i,j,k;
1455:   KSP             ksp;

1458:   PetscOptionsHead(PetscOptionsObject,"EPS Krylov-Schur Options");

1460:     PetscOptionsReal("-eps_krylovschur_restart","Proportion of vectors kept after restart","EPSKrylovSchurSetRestart",0.5,&keep,&flg);
1461:     if (flg) { EPSKrylovSchurSetRestart(eps,keep); }

1463:     PetscOptionsBool("-eps_krylovschur_locking","Choose between locking and non-locking variants","EPSKrylovSchurSetLocking",PETSC_TRUE,&lock,&flg);
1464:     if (flg) { EPSKrylovSchurSetLocking(eps,lock); }

1466:     i = ctx->npart;
1467:     PetscOptionsInt("-eps_krylovschur_partitions","Number of partitions of the communicator for spectrum slicing","EPSKrylovSchurSetPartitions",ctx->npart,&i,&flg);
1468:     if (flg) { EPSKrylovSchurSetPartitions(eps,i); }

1470:     b = ctx->detect;
1471:     PetscOptionsBool("-eps_krylovschur_detect_zeros","Check zeros during factorizations at subinterval boundaries","EPSKrylovSchurSetDetectZeros",ctx->detect,&b,&flg);
1472:     if (flg) { EPSKrylovSchurSetDetectZeros(eps,b); }

1474:     i = 1;
1475:     j = k = PETSC_DECIDE;
1476:     PetscOptionsInt("-eps_krylovschur_nev","Number of eigenvalues to compute in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",40,&i,&f1);
1477:     PetscOptionsInt("-eps_krylovschur_ncv","Number of basis vectors in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&j,&f2);
1478:     PetscOptionsInt("-eps_krylovschur_mpd","Maximum dimension of projected problem in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&k,&f3);
1479:     if (f1 || f2 || f3) { EPSKrylovSchurSetDimensions(eps,i,j,k); }

1481:   PetscOptionsTail();

1483:   /* set options of child KSP in spectrum slicing */
1484:   if (eps->which==EPS_ALL) {
1485:     if (!eps->st) { EPSGetST(eps,&eps->st); }
1486:     EPSSetDefaultST(eps);
1487:     STSetFromOptions(eps->st);  /* need to advance this to check ST type */
1488:     PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
1489:     if (!isfilt) {
1490:       EPSKrylovSchurGetKSP_KrylovSchur(eps,&ksp);
1491:       KSPSetFromOptions(ksp);
1492:     }
1493:   }
1494:   return(0);
1495: }

1497: PetscErrorCode EPSView_KrylovSchur(EPS eps,PetscViewer viewer)
1498: {
1499:   PetscErrorCode  ierr;
1500:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1501:   PetscBool       isascii,isfilt;
1502:   KSP             ksp;
1503:   PetscViewer     sviewer;

1506:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1507:   if (isascii) {
1508:     PetscViewerASCIIPrintf(viewer,"  %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
1509:     PetscViewerASCIIPrintf(viewer,"  using the %slocking variant\n",ctx->lock?"":"non-");
1510:     if (eps->which==EPS_ALL) {
1511:       PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
1512:       if (isfilt) {
1513:         PetscViewerASCIIPrintf(viewer,"  using filtering to extract all eigenvalues in an interval\n");
1514:       } else {
1515:         PetscViewerASCIIPrintf(viewer,"  doing spectrum slicing with nev=%D, ncv=%D, mpd=%D\n",ctx->nev,ctx->ncv,ctx->mpd);
1516:         if (ctx->npart>1) {
1517:           PetscViewerASCIIPrintf(viewer,"  multi-communicator spectrum slicing with %D partitions\n",ctx->npart);
1518:           if (ctx->detect) { PetscViewerASCIIPrintf(viewer,"  detecting zeros when factorizing at subinterval boundaries\n"); }
1519:         }
1520:         /* view child KSP */
1521:         EPSKrylovSchurGetKSP_KrylovSchur(eps,&ksp);
1522:         PetscViewerASCIIPushTab(viewer);
1523:         if (ctx->npart>1 && ctx->subc) {
1524:           PetscViewerGetSubViewer(viewer,ctx->subc->child,&sviewer);
1525:           if (!ctx->subc->color) {
1526:             KSPView(ksp,sviewer);
1527:           }
1528:           PetscViewerFlush(sviewer);
1529:           PetscViewerRestoreSubViewer(viewer,ctx->subc->child,&sviewer);
1530:           PetscViewerFlush(viewer);
1531:           /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1532:           PetscViewerASCIIPopSynchronized(viewer);
1533:         } else {
1534:           KSPView(ksp,viewer);
1535:         }
1536:         PetscViewerASCIIPopTab(viewer);
1537:       }
1538:     }
1539:   }
1540:   return(0);
1541: }

1543: PetscErrorCode EPSDestroy_KrylovSchur(EPS eps)
1544: {
1546:   PetscBool      isfilt;

1549:   PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
1550:   if (eps->which==EPS_ALL && !isfilt) {
1551:     EPSDestroy_KrylovSchur_Slice(eps);
1552:   }
1553:   PetscFree(eps->data);
1554:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",NULL);
1555:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",NULL);
1556:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetLocking_C",NULL);
1557:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetLocking_C",NULL);
1558:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetPartitions_C",NULL);
1559:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetPartitions_C",NULL);
1560:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDetectZeros_C",NULL);
1561:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDetectZeros_C",NULL);
1562:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",NULL);
1563:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",NULL);
1564:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetSubintervals_C",NULL);
1565:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubintervals_C",NULL);
1566:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",NULL);
1567:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommInfo_C",NULL);
1568:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommPairs_C",NULL);
1569:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommMats_C",NULL);
1570:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurUpdateSubcommMats_C",NULL);
1571:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetKSP_C",NULL);
1572:   return(0);
1573: }

1575: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
1576: {
1578:   PetscBool      isfilt;

1581:   PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
1582:   if (eps->which==EPS_ALL && !isfilt) {
1583:     EPSReset_KrylovSchur_Slice(eps);
1584:   }
1585:   return(0);
1586: }

1588: PetscErrorCode EPSSetDefaultST_KrylovSchur(EPS eps)
1589: {

1593:   if (eps->which==EPS_ALL) {
1594:     if (!((PetscObject)eps->st)->type_name) {
1595:       STSetType(eps->st,STSINVERT);
1596:     }
1597:   }
1598:   return(0);
1599: }

1601: SLEPC_EXTERN PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
1602: {
1603:   EPS_KRYLOVSCHUR *ctx;
1604:   PetscErrorCode  ierr;

1607:   PetscNewLog(eps,&ctx);
1608:   eps->data   = (void*)ctx;
1609:   ctx->lock   = PETSC_TRUE;
1610:   ctx->nev    = 1;
1611:   ctx->ncv    = PETSC_DEFAULT;
1612:   ctx->mpd    = PETSC_DEFAULT;
1613:   ctx->npart  = 1;
1614:   ctx->detect = PETSC_FALSE;
1615:   ctx->global = PETSC_TRUE;

1617:   eps->useds = PETSC_TRUE;

1619:   /* solve and computevectors determined at setup */
1620:   eps->ops->setup          = EPSSetUp_KrylovSchur;
1621:   eps->ops->setupsort      = EPSSetUpSort_KrylovSchur;
1622:   eps->ops->setfromoptions = EPSSetFromOptions_KrylovSchur;
1623:   eps->ops->destroy        = EPSDestroy_KrylovSchur;
1624:   eps->ops->reset          = EPSReset_KrylovSchur;
1625:   eps->ops->view           = EPSView_KrylovSchur;
1626:   eps->ops->backtransform  = EPSBackTransform_Default;
1627:   eps->ops->setdefaultst   = EPSSetDefaultST_KrylovSchur;

1629:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",EPSKrylovSchurSetRestart_KrylovSchur);
1630:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",EPSKrylovSchurGetRestart_KrylovSchur);
1631:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetLocking_C",EPSKrylovSchurSetLocking_KrylovSchur);
1632:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetLocking_C",EPSKrylovSchurGetLocking_KrylovSchur);
1633:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetPartitions_C",EPSKrylovSchurSetPartitions_KrylovSchur);
1634:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetPartitions_C",EPSKrylovSchurGetPartitions_KrylovSchur);
1635:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDetectZeros_C",EPSKrylovSchurSetDetectZeros_KrylovSchur);
1636:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDetectZeros_C",EPSKrylovSchurGetDetectZeros_KrylovSchur);
1637:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",EPSKrylovSchurSetDimensions_KrylovSchur);
1638:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",EPSKrylovSchurGetDimensions_KrylovSchur);
1639:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetSubintervals_C",EPSKrylovSchurSetSubintervals_KrylovSchur);
1640:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubintervals_C",EPSKrylovSchurGetSubintervals_KrylovSchur);
1641:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",EPSKrylovSchurGetInertias_KrylovSchur);
1642:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommInfo_C",EPSKrylovSchurGetSubcommInfo_KrylovSchur);
1643:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommPairs_C",EPSKrylovSchurGetSubcommPairs_KrylovSchur);
1644:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommMats_C",EPSKrylovSchurGetSubcommMats_KrylovSchur);
1645:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurUpdateSubcommMats_C",EPSKrylovSchurUpdateSubcommMats_KrylovSchur);
1646:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetKSP_C",EPSKrylovSchurGetKSP_KrylovSchur);
1647:   return(0);
1648: }