Actual source code: krylovschur.c

slepc-3.9.0 2018-04-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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:   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: {
 97:   SlepcSC        sc;
 98:   PetscReal      rleft,rright;
 99:   Mat            A;

102:   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");
103:   if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Polynomial filter only available for symmetric/Hermitian eigenproblems");
104:   if (eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Polynomial filters not available for generalized eigenproblems");
105:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs cannot be used with polynomial filters");
106:   if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL*1e-2;  /* use tighter tolerance */
107:   STFilterSetInterval(eps->st,eps->inta,eps->intb);
108:   STGetMatrix(eps->st,0,&A);
109:   STFilterGetRange(eps->st,&rleft,&rright);
110:   if (!rleft && !rright) {
111:     EstimateRange(A,&rleft,&rright);
112:     PetscInfo2(eps,"Setting eigenvalue range to [%g,%g]\n",(double)rleft,(double)rright);
113:     STFilterSetRange(eps->st,rleft,rright);
114:   }
115:   if (!eps->ncv && eps->nev==1) eps->nev = 40;  /* user did not provide nev estimation */
116:   EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
117:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
118:   if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);

120:   DSGetSlepcSC(eps->ds,&sc);
121:   sc->rg            = NULL;
122:   sc->comparison    = SlepcCompareLargestReal;
123:   sc->comparisonctx = NULL;
124:   sc->map           = NULL;
125:   sc->mapobj        = NULL;
126:   return(0);
127: }

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

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

156:   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");
157:   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");

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

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

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

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

231: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
232: {
233:   PetscErrorCode  ierr;
234:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
235:   PetscInt        i,j,*pj,k,l,nv,ld,nconv;
236:   Mat             U;
237:   PetscScalar     *S,*Q,*g;
238:   PetscReal       beta,gamma=1.0;
239:   PetscBool       breakdown,harmonic;

242:   DSGetLeadingDimension(eps->ds,&ld);
243:   harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
244:   if (harmonic) { PetscMalloc1(ld,&g); }
245:   if (eps->arbitrary) pj = &j;
246:   else pj = NULL;

248:   /* Get the starting Arnoldi vector */
249:   EPSGetStartVector(eps,0,NULL);
250:   l = 0;

252:   /* Restart loop */
253:   while (eps->reason == EPS_CONVERGED_ITERATING) {
254:     eps->its++;

256:     /* Compute an nv-step Arnoldi factorization */
257:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
258:     DSGetArray(eps->ds,DS_MAT_A,&S);
259:     EPSBasicArnoldi(eps,PETSC_FALSE,S,ld,eps->nconv+l,&nv,&beta,&breakdown);
260:     DSRestoreArray(eps->ds,DS_MAT_A,&S);
261:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
262:     if (l==0) {
263:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
264:     } else {
265:       DSSetState(eps->ds,DS_STATE_RAW);
266:     }
267:     BVSetActiveColumns(eps->V,eps->nconv,nv);

269:     /* Compute translation of Krylov decomposition if harmonic extraction used */
270:     if (harmonic) {
271:       DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,g,&gamma);
272:     }

274:     /* Solve projected problem */
275:     DSSolve(eps->ds,eps->eigr,eps->eigi);
276:     if (eps->arbitrary) {
277:       EPSGetArbitraryValues(eps,eps->rr,eps->ri);
278:       j=1;
279:     }
280:     DSSort(eps->ds,eps->eigr,eps->eigi,eps->rr,eps->ri,pj);
281:     DSSynchronize(eps->ds,eps->eigr,eps->eigi);

283:     /* Check convergence */
284:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,gamma,&k);
285:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
286:     nconv = k;

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

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

338:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
339:       BVCopyColumn(eps->V,nv,k+l);
340:     }
341:     eps->nconv = k;
342:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
343:   }

345:   if (harmonic) { PetscFree(g); }
346:   /* truncate Schur decomposition and change the state to raw so that
347:      DSVectors() computes eigenvectors from scratch */
348:   DSSetDimensions(eps->ds,eps->nconv,0,0,0);
349:   DSSetState(eps->ds,DS_STATE_RAW);
350:   return(0);
351: }

353: static PetscErrorCode EPSKrylovSchurSetRestart_KrylovSchur(EPS eps,PetscReal keep)
354: {
355:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

358:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
359:   else {
360:     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);
361:     ctx->keep = keep;
362:   }
363:   return(0);
364: }

