Actual source code: krylovschur.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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 http://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 http://slepc.upv.es.
 30: */

 32: #include <slepc/private/epsimpl.h>                /*I "slepceps.h" I*/
 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:   EPSSetOperators(eps,A,NULL);
 73:   EPSSetProblemType(eps,EPS_HEP);
 74:   EPSSetTolerances(eps,1e-3,50);
 75:   EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
 76:   EPSSolve(eps);
 77:   EPSGetConverged(eps,&nconv);
 78:   if (nconv>0) {
 79:     EPSGetEigenvalue(eps,0,&eig0,NULL);
 80:   } else eig0 = eps->eigr[0];
 81:   *left = PetscRealPart(eig0);
 82:   EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
 83:   EPSSolve(eps);
 84:   EPSGetConverged(eps,&nconv);
 85:   if (nconv>0) {
 86:     EPSGetEigenvalue(eps,0,&eig0,NULL);
 87:   } else eig0 = eps->eigr[0];
 88:   *right = PetscRealPart(eig0);
 89:   EPSDestroy(&eps);
 90:   return(0);
 91: }

 93: static PetscErrorCode EPSSetUp_KrylovSchur_Filter(EPS eps)
 94: {
 96:   SlepcSC        sc;
 97:   PetscReal      rleft,rright;
 98:   Mat            A;

101:   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");
102:   STFilterSetInterval(eps->st,eps->inta,eps->intb);
103:   STGetMatrix(eps->st,0,&A);
104:   EstimateRange(A,&rleft,&rright);
105:   STFilterSetRange(eps->st,rleft,rright);
106:   if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
107:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs cannot be used with spectrum slicing");
108:   if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL*1e-2;  /* use tighter tolerance */
109:   EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
110:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
111:   if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);

113:   DSGetSlepcSC(eps->ds,&sc);
114:   sc->rg            = NULL;
115:   sc->comparison    = SlepcCompareLargestReal;
116:   sc->comparisonctx = NULL;
117:   sc->map           = NULL;
118:   sc->mapobj        = NULL;
119:   return(0);
120: }

