Actual source code: krylovschur.c

slepc-3.15.0 2021-03-31
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:   EPS            eps;

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

 94: static PetscErrorCode EPSSetUp_KrylovSchur_Filter(EPS eps)
 95: {
 96:   PetscErrorCode  ierr;
 97:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
 98:   PetscBool       estimaterange=PETSC_TRUE;
 99:   PetscReal       rleft,rright;
100:   Mat             A;

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

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

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

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

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

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

159:   EPSAllocateSolution(eps,1);
160:   EPS_SetInnerProduct(eps);
161:   if (eps->arbitrary) {
162:     EPSSetWorkVecs(eps,2);
163:   } else if (eps->ishermitian && !eps->ispositive){
164:     EPSSetWorkVecs(eps,1);
165:   }

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

234: PetscErrorCode EPSSetUpSort_KrylovSchur(EPS eps)
235: {
236:   PetscErrorCode  ierr;
237:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
238:   SlepcSC         sc;
239:   PetscBool       isfilt;

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

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

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

284:   /* Get the starting Arnoldi vector */
285:   EPSGetStartVector(eps,0,NULL);
286:   l = 0;

288:   /* Restart loop */
289:   while (eps->reason == EPS_CONVERGED_ITERATING) {
290:     eps->its++;

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

315:     /* Compute translation of Krylov decomposition if harmonic extraction used */
316:     if (harmonic) {
317:       DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,g,&gamma);
318:     }

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

330:     /* Check convergence */
331:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,0.0,gamma,&k);
332:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
333:     nconv = k;

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

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

376:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
377:       BVCopyColumn(eps->V,nv,k+l);
378:     }
379:     eps->nconv = k;
380:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
381:   }

383:   if (harmonic) { PetscFree(g); }
384:   DSTruncate(eps->ds,eps->nconv,PETSC_TRUE);
385:   return(0);
386: }

388: static PetscErrorCode EPSKrylovSchurSetRestart_KrylovSchur(EPS eps,PetscReal keep)
389: {
390:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

393:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
394:   else {
395:     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);
396:     ctx->keep = keep;
397:   }
398:   return(0);
399: }

401: /*@
402:    EPSKrylovSchurSetRestart - Sets the restart parameter for the Krylov-Schur
403:    method, in particular the proportion of basis vectors that must be kept
404:    after restart.

406:    Logically Collective on eps

408:    Input Parameters:
409: +  eps - the eigenproblem solver context
410: -  keep - the number of vectors to be kept at restart

412:    Options Database Key:
413: .  -eps_krylovschur_restart - Sets the restart parameter

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

418:    Level: advanced

420: .seealso: EPSKrylovSchurGetRestart()
421: @*/
422: PetscErrorCode EPSKrylovSchurSetRestart(EPS eps,PetscReal keep)
423: {

429:   PetscTryMethod(eps,"EPSKrylovSchurSetRestart_C",(EPS,PetscReal),(eps,keep));
430:   return(0);
431: }

433: static PetscErrorCode EPSKrylovSchurGetRestart_KrylovSchur(EPS eps,PetscReal *keep)
434: {
435:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

438:   *keep = ctx->keep;
439:   return(0);
440: }

442: /*@
443:    EPSKrylovSchurGetRestart - Gets the restart parameter used in the
444:    Krylov-Schur method.

446:    Not Collective

448:    Input Parameter:
449: .  eps - the eigenproblem solver context

451:    Output Parameter:
452: .  keep - the restart parameter

454:    Level: advanced

456: .seealso: EPSKrylovSchurSetRestart()
457: @*/
458: PetscErrorCode EPSKrylovSchurGetRestart(EPS eps,PetscReal *keep)
459: {

465:   PetscUseMethod(eps,"EPSKrylovSchurGetRestart_C",(EPS,PetscReal*),(eps,keep));
466:   return(0);
467: }