366: /*@
367:    EPSKrylovSchurSetRestart - Sets the restart parameter for the Krylov-Schur
368:    method, in particular the proportion of basis vectors that must be kept
369:    after restart.

371:    Logically Collective on EPS

373:    Input Parameters:
374: +  eps - the eigenproblem solver context
375: -  keep - the number of vectors to be kept at restart

377:    Options Database Key:
378: .  -eps_krylovschur_restart - Sets the restart parameter

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

383:    Level: advanced

385: .seealso: EPSKrylovSchurGetRestart()
386: @*/
387: PetscErrorCode EPSKrylovSchurSetRestart(EPS eps,PetscReal keep)
388: {

394:   PetscTryMethod(eps,"EPSKrylovSchurSetRestart_C",(EPS,PetscReal),(eps,keep));
395:   return(0);
396: }

398: static PetscErrorCode EPSKrylovSchurGetRestart_KrylovSchur(EPS eps,PetscReal *keep)
399: {
400:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

403:   *keep = ctx->keep;
404:   return(0);
405: }

407: /*@
408:    EPSKrylovSchurGetRestart - Gets the restart parameter used in the
409:    Krylov-Schur method.

411:    Not Collective

413:    Input Parameter:
414: .  eps - the eigenproblem solver context

416:    Output Parameter:
417: .  keep - the restart parameter

419:    Level: advanced

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

430:   PetscUseMethod(eps,"EPSKrylovSchurGetRestart_C",(EPS,PetscReal*),(eps,keep));
431:   return(0);
432: }

434: static PetscErrorCode EPSKrylovSchurSetLocking_KrylovSchur(EPS eps,PetscBool lock)
435: {
436:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

439:   ctx->lock = lock;
440:   return(0);
441: }

443: /*@
444:    EPSKrylovSchurSetLocking - Choose between locking and non-locking variants of
445:    the Krylov-Schur method.

447:    Logically Collective on EPS

449:    Input Parameters:
450: +  eps  - the eigenproblem solver context
451: -  lock - true if the locking variant must be selected

453:    Options Database Key:
454: .  -eps_krylovschur_locking - Sets the locking flag

456:    Notes:
457:    The default is to lock converged eigenpairs when the method restarts.
458:    This behaviour can be changed so that all directions are kept in the
459:    working subspace even if already converged to working accuracy (the
460:    non-locking variant).

462:    Level: advanced

464: .seealso: EPSKrylovSchurGetLocking()
465: @*/
466: PetscErrorCode EPSKrylovSchurSetLocking(EPS eps,PetscBool lock)
467: {

473:   PetscTryMethod(eps,"EPSKrylovSchurSetLocking_C",(EPS,PetscBool),(eps,lock));
474:   return(0);
475: }

477: static PetscErrorCode EPSKrylovSchurGetLocking_KrylovSchur(EPS eps,PetscBool *lock)
478: {
479:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

482:   *lock = ctx->lock;
483:   return(0);
484: }

486: /*@
487:    EPSKrylovSchurGetLocking - Gets the locking flag used in the Krylov-Schur
488:    method.

490:    Not Collective

492:    Input Parameter:
493: .  eps - the eigenproblem solver context

495:    Output Parameter:
496: .  lock - the locking flag

498:    Level: advanced

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

509:   PetscUseMethod(eps,"EPSKrylovSchurGetLocking_C",(EPS,PetscBool*),(eps,lock));
510:   return(0);
511: }

513: static PetscErrorCode EPSKrylovSchurSetPartitions_KrylovSchur(EPS eps,PetscInt npart)
514: {
515:   PetscErrorCode  ierr;
516:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
517:   PetscMPIInt     size;

520:   if (ctx->npart!=npart) {
521:     if (ctx->commset) { PetscSubcommDestroy(&ctx->subc); }
522:     EPSDestroy(&ctx->eps);
523:   }
524:   if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
525:     ctx->npart = 1;
526:   } else {
527:     MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size);
528:     if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
529:     ctx->npart = npart;
530:   }
531:   eps->state = EPS_STATE_INITIAL;
532:   return(0);
533: }

535: /*@
536:    EPSKrylovSchurSetPartitions - Sets the number of partitions for the
537:    case of doing spectrum slicing for a computational interval with the
538:    communicator split in several sub-communicators.

540:    Logically Collective on EPS

542:    Input Parameters:
543: +  eps   - the eigenproblem solver context
544: -  npart - number of partitions

546:    Options Database Key:
547: .  -eps_krylovschur_partitions <npart> - Sets the number of partitions

549:    Notes:
550:    By default, npart=1 so all processes in the communicator participate in
551:    the processing of the whole interval. If npart>1 then the interval is
552:    divided into npart subintervals, each of them being processed by a
553:    subset of processes.

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

558:    Level: advanced

560: .seealso: EPSKrylovSchurSetSubintervals(), EPSSetInterval()
561: @*/
562: PetscErrorCode EPSKrylovSchurSetPartitions(EPS eps,PetscInt npart)
563: {

569:   PetscTryMethod(eps,"EPSKrylovSchurSetPartitions_C",(EPS,PetscInt),(eps,npart));
570:   return(0);
571: }