122: PetscErrorCode EPSSetUp_KrylovSchur(EPS eps)
123: {
124:   PetscErrorCode    ierr;
125:   PetscReal         eta;
126:   PetscBool         isfilt;
127:   BVOrthogType      otype;
128:   BVOrthogBlockType obtype;
129:   EPS_KRYLOVSCHUR   *ctx = (EPS_KRYLOVSCHUR*)eps->data;
130:   enum { EPS_KS_DEFAULT,EPS_KS_SYMM,EPS_KS_SLICE,EPS_KS_FILTER,EPS_KS_INDEF } variant;

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

149:   if (eps->isgeneralized && eps->ishermitian && !eps->ispositive && eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not implemented for indefinite problems");
150:   if (eps->ishermitian && eps->ispositive && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");

152:   if (!eps->extraction) {
153:     EPSSetExtraction(eps,EPS_RITZ);
154:   } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
155:   if (eps->extraction==EPS_HARMONIC && ctx->lock) { PetscInfo(eps,"Locking was requested but will be deactivated since is not supported with harmonic extraction\n"); }

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 {
182:     switch (eps->extraction) {
183:       case EPS_RITZ:     variant = EPS_KS_DEFAULT; break;
184:       case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
185:       default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
186:     }
187:   }
188:   switch (variant) {
189:     case EPS_KS_DEFAULT:
190:       eps->ops->solve = EPSSolve_KrylovSchur_Default;
191:       eps->ops->computevectors = EPSComputeVectors_Schur;
192:       DSSetType(eps->ds,DSNHEP);
193:       DSAllocate(eps->ds,eps->ncv+1);
194:       break;
195:     case EPS_KS_SYMM:
196:     case EPS_KS_FILTER:
197:       eps->ops->solve = EPSSolve_KrylovSchur_Symm;
198:       eps->ops->computevectors = EPSComputeVectors_Hermitian;
199:       DSSetType(eps->ds,DSHEP);
200:       DSSetCompact(eps->ds,PETSC_TRUE);
201:       DSSetExtraRow(eps->ds,PETSC_TRUE);
202:       DSAllocate(eps->ds,eps->ncv+1);
203:       break;
204:     case EPS_KS_SLICE:
205:       if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing does not support user-defined stopping test");
206:       eps->ops->solve = EPSSolve_KrylovSchur_Slice;
207:       eps->ops->computevectors = EPSComputeVectors_Slice;
208:       break;
209:     case EPS_KS_INDEF:
210:       eps->ops->solve = EPSSolve_KrylovSchur_Indefinite;
211:       eps->ops->computevectors = EPSComputeVectors_Indefinite;
212:       DSSetType(eps->ds,DSGHIEP);
213:       DSSetCompact(eps->ds,PETSC_TRUE);
214:       DSAllocate(eps->ds,eps->ncv+1);
215:       /* force reorthogonalization for pseudo-Lanczos */
216:       BVGetOrthogonalization(eps->V,&otype,NULL,&eta,&obtype);
217:       BVSetOrthogonalization(eps->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
218:       break;
219:     default: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error");
220:   }
221:   return(0);
222: }

224: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
225: {
226:   PetscErrorCode  ierr;
227:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
228:   PetscInt        i,j,*pj,k,l,nv,ld,nconv;
229:   Mat             U;
230:   PetscScalar     *S,*Q,*g;
231:   PetscReal       beta,gamma=1.0;
232:   PetscBool       breakdown,harmonic;

235:   DSGetLeadingDimension(eps->ds,&ld);
236:   harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
237:   if (harmonic) { PetscMalloc1(ld,&g); }
238:   if (eps->arbitrary) pj = &j;
239:   else pj = NULL;

241:   /* Get the starting Arnoldi vector */
242:   EPSGetStartVector(eps,0,NULL);
243:   l = 0;

245:   /* Restart loop */
246:   while (eps->reason == EPS_CONVERGED_ITERATING) {
247:     eps->its++;

249:     /* Compute an nv-step Arnoldi factorization */
250:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
251:     DSGetArray(eps->ds,DS_MAT_A,&S);
252:     EPSBasicArnoldi(eps,PETSC_FALSE,S,ld,eps->nconv+l,&nv,&beta,&breakdown);
253:     DSRestoreArray(eps->ds,DS_MAT_A,&S);
254:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
255:     if (l==0) {
256:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
257:     } else {
258:       DSSetState(eps->ds,DS_STATE_RAW);
259:     }
260:     BVSetActiveColumns(eps->V,eps->nconv,nv);

262:     /* Compute translation of Krylov decomposition if harmonic extraction used */
263:     if (harmonic) {
264:       DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,g,&gamma);
265:     }

267:     /* Solve projected problem */
268:     DSSolve(eps->ds,eps->eigr,eps->eigi);
269:     if (eps->arbitrary) {
270:       EPSGetArbitraryValues(eps,eps->rr,eps->ri);
271:       j=1;
272:     }
273:     DSSort(eps->ds,eps->eigr,eps->eigi,eps->rr,eps->ri,pj);
274:     DSSynchronize(eps->ds,eps->eigr,eps->eigi);

276:     /* Check convergence */
277:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,gamma,&k);
278:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
279:     nconv = k;

281:     /* Update l */
282:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
283:     else {
284:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
285: #if !defined(PETSC_USE_COMPLEX)
286:       DSGetArray(eps->ds,DS_MAT_A,&S);
287:       if (S[k+l+(k+l-1)*ld] != 0.0) {
288:         if (k+l<nv-1) l = l+1;
289:         else l = l-1;
290:       }
291:       DSRestoreArray(eps->ds,DS_MAT_A,&S);
292: #endif
293:     }
294:     if ((!ctx->lock || harmonic) && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */

296:     if (eps->reason == EPS_CONVERGED_ITERATING) {
297:       if (breakdown || k==nv) {
298:         /* Start a new Arnoldi factorization */
299:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
300:         if (k<eps->nev) {
301:           EPSGetStartVector(eps,k,&breakdown);
302:           if (breakdown) {
303:             eps->reason = EPS_DIVERGED_BREAKDOWN;
304:             PetscInfo(eps,"Unable to generate more start vectors\n");
305:           }
306:         }
307:       } else {
308:         /* Undo translation of Krylov decomposition */
309:         if (harmonic) {
310:           DSSetDimensions(eps->ds,nv,0,k,l);
311:           DSTranslateHarmonic(eps->ds,0.0,beta,PETSC_TRUE,g,&gamma);
312:           /* gamma u^ = u - U*g~ */
313:           BVMultColumn(eps->V,-1.0,1.0,nv,g);
314:           BVScaleColumn(eps->V,nv,1.0/gamma);
315:         }
316:         /* Prepare the Rayleigh quotient for restart */
317:         DSGetArray(eps->ds,DS_MAT_A,&S);
318:         DSGetArray(eps->ds,DS_MAT_Q,&Q);
319:         for (i=k;i<k+l;i++) {
320:           S[k+l+i*ld] = Q[nv-1+i*ld]*beta*gamma;
321:         }
322:         DSRestoreArray(eps->ds,DS_MAT_A,&S);
323:         DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
324:       }
325:     }
326:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
327:     DSGetMat(eps->ds,DS_MAT_Q,&U);
328:     BVMultInPlace(eps->V,U,eps->nconv,k+l);
329:     MatDestroy(&U);

331:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
332:       BVCopyColumn(eps->V,nv,k+l);
333:     }
334:     eps->nconv = k;
335:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
336:   }

338:   if (harmonic) { PetscFree(g); }
339:   /* truncate Schur decomposition and change the state to raw so that
340:      DSVectors() computes eigenvectors from scratch */
341:   DSSetDimensions(eps->ds,eps->nconv,0,0,0);
342:   DSSetState(eps->ds,DS_STATE_RAW);
343:   return(0);
344: }

346: static PetscErrorCode EPSKrylovSchurSetRestart_KrylovSchur(EPS eps,PetscReal keep)
347: {
348:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

351:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
352:   else {
353:     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);
354:     ctx->keep = keep;
355:   }
356:   return(0);
357: }