469: static PetscErrorCode EPSKrylovSchurSetLocking_KrylovSchur(EPS eps,PetscBool lock)
470: {
471:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

474:   ctx->lock = lock;
475:   return(0);
476: }

478: /*@
479:    EPSKrylovSchurSetLocking - Choose between locking and non-locking variants of
480:    the Krylov-Schur method.

482:    Logically Collective on eps

484:    Input Parameters:
485: +  eps  - the eigenproblem solver context
486: -  lock - true if the locking variant must be selected

488:    Options Database Key:
489: .  -eps_krylovschur_locking - Sets the locking flag

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

497:    Level: advanced

499: .seealso: EPSKrylovSchurGetLocking()
500: @*/
501: PetscErrorCode EPSKrylovSchurSetLocking(EPS eps,PetscBool lock)
502: {

508:   PetscTryMethod(eps,"EPSKrylovSchurSetLocking_C",(EPS,PetscBool),(eps,lock));
509:   return(0);
510: }

512: static PetscErrorCode EPSKrylovSchurGetLocking_KrylovSchur(EPS eps,PetscBool *lock)
513: {
514:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

517:   *lock = ctx->lock;
518:   return(0);
519: }

521: /*@
522:    EPSKrylovSchurGetLocking - Gets the locking flag used in the Krylov-Schur
523:    method.

525:    Not Collective

527:    Input Parameter:
528: .  eps - the eigenproblem solver context

530:    Output Parameter:
531: .  lock - the locking flag

533:    Level: advanced

535: .seealso: EPSKrylovSchurSetLocking()
536: @*/
537: PetscErrorCode EPSKrylovSchurGetLocking(EPS eps,PetscBool *lock)
538: {

544:   PetscUseMethod(eps,"EPSKrylovSchurGetLocking_C",(EPS,PetscBool*),(eps,lock));
545:   return(0);
546: }