573: static PetscErrorCode EPSKrylovSchurGetPartitions_KrylovSchur(EPS eps,PetscInt *npart)
574: {
575:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

578:   *npart  = ctx->npart;
579:   return(0);
580: }

582: /*@
583:    EPSKrylovSchurGetPartitions - Gets the number of partitions of the
584:    communicator in case of spectrum slicing.

586:    Not Collective

588:    Input Parameter:
589: .  eps - the eigenproblem solver context

591:    Output Parameter:
592: .  npart - number of partitions

594:    Level: advanced

596: .seealso: EPSKrylovSchurSetPartitions()
597: @*/
598: PetscErrorCode EPSKrylovSchurGetPartitions(EPS eps,PetscInt *npart)
599: {

605:   PetscUseMethod(eps,"EPSKrylovSchurGetPartitions_C",(EPS,PetscInt*),(eps,npart));
606:   return(0);
607: }

609: static PetscErrorCode EPSKrylovSchurSetDetectZeros_KrylovSchur(EPS eps,PetscBool detect)
610: {
611:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

614:   ctx->detect = detect;
615:   eps->state  = EPS_STATE_INITIAL;
616:   return(0);
617: }

619: /*@
620:    EPSKrylovSchurSetDetectZeros - Sets a flag to enforce detection of
621:    zeros during the factorizations throughout the spectrum slicing computation.

623:    Logically Collective on EPS

625:    Input Parameters:
626: +  eps    - the eigenproblem solver context
627: -  detect - check for zeros

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

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

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

641:    Level: advanced

643: .seealso: EPSKrylovSchurSetPartitions(), EPSSetInterval()
644: @*/
645: PetscErrorCode EPSKrylovSchurSetDetectZeros(EPS eps,PetscBool detect)
646: {

652:   PetscTryMethod(eps,"EPSKrylovSchurSetDetectZeros_C",(EPS,PetscBool),(eps,detect));
653:   return(0);
654: }

656: static PetscErrorCode EPSKrylovSchurGetDetectZeros_KrylovSchur(EPS eps,PetscBool *detect)
657: {
658:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

661:   *detect = ctx->detect;
662:   return(0);
663: }

665: /*@
666:    EPSKrylovSchurGetDetectZeros - Gets the flag that enforces zero detection
667:    in spectrum slicing.

669:    Not Collective

671:    Input Parameter:
672: .  eps - the eigenproblem solver context

674:    Output Parameter:
675: .  detect - whether zeros detection is enforced during factorizations

677:    Level: advanced

679: .seealso: EPSKrylovSchurSetDetectZeros()
680: @*/
681: PetscErrorCode EPSKrylovSchurGetDetectZeros(EPS eps,PetscBool *detect)
682: {

688:   PetscUseMethod(eps,"EPSKrylovSchurGetDetectZeros_C",(EPS,PetscBool*),(eps,detect));
689:   return(0);
690: }

692: static PetscErrorCode EPSKrylovSchurSetDimensions_KrylovSchur(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
693: {
694:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

697:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
698:   ctx->nev = nev;
699:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
700:     ctx->ncv = 0;
701:   } else {
702:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
703:     ctx->ncv = ncv;
704:   }
705:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
706:     ctx->mpd = 0;
707:   } else {
708:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
709:     ctx->mpd = mpd;
710:   }
711:   eps->state = EPS_STATE_INITIAL;
712:   return(0);
713: }

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

720:    Logically Collective on EPS