359: /*@
360:    EPSKrylovSchurSetRestart - Sets the restart parameter for the Krylov-Schur
361:    method, in particular the proportion of basis vectors that must be kept
362:    after restart.

364:    Logically Collective on EPS

366:    Input Parameters:
367: +  eps - the eigenproblem solver context
368: -  keep - the number of vectors to be kept at restart

370:    Options Database Key:
371: .  -eps_krylovschur_restart - Sets the restart parameter

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

376:    Level: advanced

378: .seealso: EPSKrylovSchurGetRestart()
379: @*/
380: PetscErrorCode EPSKrylovSchurSetRestart(EPS eps,PetscReal keep)
381: {

387:   PetscTryMethod(eps,"EPSKrylovSchurSetRestart_C",(EPS,PetscReal),(eps,keep));
388:   return(0);
389: }

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

396:   *keep = ctx->keep;
397:   return(0);
398: }

400: /*@
401:    EPSKrylovSchurGetRestart - Gets the restart parameter used in the
402:    Krylov-Schur method.

404:    Not Collective

406:    Input Parameter:
407: .  eps - the eigenproblem solver context

409:    Output Parameter:
410: .  keep - the restart parameter

412:    Level: advanced

414: .seealso: EPSKrylovSchurSetRestart()
415: @*/
416: PetscErrorCode EPSKrylovSchurGetRestart(EPS eps,PetscReal *keep)
417: {

423:   PetscUseMethod(eps,"EPSKrylovSchurGetRestart_C",(EPS,PetscReal*),(eps,keep));
424:   return(0);
425: }

427: static PetscErrorCode EPSKrylovSchurSetLocking_KrylovSchur(EPS eps,PetscBool lock)
428: {
429:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

432:   ctx->lock = lock;
433:   return(0);
434: }

436: /*@
437:    EPSKrylovSchurSetLocking - Choose between locking and non-locking variants of
438:    the Krylov-Schur method.

440:    Logically Collective on EPS

442:    Input Parameters:
443: +  eps  - the eigenproblem solver context
444: -  lock - true if the locking variant must be selected

446:    Options Database Key:
447: .  -eps_krylovschur_locking - Sets the locking flag

449:    Notes:
450:    The default is to lock converged eigenpairs when the method restarts.
451:    This behaviour can be changed so that all directions are kept in the
452:    working subspace even if already converged to working accuracy (the
453:    non-locking variant).

455:    Level: advanced

457: .seealso: EPSKrylovSchurGetLocking()
458: @*/
459: PetscErrorCode EPSKrylovSchurSetLocking(EPS eps,PetscBool lock)
460: {

466:   PetscTryMethod(eps,"EPSKrylovSchurSetLocking_C",(EPS,PetscBool),(eps,lock));
467:   return(0);
468: }

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

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

479: /*@
480:    EPSKrylovSchurGetLocking - Gets the locking flag used in the Krylov-Schur
481:    method.

483:    Not Collective

485:    Input Parameter:
486: .  eps - the eigenproblem solver context

488:    Output Parameter:
489: .  lock - the locking flag

491:    Level: advanced

493: .seealso: EPSKrylovSchurSetLocking()
494: @*/
495: PetscErrorCode EPSKrylovSchurGetLocking(EPS eps,PetscBool *lock)
496: {

502:   PetscUseMethod(eps,"EPSKrylovSchurGetLocking_C",(EPS,PetscBool*),(eps,lock));
503:   return(0);
504: }

506: static PetscErrorCode EPSKrylovSchurSetPartitions_KrylovSchur(EPS eps,PetscInt npart)
507: {
508:   PetscErrorCode  ierr;
509:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
510:   PetscMPIInt     size;

513:   if (ctx->npart!=npart) {
514:     if (ctx->commset) { PetscSubcommDestroy(&ctx->subc); }
515:     EPSDestroy(&ctx->eps);
516:   }
517:   if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
518:     ctx->npart = 1;
519:   } else {
520:     MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size);
521:     if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
522:     ctx->npart = npart;
523:   }
524:   eps->state = EPS_STATE_INITIAL;
525:   return(0);
526: }

528: /*@
529:    EPSKrylovSchurSetPartitions - Sets the number of partitions for the
530:    case of doing spectrum slicing for a computational interval with the
531:    communicator split in several sub-communicators.

533:    Logically Collective on EPS

535:    Input Parameters:
536: +  eps   - the eigenproblem solver context
537: -  npart - number of partitions

539:    Options Database Key:
540: .  -eps_krylovschur_partitions <npart> - Sets the number of partitions

542:    Notes:
543:    By default, npart=1 so all processes in the communicator participate in
544:    the processing of the whole interval. If npart>1 then the interval is
545:    divided into npart subintervals, each of them being processed by a
546:    subset of processes.

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

551:    Level: advanced

553: .seealso: EPSKrylovSchurSetSubintervals(), EPSSetInterval()
554: @*/
555: PetscErrorCode EPSKrylovSchurSetPartitions(EPS eps,PetscInt npart)
556: {

562:   PetscTryMethod(eps,"EPSKrylovSchurSetPartitions_C",(EPS,PetscInt),(eps,npart));
563:   return(0);
564: }

566: static PetscErrorCode EPSKrylovSchurGetPartitions_KrylovSchur(EPS eps,PetscInt *npart)
567: {
568:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

571:   *npart  = ctx->npart;
572:   return(0);
573: }