548: static PetscErrorCode EPSKrylovSchurSetPartitions_KrylovSchur(EPS eps,PetscInt npart)
549: {
550:   PetscErrorCode  ierr;
551:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
552:   PetscMPIInt     size;

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

570: /*@
571:    EPSKrylovSchurSetPartitions - Sets the number of partitions for the
572:    case of doing spectrum slicing for a computational interval with the
573:    communicator split in several sub-communicators.

575:    Logically Collective on eps

577:    Input Parameters:
578: +  eps   - the eigenproblem solver context
579: -  npart - number of partitions

581:    Options Database Key:
582: .  -eps_krylovschur_partitions <npart> - Sets the number of partitions

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

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

593:    Level: advanced

595: .seealso: EPSKrylovSchurSetSubintervals(), EPSSetInterval()
596: @*/
597: PetscErrorCode EPSKrylovSchurSetPartitions(EPS eps,PetscInt npart)
598: {

604:   PetscTryMethod(eps,"EPSKrylovSchurSetPartitions_C",(EPS,PetscInt),(eps,npart));
605:   return(0);
606: }

608: static PetscErrorCode EPSKrylovSchurGetPartitions_KrylovSchur(EPS eps,PetscInt *npart)
609: {
610:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

613:   *npart  = ctx->npart;
614:   return(0);
615: }

617: /*@
618:    EPSKrylovSchurGetPartitions - Gets the number of partitions of the
619:    communicator in case of spectrum slicing.

621:    Not Collective

623:    Input Parameter:
624: .  eps - the eigenproblem solver context

626:    Output Parameter:
627: .  npart - number of partitions

629:    Level: advanced

631: .seealso: EPSKrylovSchurSetPartitions()
632: @*/
633: PetscErrorCode EPSKrylovSchurGetPartitions(EPS eps,PetscInt *npart)
634: {

640:   PetscUseMethod(eps,"EPSKrylovSchurGetPartitions_C",(EPS,PetscInt*),(eps,npart));
641:   return(0);
642: }

644: static PetscErrorCode EPSKrylovSchurSetDetectZeros_KrylovSchur(EPS eps,PetscBool detect)
645: {
646:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

649:   ctx->detect = detect;
650:   eps->state  = EPS_STATE_INITIAL;
651:   return(0);
652: }

654: /*@
655:    EPSKrylovSchurSetDetectZeros - Sets a flag to enforce detection of
656:    zeros during the factorizations throughout the spectrum slicing computation.

658:    Logically Collective on eps

660:    Input Parameters:
661: +  eps    - the eigenproblem solver context
662: -  detect - check for zeros

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

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

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

676:    Level: advanced

678: .seealso: EPSKrylovSchurSetPartitions(), EPSSetInterval()
679: @*/
680: PetscErrorCode EPSKrylovSchurSetDetectZeros(EPS eps,PetscBool detect)
681: {

687:   PetscTryMethod(eps,"EPSKrylovSchurSetDetectZeros_C",(EPS,PetscBool),(eps,detect));
688:   return(0);
689: }

691: static PetscErrorCode EPSKrylovSchurGetDetectZeros_KrylovSchur(EPS eps,PetscBool *detect)
692: {
693:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

696:   *detect = ctx->detect;
697:   return(0);
698: }

700: /*@
701:    EPSKrylovSchurGetDetectZeros - Gets the flag that enforces zero detection
702:    in spectrum slicing.

704:    Not Collective

706:    Input Parameter:
707: .  eps - the eigenproblem solver context

709:    Output Parameter:
710: .  detect - whether zeros detection is enforced during factorizations

712:    Level: advanced

714: .seealso: EPSKrylovSchurSetDetectZeros()
715: @*/
716: PetscErrorCode EPSKrylovSchurGetDetectZeros(EPS eps,PetscBool *detect)
717: {

723:   PetscUseMethod(eps,"EPSKrylovSchurGetDetectZeros_C",(EPS,PetscBool*),(eps,detect));
724:   return(0);
725: }

727: static PetscErrorCode EPSKrylovSchurSetDimensions_KrylovSchur(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
728: {
729:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

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

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

755:    Logically Collective on eps

757:    Input Parameters:
758: +  eps - the eigenproblem solver context
759: .  nev - number of eigenvalues to compute
760: .  ncv - the maximum dimension of the subspace to be used by the subsolve
761: -  mpd - the maximum dimension allowed for the projected problem

763:    Options Database Key:
764: +  -eps_krylovschur_nev <nev> - Sets the number of eigenvalues
765: .  -eps_krylovschur_ncv <ncv> - Sets the dimension of the subspace
766: -  -eps_krylovschur_mpd <mpd> - Sets the maximum projected dimension

768:    Level: advanced

770: .seealso: EPSKrylovSchurGetDimensions(), EPSSetDimensions(), EPSSetInterval()
771: @*/
772: PetscErrorCode EPSKrylovSchurSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
773: {

781:   PetscTryMethod(eps,"EPSKrylovSchurSetDimensions_C",(EPS,PetscInt,PetscInt,PetscInt),(eps,nev,ncv,mpd));
782:   return(0);
783: }

785: static PetscErrorCode EPSKrylovSchurGetDimensions_KrylovSchur(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
786: {
787:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

790:   if (nev) *nev = ctx->nev;
791:   if (ncv) *ncv = ctx->ncv;
792:   if (mpd) *mpd = ctx->mpd;
793:   return(0);
794: }

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

800:    Not Collective

802:    Input Parameter:
803: .  eps - the eigenproblem solver context

805:    Output Parameters:
806: +  nev - number of eigenvalues to compute
807: .  ncv - the maximum dimension of the subspace to be used by the subsolve
808: -  mpd - the maximum dimension allowed for the projected problem

810:    Level: advanced

812: .seealso: EPSKrylovSchurSetDimensions()
813: @*/
814: PetscErrorCode EPSKrylovSchurGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
815: {

820:   PetscUseMethod(eps,"EPSKrylovSchurGetDimensions_C",(EPS,PetscInt*,PetscInt*,PetscInt*),(eps,nev,ncv,mpd));
821:   return(0);
822: }

824: static PetscErrorCode EPSKrylovSchurSetSubintervals_KrylovSchur(EPS eps,PetscReal* subint)
825: {
826:   PetscErrorCode  ierr;
827:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
828:   PetscInt        i;

831:   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()");
832:   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");
833:   if (ctx->subintervals) { PetscFree(ctx->subintervals); }
834:   PetscMalloc1(ctx->npart+1,&ctx->subintervals);
835:   for (i=0;i<ctx->npart+1;i++) ctx->subintervals[i] = subint[i];
836:   ctx->subintset = PETSC_TRUE;
837:   eps->state = EPS_STATE_INITIAL;
838:   return(0);
839: }

841: /*@C
842:    EPSKrylovSchurSetSubintervals - Sets the points that delimit the
843:    subintervals to be used in spectrum slicing with several partitions.

845:    Logically Collective on eps

847:    Input Parameters:
848: +  eps    - the eigenproblem solver context
849: -  subint - array of real values specifying subintervals

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

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

861:    Level: advanced

863: .seealso: EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubintervals(), EPSSetInterval()
864: @*/
865: PetscErrorCode EPSKrylovSchurSetSubintervals(EPS eps,PetscReal *subint)
866: {

871:   PetscTryMethod(eps,"EPSKrylovSchurSetSubintervals_C",(EPS,PetscReal*),(eps,subint));
872:   return(0);
873: }

875: static PetscErrorCode EPSKrylovSchurGetSubintervals_KrylovSchur(EPS eps,PetscReal **subint)
876: {
877:   PetscErrorCode  ierr;
878:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
879:   PetscInt        i;

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

891: /*@C
892:    EPSKrylovSchurGetSubintervals - Returns the points that delimit the
893:    subintervals used in spectrum slicing with several partitions.

895:    Logically Collective on eps

897:    Input Parameter:
898: .  eps    - the eigenproblem solver context

900:    Output Parameter:
901: .  subint - array of real values specifying subintervals

903:    Notes:
904:    If the user passed values with EPSKrylovSchurSetSubintervals(), then the
905:    same values are returned. Otherwise, the values computed internally are
906:    obtained.

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

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

913:    Fortran Notes:
914:    The calling sequence from Fortran is
915: .vb
916:    EPSKrylovSchurGetSubintervals(eps,subint,ierr)
917:    double precision subint(npart+1) output
918: .ve

920:    Level: advanced

922: .seealso: EPSKrylovSchurSetSubintervals(), EPSKrylovSchurGetPartitions(), EPSSetInterval()
923: @*/
924: PetscErrorCode EPSKrylovSchurGetSubintervals(EPS eps,PetscReal **subint)
925: {

931:   PetscUseMethod(eps,"EPSKrylovSchurGetSubintervals_C",(EPS,PetscReal**),(eps,subint));
932:   return(0);
933: }

935: static PetscErrorCode EPSKrylovSchurGetInertias_KrylovSchur(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
936: {
937:   PetscErrorCode  ierr;
938:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
939:   PetscInt        i,numsh;
940:   EPS_SR          sr = ctx->sr;

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

981: /*@C
982:    EPSKrylovSchurGetInertias - Gets the values of the shifts and their
983:    corresponding inertias in case of doing spectrum slicing for a
984:    computational interval.

986:    Not Collective

988:    Input Parameter:
989: .  eps - the eigenproblem solver context

991:    Output Parameters:
992: +  n        - number of shifts, including the endpoints of the interval
993: .  shifts   - the values of the shifts used internally in the solver
994: -  inertias - the values of the inertia in each shift

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

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

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

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

1021:    Level: advanced

1023: .seealso: EPSSetInterval(), EPSKrylovSchurSetSubintervals()
1024: @*/
1025: PetscErrorCode EPSKrylovSchurGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
1026: {

1032:   PetscUseMethod(eps,"EPSKrylovSchurGetInertias_C",(EPS,PetscInt*,PetscReal**,PetscInt**),(eps,n,shifts,inertias));
1033:   return(0);
1034: }

1036: static PetscErrorCode EPSKrylovSchurGetSubcommInfo_KrylovSchur(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
1037: {
1038:   PetscErrorCode  ierr;
1039:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1040:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

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

1053: /*@C
1054:    EPSKrylovSchurGetSubcommInfo - Gets information related to the case of
1055:    doing spectrum slicing for a computational interval with multiple
1056:    communicators.

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

1060:    Input Parameter:
1061: .  eps - the eigenproblem solver context

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

1069:    Notes:
1070:    This function is only available for spectrum slicing runs.

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

1074:    Level: advanced

1076: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommPairs()
1077: @*/
1078: PetscErrorCode EPSKrylovSchurGetSubcommInfo(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
1079: {

1084:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommInfo_C",(EPS,PetscInt*,PetscInt*,Vec*),(eps,k,n,v));
1085:   return(0);
1086: }

1088: static PetscErrorCode EPSKrylovSchurGetSubcommPairs_KrylovSchur(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1089: {
1090:   PetscErrorCode  ierr;
1091:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1092:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

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

1103: /*@
1104:    EPSKrylovSchurGetSubcommPairs - Gets the i-th eigenpair stored
1105:    internally in the subcommunicator to which the calling process belongs.

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

1109:    Input Parameter:
1110: +  eps - the eigenproblem solver context
1111: -  i   - index of the solution

1113:    Output Parameters:
1114: +  eig - the eigenvalue
1115: -  v   - the eigenvector

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

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

1125:    Level: advanced

1127: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo(), EPSKrylovSchurGetSubcommMats()
1128: @*/
1129: PetscErrorCode EPSKrylovSchurGetSubcommPairs(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1130: {

1136:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommPairs_C",(EPS,PetscInt,PetscScalar*,Vec),(eps,i,eig,v));
1137:   return(0);
1138: }

1140: static PetscErrorCode EPSKrylovSchurGetSubcommMats_KrylovSchur(EPS eps,Mat *A,Mat *B)
1141: {
1142:   PetscErrorCode  ierr;
1143:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

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

1152: /*@C
1153:    EPSKrylovSchurGetSubcommMats - Gets the eigenproblem matrices stored
1154:    internally in the subcommunicator to which the calling process belongs.

1156:    Collective on the subcommunicator

1158:    Input Parameter:
1159: .  eps - the eigenproblem solver context

1161:    Output Parameters:
1162: +  A  - the matrix associated with the eigensystem
1163: -  B  - the second matrix in the case of generalized eigenproblems

1165:    Notes:
1166:    This is the analog of EPSGetOperators(), but returns the matrices distributed
1167:    differently (in the subcommunicator rather than in the parent communicator).

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

1171:    Level: advanced

1173: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo()
1174: @*/
1175: PetscErrorCode EPSKrylovSchurGetSubcommMats(EPS eps,Mat *A,Mat *B)
1176: {

1181:   PetscTryMethod(eps,"EPSKrylovSchurGetSubcommMats_C",(EPS,Mat*,Mat*),(eps,A,B));
1182:   return(0);
1183: }

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

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

1198:   MatScale(A,a);
1199:   if (Au) {
1200:     MatAXPY(A,ap,Au,str);
1201:   }
1202:   if (B) MatScale(B,b);
1203:   if (Bu) {
1204:     MatAXPY(B,bp,Bu,str);
1205:   }
1206:   EPSSetOperators(ctx->eps,A,B);

1208:   /* Update stored matrix state */
1209:   subctx = (EPS_KRYLOVSCHUR*)ctx->eps->data;
1210:   PetscObjectStateGet((PetscObject)A,&subctx->Astate);
1211:   if (B) { PetscObjectStateGet((PetscObject)B,&subctx->Bstate); }

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

1241: /*@
1242:    EPSKrylovSchurUpdateSubcommMats - Update the eigenproblem matrices stored
1243:    internally in the subcommunicator to which the calling process belongs.

1245:    Collective on eps

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

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

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

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

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

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

1275:    Level: advanced

1277: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommMats()
1278: @*/
1279: PetscErrorCode EPSKrylovSchurUpdateSubcommMats(EPS eps,PetscScalar s,PetscScalar a,Mat Au,PetscScalar t,PetscScalar b,Mat Bu,MatStructure str,PetscBool globalup)
1280: {

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

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

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

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

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

1380:   ctx_local = (EPS_KRYLOVSCHUR*)ctx->eps->data;
1381:   ctx_local->npart = ctx->npart;
1382:   ctx_local->global = PETSC_FALSE;
1383:   ctx_local->eps = eps;
1384:   ctx_local->subc = ctx->subc;
1385:   ctx_local->commrank = ctx->commrank;
1386:   *child = ctx->eps;
1387:   return(0);
1388: }

1390: static PetscErrorCode EPSKrylovSchurGetKSP_KrylovSchur(EPS eps,KSP *ksp)
1391: {
1392:   PetscErrorCode  ierr;
1393:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1394:   ST              st;
1395:   PetscBool       isfilt;

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

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

1411:    Collective on eps

1413:    Input Parameter:
1414: .  eps - the eigenproblem solver context

1416:    Output Parameters:
1417: -  ksp - the internal KSP object

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

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

1430:    Level: advanced

1432: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions()
1433: @*/
1434: PetscErrorCode EPSKrylovSchurGetKSP(EPS eps,KSP *ksp)
1435: {

1440:   PetscUseMethod(eps,"EPSKrylovSchurGetKSP_C",(EPS,KSP*),(eps,ksp));
1441:   return(0);
1442: }

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

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

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

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

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

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

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

1477:   PetscOptionsTail();

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

1493: PetscErrorCode EPSView_KrylovSchur(EPS eps,PetscViewer viewer)
1494: {
1495:   PetscErrorCode  ierr;
1496:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1497:   PetscBool       isascii,isfilt;
1498:   KSP             ksp;
1499:   PetscViewer     sviewer;

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

1539: PetscErrorCode EPSDestroy_KrylovSchur(EPS eps)
1540: {
1542:   PetscBool      isfilt;

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

1571: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
1572: {
1574:   PetscBool      isfilt;

1577:   PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
1578:   if (eps->which==EPS_ALL && !isfilt) {
1579:     EPSReset_KrylovSchur_Slice(eps);
1580:   }
1581:   return(0);
1582: }

1584: PetscErrorCode EPSSetDefaultST_KrylovSchur(EPS eps)
1585: {

1589:   if (eps->which==EPS_ALL) {
1590:     if (!((PetscObject)eps->st)->type_name) {
1591:       STSetType(eps->st,STSINVERT);
1592:     }
1593:   }
1594:   return(0);
1595: }

1597: SLEPC_EXTERN PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
1598: {
1599:   EPS_KRYLOVSCHUR *ctx;
1600:   PetscErrorCode  ierr;

1603:   PetscNewLog(eps,&ctx);
1604:   eps->data   = (void*)ctx;
1605:   ctx->lock   = PETSC_TRUE;
1606:   ctx->nev    = 1;
1607:   ctx->ncv    = PETSC_DEFAULT;
1608:   ctx->mpd    = PETSC_DEFAULT;
1609:   ctx->npart  = 1;
1610:   ctx->detect = PETSC_FALSE;
1611:   ctx->global = PETSC_TRUE;

1613:   eps->useds = PETSC_TRUE;

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

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