722:    Input Parameters:
723: +  eps - the eigenproblem solver context
724: .  nev - number of eigenvalues to compute
725: .  ncv - the maximum dimension of the subspace to be used by the subsolve
726: -  mpd - the maximum dimension allowed for the projected problem

728:    Options Database Key:
729: +  -eps_krylovschur_nev <nev> - Sets the number of eigenvalues
730: .  -eps_krylovschur_ncv <ncv> - Sets the dimension of the subspace
731: -  -eps_krylovschur_mpd <mpd> - Sets the maximum projected dimension

733:    Level: advanced

735: .seealso: EPSKrylovSchurGetDimensions(), EPSSetDimensions(), EPSSetInterval()
736: @*/
737: PetscErrorCode EPSKrylovSchurSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
738: {

746:   PetscTryMethod(eps,"EPSKrylovSchurSetDimensions_C",(EPS,PetscInt,PetscInt,PetscInt),(eps,nev,ncv,mpd));
747:   return(0);
748: }

750: static PetscErrorCode EPSKrylovSchurGetDimensions_KrylovSchur(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
751: {
752:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

755:   if (nev) *nev = ctx->nev;
756:   if (ncv) *ncv = ctx->ncv;
757:   if (mpd) *mpd = ctx->mpd;
758:   return(0);
759: }

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

765:    Not Collective

767:    Input Parameter:
768: .  eps - the eigenproblem solver context

770:    Output Parameters:
771: +  nev - number of eigenvalues to compute
772: .  ncv - the maximum dimension of the subspace to be used by the subsolve
773: -  mpd - the maximum dimension allowed for the projected problem

775:    Level: advanced

777: .seealso: EPSKrylovSchurSetDimensions()
778: @*/
779: PetscErrorCode EPSKrylovSchurGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
780: {

785:   PetscUseMethod(eps,"EPSKrylovSchurGetDimensions_C",(EPS,PetscInt*,PetscInt*,PetscInt*),(eps,nev,ncv,mpd));
786:   return(0);
787: }

789: static PetscErrorCode EPSKrylovSchurSetSubintervals_KrylovSchur(EPS eps,PetscReal* subint)
790: {
791:   PetscErrorCode  ierr;
792:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
793:   PetscInt        i;

796:   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()");
797:   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");
798:   if (ctx->subintervals) { PetscFree(ctx->subintervals); }
799:   PetscMalloc1(ctx->npart+1,&ctx->subintervals);
800:   for (i=0;i<ctx->npart+1;i++) ctx->subintervals[i] = subint[i];
801:   ctx->subintset = PETSC_TRUE;
802:   eps->state = EPS_STATE_INITIAL;
803:   return(0);
804: }

806: /*@C
807:    EPSKrylovSchurSetSubintervals - Sets the points that delimit the
808:    subintervals to be used in spectrum slicing with several partitions.

810:    Logically Collective on EPS

812:    Input Parameters:
813: +  eps    - the eigenproblem solver context
814: -  subint - array of real values specifying subintervals

816:    Notes:
817:    This function must be called after EPSKrylovSchurSetPartitions(). For npart
818:    partitions, the argument subint must contain npart+1 real values sorted in
819:    ascending order: subint_0, subint_1, ..., subint_npart, where the first
820:    and last values must coincide with the interval endpoints set with
821:    EPSSetInterval().

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

826:    Level: advanced

828: .seealso: EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubintervals(), EPSSetInterval()
829: @*/
830: PetscErrorCode EPSKrylovSchurSetSubintervals(EPS eps,PetscReal *subint)
831: {

836:   PetscTryMethod(eps,"EPSKrylovSchurSetSubintervals_C",(EPS,PetscReal*),(eps,subint));
837:   return(0);
838: }

840: static PetscErrorCode EPSKrylovSchurGetSubintervals_KrylovSchur(EPS eps,PetscReal **subint)
841: {
842:   PetscErrorCode  ierr;
843:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
844:   PetscInt        i;

847:   if (!ctx->subintset) {
848:     if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
849:     if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
850:   }
851:   PetscMalloc1(ctx->npart+1,subint);
852:   for (i=0;i<=ctx->npart;i++) (*subint)[i] = ctx->subintervals[i];
853:   return(0);
854: }

856: /*@C
857:    EPSKrylovSchurGetSubintervals - Returns the points that delimit the
858:    subintervals used in spectrum slicing with several partitions.

860:    Logically Collective on EPS

862:    Input Parameter:
863: .  eps    - the eigenproblem solver context

865:    Output Parameter:
866: .  subint - array of real values specifying subintervals

868:    Notes:
869:    If the user passed values with EPSKrylovSchurSetSubintervals(), then the
870:    same values are returned. Otherwise, the values computed internally are
871:    obtained.

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

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

878:    Fortran Notes:
879:    The calling sequence from Fortran is
880: .vb
881:    EPSKrylovSchurGetSubintervals(eps,subint,ierr)
882:    double precision subint(npart+1) output
883: .ve

885:    Level: advanced

887: .seealso: EPSKrylovSchurSetSubintervals(), EPSKrylovSchurGetPartitions(), EPSSetInterval()
888: @*/
889: PetscErrorCode EPSKrylovSchurGetSubintervals(EPS eps,PetscReal **subint)
890: {

896:   PetscUseMethod(eps,"EPSKrylovSchurGetSubintervals_C",(EPS,PetscReal**),(eps,subint));
897:   return(0);
898: }

900: static PetscErrorCode EPSKrylovSchurGetInertias_KrylovSchur(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
901: {
902:   PetscErrorCode  ierr;
903:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
904:   PetscInt        i,numsh;
905:   EPS_SR          sr = ctx->sr;

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

946: /*@C
947:    EPSKrylovSchurGetInertias - Gets the values of the shifts and their
948:    corresponding inertias in case of doing spectrum slicing for a
949:    computational interval.

951:    Not Collective

953:    Input Parameter:
954: .  eps - the eigenproblem solver context

956:    Output Parameters:
957: +  n        - number of shifts, including the endpoints of the interval
958: .  shifts   - the values of the shifts used internally in the solver
959: -  inertias - the values of the inertia in each shift

961:    Notes:
962:    If called after EPSSolve(), all shifts used internally by the solver are
963:    returned (including both endpoints and any intermediate ones). If called
964:    before EPSSolve() and after EPSSetUp() then only the information of the
965:    endpoints of subintervals is available.

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

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

972:    Fortran Notes:
973:    The calling sequence from Fortran is
974: .vb
975:    EPSKrylovSchurGetInertias(eps,n,shifts,inertias,ierr)
976:    integer n
977:    double precision shifts(*)
978:    integer inertias(*)
979: .ve
980:    The arrays should be at least of length n. The value of n can be determined
981:    by an initial call
982: .vb
983:    EPSKrylovSchurGetInertias(eps,n,PETSC_NULL_REAL,PETSC_NULL_INTEGER,ierr)
984: .ve

986:    Level: advanced

988: .seealso: EPSSetInterval(), EPSKrylovSchurSetSubintervals()
989: @*/
990: PetscErrorCode EPSKrylovSchurGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
991: {

997:   PetscUseMethod(eps,"EPSKrylovSchurGetInertias_C",(EPS,PetscInt*,PetscReal**,PetscInt**),(eps,n,shifts,inertias));
998:   return(0);
999: }

1001: static PetscErrorCode EPSKrylovSchurGetSubcommInfo_KrylovSchur(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
1002: {
1003:   PetscErrorCode  ierr;
1004:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1005:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

1008:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1009:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1010:   if (k) *k = (ctx->npart==1)? 0: ctx->subc->color;
1011:   if (n) *n = sr->numEigs;
1012:   if (v) {
1013:     BVCreateVec(sr->V,v);
1014:   }
1015:   return(0);
1016: }

1018: /*@C
1019:    EPSKrylovSchurGetSubcommInfo - Gets information related to the case of
1020:    doing spectrum slicing for a computational interval with multiple
1021:    communicators.

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

1025:    Input Parameter:
1026: .  eps - the eigenproblem solver context

1028:    Output Parameters:
1029: +  k - index of the subinterval for the calling process
1030: .  n - number of eigenvalues found in the k-th subinterval
1031: -  v - a vector owned by processes in the subcommunicator with dimensions
1032:        compatible for locally computed eigenvectors (or NULL)

1034:    Notes:
1035:    This function is only available for spectrum slicing runs.

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

1039:    Level: advanced

1041: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommPairs()
1042: @*/
1043: PetscErrorCode EPSKrylovSchurGetSubcommInfo(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
1044: {

1049:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommInfo_C",(EPS,PetscInt*,PetscInt*,Vec*),(eps,k,n,v));
1050:   return(0);
1051: }

1053: static PetscErrorCode EPSKrylovSchurGetSubcommPairs_KrylovSchur(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1054: {
1055:   PetscErrorCode  ierr;
1056:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1057:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

1060:   EPSCheckSolved(eps,1);
1061:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1062:   if (i<0 || i>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
1063:   if (eig) *eig = sr->eigr[sr->perm[i]];
1064:   if (v) { BVCopyVec(sr->V,sr->perm[i],v); }
1065:   return(0);
1066: }

1068: /*@C
1069:    EPSKrylovSchurGetSubcommPairs - Gets the i-th eigenpair stored
1070:    internally in the subcommunicator to which the calling process belongs.

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

1074:    Input Parameter:
1075: +  eps - the eigenproblem solver context
1076: -  i   - index of the solution

1078:    Output Parameters:
1079: +  eig - the eigenvalue
1080: -  v   - the eigenvector

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

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

1090:    Level: advanced

1092: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo(), EPSKrylovSchurGetSubcommMats()
1093: @*/
1094: PetscErrorCode EPSKrylovSchurGetSubcommPairs(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1095: {

1101:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommPairs_C",(EPS,PetscInt,PetscScalar*,Vec),(eps,i,eig,v));
1102:   return(0);
1103: }

1105: static PetscErrorCode EPSKrylovSchurGetSubcommMats_KrylovSchur(EPS eps,Mat *A,Mat *B)
1106: {
1107:   PetscErrorCode  ierr;
1108:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

1111:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1112:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1113:   EPSGetOperators(ctx->eps,A,B);
1114:   return(0);
1115: }

1117: /*@C
1118:    EPSKrylovSchurGetSubcommMats - Gets the eigenproblem matrices stored
1119:    internally in the subcommunicator to which the calling process belongs.

1121:    Collective on the subcommunicator

1123:    Input Parameter:
1124: .  eps - the eigenproblem solver context

1126:    Output Parameters:
1127: +  A  - the matrix associated with the eigensystem
1128: -  B  - the second matrix in the case of generalized eigenproblems

1130:    Notes:
1131:    This is the analog of EPSGetOperators(), but returns the matrices distributed
1132:    differently (in the subcommunicator rather than in the parent communicator).

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

1136:    Level: advanced

1138: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo()
1139: @*/
1140: PetscErrorCode EPSKrylovSchurGetSubcommMats(EPS eps,Mat *A,Mat *B)
1141: {

1146:   PetscTryMethod(eps,"EPSKrylovSchurGetSubcommMats_C",(EPS,Mat*,Mat*),(eps,A,B));
1147:   return(0);
1148: }

1150: static PetscErrorCode EPSKrylovSchurUpdateSubcommMats_KrylovSchur(EPS eps,PetscScalar a,PetscScalar ap,Mat Au,PetscScalar b,PetscScalar bp, Mat Bu,MatStructure str,PetscBool globalup)
1151: {
1152:   PetscErrorCode  ierr;
1153:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data,*subctx;
1154:   Mat             A,B=NULL,Ag,Bg=NULL;
1155:   PetscBool       reuse=PETSC_TRUE;

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

1163:   MatScale(A,a);
1164:   if (Au) {
1165:     MatAXPY(A,ap,Au,str);
1166:   }
1167:   if (B) MatScale(B,b);
1168:   if (Bu) {
1169:     MatAXPY(B,bp,Bu,str);
1170:   }
1171:   EPSSetOperators(ctx->eps,A,B);

1173:   /* Update stored matrix state */
1174:   subctx = (EPS_KRYLOVSCHUR*)ctx->eps->data;
1175:   PetscObjectStateGet((PetscObject)A,&subctx->Astate);
1176:   if (B) { PetscObjectStateGet((PetscObject)B,&subctx->Bstate); }

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

1206: /*@C
1207:    EPSKrylovSchurUpdateSubcommMats - Update the eigenproblem matrices stored
1208:    internally in the subcommunicator to which the calling process belongs.

1210:    Collective on EPS

1212:    Input Parameters:
1213: +  eps - the eigenproblem solver context
1214: .  s   - scalar that multiplies the existing A matrix
1215: .  a   - scalar used in the axpy operation on A
1216: .  Au  - matrix used in the axpy operation on A
1217: .  t   - scalar that multiplies the existing B matrix
1218: .  b   - scalar used in the axpy operation on B
1219: .  Bu  - matrix used in the axpy operation on B
1220: .  str - structure flag
1221: -  globalup - flag indicating if global matrices must be updated

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

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

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

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

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

1240:    Level: advanced

1242: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommMats()
1243: @*/
1244: PetscErrorCode EPSKrylovSchurUpdateSubcommMats(EPS eps,PetscScalar s,PetscScalar a,Mat Au,PetscScalar t,PetscScalar b,Mat Bu,MatStructure str,PetscBool globalup)
1245: {

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

1262: PetscErrorCode EPSSetFromOptions_KrylovSchur(PetscOptionItems *PetscOptionsObject,EPS eps)
1263: {
1264:   PetscErrorCode  ierr;
1265:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1266:   PetscBool       flg,lock,b,f1,f2,f3;
1267:   PetscReal       keep;
1268:   PetscInt        i,j,k;

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

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

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

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

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

1287:     i = 1;
1288:     j = k = PETSC_DECIDE;
1289:     PetscOptionsInt("-eps_krylovschur_nev","Number of eigenvalues to compute in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",40,&i,&f1);
1290:     PetscOptionsInt("-eps_krylovschur_ncv","Number of basis vectors in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&j,&f2);
1291:     PetscOptionsInt("-eps_krylovschur_mpd","Maximum dimension of projected problem in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&k,&f3);
1292:     if (f1 || f2 || f3) { EPSKrylovSchurSetDimensions(eps,i,j,k); }

1294:   PetscOptionsTail();
1295:   return(0);
1296: }

1298: PetscErrorCode EPSView_KrylovSchur(EPS eps,PetscViewer viewer)
1299: {
1300:   PetscErrorCode  ierr;
1301:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1302:   PetscBool       isascii,isfilt;

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

1325: PetscErrorCode EPSDestroy_KrylovSchur(EPS eps)
1326: {

1330:   PetscFree(eps->data);
1331:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",NULL);
1332:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",NULL);
1333:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetLocking_C",NULL);
1334:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetLocking_C",NULL);
1335:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetPartitions_C",NULL);
1336:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetPartitions_C",NULL);
1337:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDetectZeros_C",NULL);
1338:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDetectZeros_C",NULL);
1339:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",NULL);
1340:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",NULL);
1341:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetSubintervals_C",NULL);
1342:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubintervals_C",NULL);
1343:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",NULL);
1344:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommInfo_C",NULL);
1345:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommPairs_C",NULL);
1346:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommMats_C",NULL);
1347:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurUpdateSubcommMats_C",NULL);
1348:   return(0);
1349: }

1351: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
1352: {
1354:   PetscBool      isfilt;

1357:   PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
1358:   if (eps->which==EPS_ALL && !isfilt) {
1359:     EPSReset_KrylovSchur_Slice(eps);
1360:   }
1361:   return(0);
1362: }

1364: PetscErrorCode EPSSetDefaultST_KrylovSchur(EPS eps)
1365: {

1369:   if (eps->which==EPS_ALL) {
1370:     if (!((PetscObject)eps->st)->type_name) {
1371:       STSetType(eps->st,STSINVERT);
1372:     }
1373:   }
1374:   return(0);
1375: }

1377: PETSC_EXTERN PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
1378: {
1379:   EPS_KRYLOVSCHUR *ctx;
1380:   PetscErrorCode  ierr;

1383:   PetscNewLog(eps,&ctx);
1384:   eps->data   = (void*)ctx;
1385:   ctx->lock   = PETSC_TRUE;
1386:   ctx->nev    = 1;
1387:   ctx->npart  = 1;
1388:   ctx->detect = PETSC_FALSE;
1389:   ctx->global = PETSC_TRUE;

1391:   eps->useds = PETSC_TRUE;

1393:   /* solve and computevectors determined at setup */
1394:   eps->ops->setup          = EPSSetUp_KrylovSchur;
1395:   eps->ops->setfromoptions = EPSSetFromOptions_KrylovSchur;
1396:   eps->ops->destroy        = EPSDestroy_KrylovSchur;
1397:   eps->ops->reset          = EPSReset_KrylovSchur;
1398:   eps->ops->view           = EPSView_KrylovSchur;
1399:   eps->ops->backtransform  = EPSBackTransform_Default;
1400:   eps->ops->setdefaultst   = EPSSetDefaultST_KrylovSchur;

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