575: /*@
576:    EPSKrylovSchurGetPartitions - Gets the number of partitions of the
577:    communicator in case of spectrum slicing.

579:    Not Collective

581:    Input Parameter:
582: .  eps - the eigenproblem solver context

584:    Output Parameter:
585: .  npart - number of partitions

587:    Level: advanced

589: .seealso: EPSKrylovSchurSetPartitions()
590: @*/
591: PetscErrorCode EPSKrylovSchurGetPartitions(EPS eps,PetscInt *npart)
592: {

598:   PetscUseMethod(eps,"EPSKrylovSchurGetPartitions_C",(EPS,PetscInt*),(eps,npart));
599:   return(0);
600: }

602: static PetscErrorCode EPSKrylovSchurSetDetectZeros_KrylovSchur(EPS eps,PetscBool detect)
603: {
604:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

607:   ctx->detect = detect;
608:   eps->state  = EPS_STATE_INITIAL;
609:   return(0);
610: }

612: /*@
613:    EPSKrylovSchurSetDetectZeros - Sets a flag to enforce detection of
614:    zeros during the factorizations throughout the spectrum slicing computation.

616:    Logically Collective on EPS

618:    Input Parameters:
619: +  eps    - the eigenproblem solver context
620: -  detect - check for zeros

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

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

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

634:    Level: advanced

636: .seealso: EPSKrylovSchurSetPartitions(), EPSSetInterval()
637: @*/
638: PetscErrorCode EPSKrylovSchurSetDetectZeros(EPS eps,PetscBool detect)
639: {

645:   PetscTryMethod(eps,"EPSKrylovSchurSetDetectZeros_C",(EPS,PetscBool),(eps,detect));
646:   return(0);
647: }

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

654:   *detect = ctx->detect;
655:   return(0);
656: }

658: /*@
659:    EPSKrylovSchurGetDetectZeros - Gets the flag that enforces zero detection
660:    in spectrum slicing.

662:    Not Collective

664:    Input Parameter:
665: .  eps - the eigenproblem solver context

667:    Output Parameter:
668: .  detect - whether zeros detection is enforced during factorizations

670:    Level: advanced

672: .seealso: EPSKrylovSchurSetDetectZeros()
673: @*/
674: PetscErrorCode EPSKrylovSchurGetDetectZeros(EPS eps,PetscBool *detect)
675: {

681:   PetscUseMethod(eps,"EPSKrylovSchurGetDetectZeros_C",(EPS,PetscBool*),(eps,detect));
682:   return(0);
683: }

685: static PetscErrorCode EPSKrylovSchurSetDimensions_KrylovSchur(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
686: {
687:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

690:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
691:   ctx->nev = nev;
692:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
693:     ctx->ncv = 0;
694:   } else {
695:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
696:     ctx->ncv = ncv;
697:   }
698:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
699:     ctx->mpd = 0;
700:   } else {
701:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
702:     ctx->mpd = mpd;
703:   }
704:   eps->state = EPS_STATE_INITIAL;
705:   return(0);
706: }

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

713:    Logically Collective on EPS

715:    Input Parameters:
716: +  eps - the eigenproblem solver context
717: .  nev - number of eigenvalues to compute
718: .  ncv - the maximum dimension of the subspace to be used by the subsolve
719: -  mpd - the maximum dimension allowed for the projected problem

721:    Options Database Key:
722: +  -eps_krylovschur_nev <nev> - Sets the number of eigenvalues
723: .  -eps_krylovschur_ncv <ncv> - Sets the dimension of the subspace
724: -  -eps_krylovschur_mpd <mpd> - Sets the maximum projected dimension

726:    Level: advanced

728: .seealso: EPSKrylovSchurGetDimensions(), EPSSetDimensions(), EPSSetInterval()
729: @*/
730: PetscErrorCode EPSKrylovSchurSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
731: {

739:   PetscTryMethod(eps,"EPSKrylovSchurSetDimensions_C",(EPS,PetscInt,PetscInt,PetscInt),(eps,nev,ncv,mpd));
740:   return(0);
741: }

743: static PetscErrorCode EPSKrylovSchurGetDimensions_KrylovSchur(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
744: {
745:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

748:   if (nev) *nev = ctx->nev;
749:   if (ncv) *ncv = ctx->ncv;
750:   if (mpd) *mpd = ctx->mpd;
751:   return(0);
752: }

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

758:    Not Collective

760:    Input Parameter:
761: .  eps - the eigenproblem solver context

763:    Output Parameters:
764: +  nev - number of eigenvalues to compute
765: .  ncv - the maximum dimension of the subspace to be used by the subsolve
766: -  mpd - the maximum dimension allowed for the projected problem

768:    Level: advanced

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

778:   PetscUseMethod(eps,"EPSKrylovSchurGetDimensions_C",(EPS,PetscInt*,PetscInt*,PetscInt*),(eps,nev,ncv,mpd));
779:   return(0);
780: }

782: static PetscErrorCode EPSKrylovSchurSetSubintervals_KrylovSchur(EPS eps,PetscReal* subint)
783: {
784:   PetscErrorCode  ierr;
785:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
786:   PetscInt        i;

789:   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()");
790:   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");
791:   if (ctx->subintervals) { PetscFree(ctx->subintervals); }
792:   PetscMalloc1(ctx->npart+1,&ctx->subintervals);
793:   for (i=0;i<ctx->npart+1;i++) ctx->subintervals[i] = subint[i];
794:   ctx->subintset = PETSC_TRUE;
795:   eps->state = EPS_STATE_INITIAL;
796:   return(0);
797: }

799: /*@C
800:    EPSKrylovSchurSetSubintervals - Sets the points that delimit the
801:    subintervals to be used in spectrum slicing with several partitions.

803:    Logically Collective on EPS

805:    Input Parameters:
806: +  eps    - the eigenproblem solver context
807: -  subint - array of real values specifying subintervals

809:    Notes:
810:    This function must be called after EPSKrylovSchurSetPartitions(). For npart
811:    partitions, the argument subint must contain npart+1 real values sorted in
812:    ascending order: subint_0, subint_1, ..., subint_npart, where the first
813:    and last values must coincide with the interval endpoints set with
814:    EPSSetInterval().

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

819:    Level: advanced

821: .seealso: EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubintervals(), EPSSetInterval()
822: @*/
823: PetscErrorCode EPSKrylovSchurSetSubintervals(EPS eps,PetscReal *subint)
824: {

829:   PetscTryMethod(eps,"EPSKrylovSchurSetSubintervals_C",(EPS,PetscReal*),(eps,subint));
830:   return(0);
831: }

833: static PetscErrorCode EPSKrylovSchurGetSubintervals_KrylovSchur(EPS eps,PetscReal **subint)
834: {
835:   PetscErrorCode  ierr;
836:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
837:   PetscInt        i;

840:   if (!ctx->subintset) {
841:     if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
842:     if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
843:   }
844:   PetscMalloc1(ctx->npart+1,subint);
845:   for (i=0;i<=ctx->npart;i++) (*subint)[i] = ctx->subintervals[i];
846:   return(0);
847: }

849: /*@C
850:    EPSKrylovSchurGetSubintervals - Returns the points that delimit the
851:    subintervals used in spectrum slicing with several partitions.

853:    Logically Collective on EPS

855:    Input Parameter:
856: .  eps    - the eigenproblem solver context

858:    Output Parameter:
859: .  subint - array of real values specifying subintervals

861:    Notes:
862:    If the user passed values with EPSKrylovSchurSetSubintervals(), then the
863:    same values are returned. Otherwise, the values computed internally are
864:    obtained.

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

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

871:    Fortran Notes:
872:    The calling sequence from Fortran is
873: .vb
874:    EPSKrylovSchurGetSubintervals(eps,subint,ierr)
875:    double precision subint(npart+1) output
876: .ve

878:    Level: advanced

880: .seealso: EPSKrylovSchurSetSubintervals(), EPSKrylovSchurGetPartitions(), EPSSetInterval()
881: @*/
882: PetscErrorCode EPSKrylovSchurGetSubintervals(EPS eps,PetscReal **subint)
883: {

889:   PetscUseMethod(eps,"EPSKrylovSchurGetSubintervals_C",(EPS,PetscReal**),(eps,subint));
890:   return(0);
891: }

893: static PetscErrorCode EPSKrylovSchurGetInertias_KrylovSchur(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
894: {
895:   PetscErrorCode  ierr;
896:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
897:   PetscInt        i,numsh;
898:   EPS_SR          sr = ctx->sr;

901:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
902:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
903:   switch (eps->state) {
904:   case EPS_STATE_INITIAL:
905:     break;
906:   case EPS_STATE_SETUP:
907:     numsh = ctx->npart+1;
908:     if (n) *n = numsh;
909:     if (shifts) {
910:       PetscMalloc1(numsh,shifts);
911:       (*shifts)[0] = eps->inta;
912:       if (ctx->npart==1) (*shifts)[1] = eps->intb;
913:       else for (i=1;i<numsh;i++) (*shifts)[i] = ctx->subintervals[i];
914:     }
915:     if (inertias) {
916:       PetscMalloc1(numsh,inertias);
917:       (*inertias)[0] = (sr->dir==1)?sr->inertia0:sr->inertia1;
918:       if (ctx->npart==1) (*inertias)[1] = (sr->dir==1)?sr->inertia1:sr->inertia0;
919:       else for (i=1;i<numsh;i++) (*inertias)[i] = (*inertias)[i-1]+ctx->nconv_loc[i-1];
920:     }
921:     break;
922:   case EPS_STATE_SOLVED:
923:   case EPS_STATE_EIGENVECTORS:
924:     numsh = ctx->nshifts;
925:     if (n) *n = numsh;
926:     if (shifts) {
927:       PetscMalloc1(numsh,shifts);
928:       for (i=0;i<numsh;i++) (*shifts)[i] = ctx->shifts[i];
929:     }
930:     if (inertias) {
931:       PetscMalloc1(numsh,inertias);
932:       for (i=0;i<numsh;i++) (*inertias)[i] = ctx->inertias[i];
933:     }
934:     break;
935:   }
936:   return(0);
937: }

939: /*@C
940:    EPSKrylovSchurGetInertias - Gets the values of the shifts and their
941:    corresponding inertias in case of doing spectrum slicing for a
942:    computational interval.

944:    Not Collective

946:    Input Parameter:
947: .  eps - the eigenproblem solver context

949:    Output Parameters:
950: +  n        - number of shifts, including the endpoints of the interval
951: .  shifts   - the values of the shifts used internally in the solver
952: -  inertias - the values of the inertia in each shift

954:    Notes:
955:    If called after EPSSolve(), all shifts used internally by the solver are
956:    returned (including both endpoints and any intermediate ones). If called
957:    before EPSSolve() and after EPSSetUp() then only the information of the
958:    endpoints of subintervals is available.

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

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

965:    Fortran Notes:
966:    The calling sequence from Fortran is
967: .vb
968:    EPSKrylovSchurGetInertias(eps,n,shifts,inertias,ierr)
969:    integer n
970:    double precision shifts(*)
971:    integer inertias(*)
972: .ve
973:    The arrays should be at least of length n. The value of n can be determined
974:    by an initial call
975: .vb
976:    EPSKrylovSchurGetInertias(eps,n,PETSC_NULL_REAL,PETSC_NULL_INTEGER,ierr)
977: .ve

979:    Level: advanced

981: .seealso: EPSSetInterval(), EPSKrylovSchurSetSubintervals()
982: @*/
983: PetscErrorCode EPSKrylovSchurGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
984: {

990:   PetscUseMethod(eps,"EPSKrylovSchurGetInertias_C",(EPS,PetscInt*,PetscReal**,PetscInt**),(eps,n,shifts,inertias));
991:   return(0);
992: }

994: static PetscErrorCode EPSKrylovSchurGetSubcommInfo_KrylovSchur(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
995: {
996:   PetscErrorCode  ierr;
997:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
998:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

1001:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1002:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1003:   if (k) *k = (ctx->npart==1)? 0: ctx->subc->color;
1004:   if (n) *n = sr->numEigs;
1005:   if (v) {
1006:     BVCreateVec(sr->V,v);
1007:   }
1008:   return(0);
1009: }

1011: /*@C
1012:    EPSKrylovSchurGetSubcommInfo - Gets information related to the case of
1013:    doing spectrum slicing for a computational interval with multiple
1014:    communicators.

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

1018:    Input Parameter:
1019: .  eps - the eigenproblem solver context

1021:    Output Parameters:
1022: +  k - index of the subinterval for the calling process
1023: .  n - number of eigenvalues found in the k-th subinterval
1024: -  v - a vector owned by processes in the subcommunicator with dimensions
1025:        compatible for locally computed eigenvectors (or NULL)

1027:    Notes:
1028:    This function is only available for spectrum slicing runs.

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

1032:    Level: advanced

1034: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommPairs()
1035: @*/
1036: PetscErrorCode EPSKrylovSchurGetSubcommInfo(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
1037: {

1042:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommInfo_C",(EPS,PetscInt*,PetscInt*,Vec*),(eps,k,n,v));
1043:   return(0);
1044: }

1046: static PetscErrorCode EPSKrylovSchurGetSubcommPairs_KrylovSchur(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1047: {
1048:   PetscErrorCode  ierr;
1049:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1050:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

1053:   EPSCheckSolved(eps,1);
1054:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1055:   if (i<0 || i>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
1056:   if (eig) *eig = sr->eigr[sr->perm[i]];
1057:   if (v) { BVCopyVec(sr->V,sr->perm[i],v); }
1058:   return(0);
1059: }

1061: /*@C
1062:    EPSKrylovSchurGetSubcommPairs - Gets the i-th eigenpair stored
1063:    internally in the subcommunicator to which the calling process belongs.

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

1067:    Input Parameter:
1068: +  eps - the eigenproblem solver context
1069: -  i   - index of the solution

1071:    Output Parameters:
1072: +  eig - the eigenvalue
1073: -  v   - the eigenvector

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

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

1083:    Level: advanced

1085: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo(), EPSKrylovSchurGetSubcommMats()
1086: @*/
1087: PetscErrorCode EPSKrylovSchurGetSubcommPairs(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1088: {

1094:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommPairs_C",(EPS,PetscInt,PetscScalar*,Vec),(eps,i,eig,v));
1095:   return(0);
1096: }

1098: static PetscErrorCode EPSKrylovSchurGetSubcommMats_KrylovSchur(EPS eps,Mat *A,Mat *B)
1099: {
1100:   PetscErrorCode  ierr;
1101:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

1104:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1105:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1106:   EPSGetOperators(ctx->eps,A,B);
1107:   return(0);
1108: }

1110: /*@C
1111:    EPSKrylovSchurGetSubcommMats - Gets the eigenproblem matrices stored
1112:    internally in the subcommunicator to which the calling process belongs.

1114:    Collective on the subcommunicator

1116:    Input Parameter:
1117: .  eps - the eigenproblem solver context

1119:    Output Parameters:
1120: +  A  - the matrix associated with the eigensystem
1121: -  B  - the second matrix in the case of generalized eigenproblems

1123:    Notes:
1124:    This is the analog of EPSGetOperators(), but returns the matrices distributed
1125:    differently (in the subcommunicator rather than in the parent communicator).

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

1129:    Level: advanced

1131: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo()
1132: @*/
1133: PetscErrorCode EPSKrylovSchurGetSubcommMats(EPS eps,Mat *A,Mat *B)
1134: {

1139:   PetscTryMethod(eps,"EPSKrylovSchurGetSubcommMats_C",(EPS,Mat*,Mat*),(eps,A,B));
1140:   return(0);
1141: }

1143: static PetscErrorCode EPSKrylovSchurUpdateSubcommMats_KrylovSchur(EPS eps,PetscScalar a,PetscScalar ap,Mat Au,PetscScalar b,PetscScalar bp, Mat Bu,MatStructure str,PetscBool globalup)
1144: {
1145:   PetscErrorCode  ierr;
1146:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data,*subctx;
1147:   Mat             A,B=NULL,Ag,Bg=NULL;
1148:   PetscBool       reuse=PETSC_TRUE;

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

1156:   MatScale(A,a);
1157:   if (Au) {
1158:     MatAXPY(A,ap,Au,str);
1159:   }
1160:   if (B) MatScale(B,b);
1161:   if (Bu) {
1162:     MatAXPY(B,bp,Bu,str);
1163:   }
1164:   EPSSetOperators(ctx->eps,A,B);

1166:   /* Update stored matrix state */
1167:   subctx = (EPS_KRYLOVSCHUR*)ctx->eps->data;
1168:   PetscObjectStateGet((PetscObject)A,&subctx->Astate);
1169:   if (B) { PetscObjectStateGet((PetscObject)B,&subctx->Bstate); }

1171:   /* Update matrices in the parent communicator if requested by user */
1172:   if (globalup) {
1173:     if (ctx->npart>1) {
1174:       if (!ctx->isrow) {
1175:         MatGetOwnershipIS(Ag,&ctx->isrow,&ctx->iscol);
1176:         reuse = PETSC_FALSE;
1177:       }
1178:       if (str==DIFFERENT_NONZERO_PATTERN) reuse = PETSC_FALSE;
1179:       if (ctx->submata && !reuse) {
1180:         MatDestroyMatrices(1,&ctx->submata);
1181:       }
1182:       MatCreateSubMatrices(A,1,&ctx->isrow,&ctx->iscol,(reuse)?MAT_REUSE_MATRIX:MAT_INITIAL_MATRIX,&ctx->submata);
1183:       MatCreateMPIMatConcatenateSeqMat(((PetscObject)Ag)->comm,ctx->submata[0],PETSC_DECIDE,MAT_REUSE_MATRIX,&Ag);
1184:       if (B) {
1185:         if (ctx->submatb && !reuse) {
1186:           MatDestroyMatrices(1,&ctx->submatb);
1187:         }
1188:         MatCreateSubMatrices(B,1,&ctx->isrow,&ctx->iscol,(reuse)?MAT_REUSE_MATRIX:MAT_INITIAL_MATRIX,&ctx->submatb);
1189:         MatCreateMPIMatConcatenateSeqMat(((PetscObject)Bg)->comm,ctx->submatb[0],PETSC_DECIDE,MAT_REUSE_MATRIX,&Bg);
1190:       }
1191:     }
1192:     PetscObjectStateGet((PetscObject)Ag,&ctx->Astate);
1193:     if (Bg) { PetscObjectStateGet((PetscObject)Bg,&ctx->Bstate); }
1194:   }
1195:   EPSSetOperators(eps,Ag,Bg);
1196:   return(0);
1197: }

1199: /*@C
1200:    EPSKrylovSchurUpdateSubcommMats - Update the eigenproblem matrices stored
1201:    internally in the subcommunicator to which the calling process belongs.

1203:    Collective on EPS

1205:    Input Parameters:
1206: +  eps - the eigenproblem solver context
1207: .  s   - scalar that multiplies the existing A matrix
1208: .  a   - scalar used in the axpy operation on A
1209: .  Au  - matrix used in the axpy operation on A
1210: .  t   - scalar that multiplies the existing B matrix
1211: .  b   - scalar used in the axpy operation on B
1212: .  Bu  - matrix used in the axpy operation on B
1213: .  str - structure flag
1214: -  globalup - flag indicating if global matrices must be updated

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

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

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

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

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

1233:    Level: advanced

1235: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommMats()
1236: @*/
1237: PetscErrorCode EPSKrylovSchurUpdateSubcommMats(EPS eps,PetscScalar s,PetscScalar a,Mat Au,PetscScalar t,PetscScalar b,Mat Bu,MatStructure str,PetscBool globalup)
1238: {

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

1255: PetscErrorCode EPSSetFromOptions_KrylovSchur(PetscOptionItems *PetscOptionsObject,EPS eps)
1256: {
1257:   PetscErrorCode  ierr;
1258:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1259:   PetscBool       flg,lock,b,f1,f2,f3;
1260:   PetscReal       keep;
1261:   PetscInt        i,j,k;

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

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

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

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

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

1280:     i = 1;
1281:     j = k = PETSC_DECIDE;
1282:     PetscOptionsInt("-eps_krylovschur_nev","Number of eigenvalues to compute in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",40,&i,&f1);
1283:     PetscOptionsInt("-eps_krylovschur_ncv","Number of basis vectors in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&j,&f2);
1284:     PetscOptionsInt("-eps_krylovschur_mpd","Maximum dimension of projected problem in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&k,&f3);
1285:     if (f1 || f2 || f3) { EPSKrylovSchurSetDimensions(eps,i,j,k); }

1287:   PetscOptionsTail();
1288:   return(0);
1289: }

1291: PetscErrorCode EPSView_KrylovSchur(EPS eps,PetscViewer viewer)
1292: {
1293:   PetscErrorCode  ierr;
1294:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1295:   PetscBool       isascii,isfilt;

1298:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1299:   if (isascii) {
1300:     PetscViewerASCIIPrintf(viewer,"  %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
1301:     PetscViewerASCIIPrintf(viewer,"  using the %slocking variant\n",ctx->lock?"":"non-");
1302:     if (eps->which==EPS_ALL) {
1303:       PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
1304:       if (isfilt) {
1305:         PetscViewerASCIIPrintf(viewer,"  using filtering to extract all eigenvalues in an interval\n");
1306:       } else {
1307:         PetscViewerASCIIPrintf(viewer,"  doing spectrum slicing with nev=%D, ncv=%D, mpd=%D\n",ctx->nev,ctx->ncv,ctx->mpd);
1308:         if (ctx->npart>1) {
1309:           PetscViewerASCIIPrintf(viewer,"  multi-communicator spectrum slicing with %D partitions\n",ctx->npart);
1310:           if (ctx->detect) { PetscViewerASCIIPrintf(viewer,"  detecting zeros when factorizing at subinterval boundaries\n"); }
1311:         }
1312:       }
1313:     }
1314:   }
1315:   return(0);
1316: }

1318: PetscErrorCode EPSDestroy_KrylovSchur(EPS eps)
1319: {

1323:   PetscFree(eps->data);
1324:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",NULL);
1325:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",NULL);
1326:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetLocking_C",NULL);
1327:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetLocking_C",NULL);
1328:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetPartitions_C",NULL);
1329:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetPartitions_C",NULL);
1330:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDetectZeros_C",NULL);
1331:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDetectZeros_C",NULL);
1332:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",NULL);
1333:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",NULL);
1334:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetSubintervals_C",NULL);
1335:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubintervals_C",NULL);
1336:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",NULL);
1337:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommInfo_C",NULL);
1338:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommPairs_C",NULL);
1339:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommMats_C",NULL);
1340:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurUpdateSubcommMats_C",NULL);
1341:   return(0);
1342: }

1344: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
1345: {
1347:   PetscBool      isfilt;

1350:   PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
1351:   if (eps->which==EPS_ALL && !isfilt) {
1352:     EPSReset_KrylovSchur_Slice(eps);
1353:   }
1354:   return(0);
1355: }

1357: PetscErrorCode EPSSetDefaultST_KrylovSchur(EPS eps)
1358: {

1362:   if (eps->which==EPS_ALL) {
1363:     if (!((PetscObject)eps->st)->type_name) {
1364:       STSetType(eps->st,STSINVERT);
1365:     }
1366:   }
1367:   return(0);
1368: }

1370: PETSC_EXTERN PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
1371: {
1372:   EPS_KRYLOVSCHUR *ctx;
1373:   PetscErrorCode  ierr;

1376:   PetscNewLog(eps,&ctx);
1377:   eps->data   = (void*)ctx;
1378:   ctx->lock   = PETSC_TRUE;
1379:   ctx->nev    = 1;
1380:   ctx->npart  = 1;
1381:   ctx->detect = PETSC_FALSE;
1382:   ctx->global = PETSC_TRUE;

1384:   eps->useds = PETSC_TRUE;

1386:   /* solve and computevectors determined at setup */
1387:   eps->ops->setup          = EPSSetUp_KrylovSchur;
1388:   eps->ops->setfromoptions = EPSSetFromOptions_KrylovSchur;
1389:   eps->ops->destroy        = EPSDestroy_KrylovSchur;
1390:   eps->ops->reset          = EPSReset_KrylovSchur;
1391:   eps->ops->view           = EPSView_KrylovSchur;
1392:   eps->ops->backtransform  = EPSBackTransform_Default;
1393:   eps->ops->setdefaultst   = EPSSetDefaultST_KrylovSchur;

1395:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",EPSKrylovSchurSetRestart_KrylovSchur);
1396:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",EPSKrylovSchurGetRestart_KrylovSchur);
1397:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetLocking_C",EPSKrylovSchurSetLocking_KrylovSchur);
1398:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetLocking_C",EPSKrylovSchurGetLocking_KrylovSchur);
1399:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetPartitions_C",EPSKrylovSchurSetPartitions_KrylovSchur);
1400:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetPartitions_C",EPSKrylovSchurGetPartitions_KrylovSchur);
1401:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDetectZeros_C",EPSKrylovSchurSetDetectZeros_KrylovSchur);
1402:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDetectZeros_C",EPSKrylovSchurGetDetectZeros_KrylovSchur);
1403:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",EPSKrylovSchurSetDimensions_KrylovSchur);
1404:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",EPSKrylovSchurGetDimensions_KrylovSchur);
1405:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetSubintervals_C",EPSKrylovSchurSetSubintervals_KrylovSchur);
1406:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubintervals_C",EPSKrylovSchurGetSubintervals_KrylovSchur);
1407:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",EPSKrylovSchurGetInertias_KrylovSchur);
1408:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommInfo_C",EPSKrylovSchurGetSubcommInfo_KrylovSchur);
1409:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommPairs_C",EPSKrylovSchurGetSubcommPairs_KrylovSchur);
1410:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommMats_C",EPSKrylovSchurGetSubcommMats_KrylovSchur);
1411:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurUpdateSubcommMats_C",EPSKrylovSchurUpdateSubcommMats_KrylovSchur);
1412:   return(0);
1413: }