Actual source code: krylovschur.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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,EPS_KS_TWOSIDED } 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 if (eps->twosided) {
189:     variant = EPS_KS_TWOSIDED;
190:   } else {
191:     switch (eps->extraction) {
192:       case EPS_RITZ:     variant = EPS_KS_DEFAULT; break;
193:       case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
194:       default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
195:     }
196:   }
197:   switch (variant) {
198:     case EPS_KS_DEFAULT:
199:       eps->ops->solve = EPSSolve_KrylovSchur_Default;
200:       eps->ops->computevectors = EPSComputeVectors_Schur;
201:       DSSetType(eps->ds,DSNHEP);
202:       DSAllocate(eps->ds,eps->ncv+1);
203:       break;
204:     case EPS_KS_SYMM:
205:     case EPS_KS_FILTER:
206:       eps->ops->solve = EPSSolve_KrylovSchur_Symm;
207:       eps->ops->computevectors = EPSComputeVectors_Hermitian;
208:       DSSetType(eps->ds,DSHEP);
209:       DSSetCompact(eps->ds,PETSC_TRUE);
210:       DSSetExtraRow(eps->ds,PETSC_TRUE);
211:       DSAllocate(eps->ds,eps->ncv+1);
212:       break;
213:     case EPS_KS_SLICE:
214:       if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing does not support user-defined stopping test");
215:       eps->ops->solve = EPSSolve_KrylovSchur_Slice;
216:       eps->ops->computevectors = EPSComputeVectors_Slice;
217:       break;
218:     case EPS_KS_INDEF:
219:       eps->ops->solve = EPSSolve_KrylovSchur_Indefinite;
220:       eps->ops->computevectors = EPSComputeVectors_Indefinite;
221:       DSSetType(eps->ds,DSGHIEP);
222:       DSSetCompact(eps->ds,PETSC_TRUE);
223:       DSAllocate(eps->ds,eps->ncv+1);
224:       /* force reorthogonalization for pseudo-Lanczos */
225:       BVGetOrthogonalization(eps->V,&otype,NULL,&eta,&obtype);
226:       BVSetOrthogonalization(eps->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
227:       break;
228:     case EPS_KS_TWOSIDED:
229:       eps->ops->solve = EPSSolve_KrylovSchur_TwoSided;
230:       eps->ops->computevectors = EPSComputeVectors_Schur;
231:       DSSetType(eps->ds,DSNHEP);
232:       DSAllocate(eps->ds,eps->ncv+1);
233:       DSSetType(eps->dsts,DSNHEP);
234:       DSAllocate(eps->dsts,eps->ncv+1);
235:       break;
236:     default: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error");
237:   }
238:   return(0);
239: }

241: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
242: {
243:   PetscErrorCode  ierr;
244:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
245:   PetscInt        i,j,*pj,k,l,nv,ld,nconv;
246:   Mat             U;
247:   PetscScalar     *S,*Q,*g;
248:   PetscReal       beta,gamma=1.0;
249:   PetscBool       breakdown,harmonic;

252:   DSGetLeadingDimension(eps->ds,&ld);
253:   harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
254:   if (harmonic) { PetscMalloc1(ld,&g); }
255:   if (eps->arbitrary) pj = &j;
256:   else pj = NULL;

258:   /* Get the starting Arnoldi vector */
259:   EPSGetStartVector(eps,0,NULL);
260:   l = 0;

262:   /* Restart loop */
263:   while (eps->reason == EPS_CONVERGED_ITERATING) {
264:     eps->its++;

266:     /* Compute an nv-step Arnoldi factorization */
267:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
268:     DSGetArray(eps->ds,DS_MAT_A,&S);
269:     EPSBasicArnoldi(eps,PETSC_FALSE,S,ld,eps->nconv+l,&nv,&beta,&breakdown);
270:     DSRestoreArray(eps->ds,DS_MAT_A,&S);
271:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
272:     if (l==0) {
273:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
274:     } else {
275:       DSSetState(eps->ds,DS_STATE_RAW);
276:     }
277:     BVSetActiveColumns(eps->V,eps->nconv,nv);

279:     /* Compute translation of Krylov decomposition if harmonic extraction used */
280:     if (harmonic) {
281:       DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,g,&gamma);
282:     }

284:     /* Solve projected problem */
285:     DSSolve(eps->ds,eps->eigr,eps->eigi);
286:     if (eps->arbitrary) {
287:       EPSGetArbitraryValues(eps,eps->rr,eps->ri);
288:       j=1;
289:     }
290:     DSSort(eps->ds,eps->eigr,eps->eigi,eps->rr,eps->ri,pj);
291:     DSSynchronize(eps->ds,eps->eigr,eps->eigi);

293:     /* Check convergence */
294:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,0.0,gamma,&k);
295:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
296:     nconv = k;

298:     /* Update l */
299:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
300:     else {
301:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
302: #if !defined(PETSC_USE_COMPLEX)
303:       DSGetArray(eps->ds,DS_MAT_A,&S);
304:       if (S[k+l+(k+l-1)*ld] != 0.0) {
305:         if (k+l<nv-1) l = l+1;
306:         else l = l-1;
307:       }
308:       DSRestoreArray(eps->ds,DS_MAT_A,&S);
309: #endif
310:     }
311:     if ((!ctx->lock || harmonic) && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */

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

348:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
349:       BVCopyColumn(eps->V,nv,k+l);
350:     }
351:     eps->nconv = k;
352:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
353:   }

355:   if (harmonic) { PetscFree(g); }
356:   /* truncate Schur decomposition and change the state to raw so that
357:      DSVectors() computes eigenvectors from scratch */
358:   DSSetDimensions(eps->ds,eps->nconv,0,0,0);
359:   DSSetState(eps->ds,DS_STATE_RAW);
360:   return(0);
361: }

363: static PetscErrorCode EPSKrylovSchurSetRestart_KrylovSchur(EPS eps,PetscReal keep)
364: {
365:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

368:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
369:   else {
370:     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);
371:     ctx->keep = keep;
372:   }
373:   return(0);
374: }

376: /*@
377:    EPSKrylovSchurSetRestart - Sets the restart parameter for the Krylov-Schur
378:    method, in particular the proportion of basis vectors that must be kept
379:    after restart.

381:    Logically Collective on EPS

383:    Input Parameters:
384: +  eps - the eigenproblem solver context
385: -  keep - the number of vectors to be kept at restart

387:    Options Database Key:
388: .  -eps_krylovschur_restart - Sets the restart parameter

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

393:    Level: advanced

395: .seealso: EPSKrylovSchurGetRestart()
396: @*/
397: PetscErrorCode EPSKrylovSchurSetRestart(EPS eps,PetscReal keep)
398: {

404:   PetscTryMethod(eps,"EPSKrylovSchurSetRestart_C",(EPS,PetscReal),(eps,keep));
405:   return(0);
406: }

408: static PetscErrorCode EPSKrylovSchurGetRestart_KrylovSchur(EPS eps,PetscReal *keep)
409: {
410:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

413:   *keep = ctx->keep;
414:   return(0);
415: }

417: /*@
418:    EPSKrylovSchurGetRestart - Gets the restart parameter used in the
419:    Krylov-Schur method.

421:    Not Collective

423:    Input Parameter:
424: .  eps - the eigenproblem solver context

426:    Output Parameter:
427: .  keep - the restart parameter

429:    Level: advanced

431: .seealso: EPSKrylovSchurSetRestart()
432: @*/
433: PetscErrorCode EPSKrylovSchurGetRestart(EPS eps,PetscReal *keep)
434: {

440:   PetscUseMethod(eps,"EPSKrylovSchurGetRestart_C",(EPS,PetscReal*),(eps,keep));
441:   return(0);
442: }

444: static PetscErrorCode EPSKrylovSchurSetLocking_KrylovSchur(EPS eps,PetscBool lock)
445: {
446:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

449:   ctx->lock = lock;
450:   return(0);
451: }

453: /*@
454:    EPSKrylovSchurSetLocking - Choose between locking and non-locking variants of
455:    the Krylov-Schur method.

457:    Logically Collective on EPS

459:    Input Parameters:
460: +  eps  - the eigenproblem solver context
461: -  lock - true if the locking variant must be selected

463:    Options Database Key:
464: .  -eps_krylovschur_locking - Sets the locking flag

466:    Notes:
467:    The default is to lock converged eigenpairs when the method restarts.
468:    This behaviour can be changed so that all directions are kept in the
469:    working subspace even if already converged to working accuracy (the
470:    non-locking variant).

472:    Level: advanced

474: .seealso: EPSKrylovSchurGetLocking()
475: @*/
476: PetscErrorCode EPSKrylovSchurSetLocking(EPS eps,PetscBool lock)
477: {

483:   PetscTryMethod(eps,"EPSKrylovSchurSetLocking_C",(EPS,PetscBool),(eps,lock));
484:   return(0);
485: }

487: static PetscErrorCode EPSKrylovSchurGetLocking_KrylovSchur(EPS eps,PetscBool *lock)
488: {
489:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

492:   *lock = ctx->lock;
493:   return(0);
494: }

496: /*@
497:    EPSKrylovSchurGetLocking - Gets the locking flag used in the Krylov-Schur
498:    method.

500:    Not Collective

502:    Input Parameter:
503: .  eps - the eigenproblem solver context

505:    Output Parameter:
506: .  lock - the locking flag

508:    Level: advanced

510: .seealso: EPSKrylovSchurSetLocking()
511: @*/
512: PetscErrorCode EPSKrylovSchurGetLocking(EPS eps,PetscBool *lock)
513: {

519:   PetscUseMethod(eps,"EPSKrylovSchurGetLocking_C",(EPS,PetscBool*),(eps,lock));
520:   return(0);
521: }

523: static PetscErrorCode EPSKrylovSchurSetPartitions_KrylovSchur(EPS eps,PetscInt npart)
524: {
525:   PetscErrorCode  ierr;
526:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
527:   PetscMPIInt     size;

530:   if (ctx->npart!=npart) {
531:     if (ctx->commset) { PetscSubcommDestroy(&ctx->subc); }
532:     EPSDestroy(&ctx->eps);
533:   }
534:   if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
535:     ctx->npart = 1;
536:   } else {
537:     MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size);
538:     if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
539:     ctx->npart = npart;
540:   }
541:   eps->state = EPS_STATE_INITIAL;
542:   return(0);
543: }

545: /*@
546:    EPSKrylovSchurSetPartitions - Sets the number of partitions for the
547:    case of doing spectrum slicing for a computational interval with the
548:    communicator split in several sub-communicators.

550:    Logically Collective on EPS

552:    Input Parameters:
553: +  eps   - the eigenproblem solver context
554: -  npart - number of partitions

556:    Options Database Key:
557: .  -eps_krylovschur_partitions <npart> - Sets the number of partitions

559:    Notes:
560:    By default, npart=1 so all processes in the communicator participate in
561:    the processing of the whole interval. If npart>1 then the interval is
562:    divided into npart subintervals, each of them being processed by a
563:    subset of processes.

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

568:    Level: advanced

570: .seealso: EPSKrylovSchurSetSubintervals(), EPSSetInterval()
571: @*/
572: PetscErrorCode EPSKrylovSchurSetPartitions(EPS eps,PetscInt npart)
573: {

579:   PetscTryMethod(eps,"EPSKrylovSchurSetPartitions_C",(EPS,PetscInt),(eps,npart));
580:   return(0);
581: }

583: static PetscErrorCode EPSKrylovSchurGetPartitions_KrylovSchur(EPS eps,PetscInt *npart)
584: {
585:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

588:   *npart  = ctx->npart;
589:   return(0);
590: }

592: /*@
593:    EPSKrylovSchurGetPartitions - Gets the number of partitions of the
594:    communicator in case of spectrum slicing.

596:    Not Collective

598:    Input Parameter:
599: .  eps - the eigenproblem solver context

601:    Output Parameter:
602: .  npart - number of partitions

604:    Level: advanced

606: .seealso: EPSKrylovSchurSetPartitions()
607: @*/
608: PetscErrorCode EPSKrylovSchurGetPartitions(EPS eps,PetscInt *npart)
609: {

615:   PetscUseMethod(eps,"EPSKrylovSchurGetPartitions_C",(EPS,PetscInt*),(eps,npart));
616:   return(0);
617: }

619: static PetscErrorCode EPSKrylovSchurSetDetectZeros_KrylovSchur(EPS eps,PetscBool detect)
620: {
621:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

624:   ctx->detect = detect;
625:   eps->state  = EPS_STATE_INITIAL;
626:   return(0);
627: }

629: /*@
630:    EPSKrylovSchurSetDetectZeros - Sets a flag to enforce detection of
631:    zeros during the factorizations throughout the spectrum slicing computation.

633:    Logically Collective on EPS

635:    Input Parameters:
636: +  eps    - the eigenproblem solver context
637: -  detect - check for zeros

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

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

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

651:    Level: advanced

653: .seealso: EPSKrylovSchurSetPartitions(), EPSSetInterval()
654: @*/
655: PetscErrorCode EPSKrylovSchurSetDetectZeros(EPS eps,PetscBool detect)
656: {

662:   PetscTryMethod(eps,"EPSKrylovSchurSetDetectZeros_C",(EPS,PetscBool),(eps,detect));
663:   return(0);
664: }

666: static PetscErrorCode EPSKrylovSchurGetDetectZeros_KrylovSchur(EPS eps,PetscBool *detect)
667: {
668:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

671:   *detect = ctx->detect;
672:   return(0);
673: }

675: /*@
676:    EPSKrylovSchurGetDetectZeros - Gets the flag that enforces zero detection
677:    in spectrum slicing.

679:    Not Collective

681:    Input Parameter:
682: .  eps - the eigenproblem solver context

684:    Output Parameter:
685: .  detect - whether zeros detection is enforced during factorizations

687:    Level: advanced

689: .seealso: EPSKrylovSchurSetDetectZeros()
690: @*/
691: PetscErrorCode EPSKrylovSchurGetDetectZeros(EPS eps,PetscBool *detect)
692: {

698:   PetscUseMethod(eps,"EPSKrylovSchurGetDetectZeros_C",(EPS,PetscBool*),(eps,detect));
699:   return(0);
700: }

702: static PetscErrorCode EPSKrylovSchurSetDimensions_KrylovSchur(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
703: {
704:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

707:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
708:   ctx->nev = nev;
709:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
710:     ctx->ncv = 0;
711:   } else {
712:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
713:     ctx->ncv = ncv;
714:   }
715:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
716:     ctx->mpd = 0;
717:   } else {
718:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
719:     ctx->mpd = mpd;
720:   }
721:   eps->state = EPS_STATE_INITIAL;
722:   return(0);
723: }

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

730:    Logically Collective on EPS

732:    Input Parameters:
733: +  eps - the eigenproblem solver context
734: .  nev - number of eigenvalues to compute
735: .  ncv - the maximum dimension of the subspace to be used by the subsolve
736: -  mpd - the maximum dimension allowed for the projected problem

738:    Options Database Key:
739: +  -eps_krylovschur_nev <nev> - Sets the number of eigenvalues
740: .  -eps_krylovschur_ncv <ncv> - Sets the dimension of the subspace
741: -  -eps_krylovschur_mpd <mpd> - Sets the maximum projected dimension

743:    Level: advanced

745: .seealso: EPSKrylovSchurGetDimensions(), EPSSetDimensions(), EPSSetInterval()
746: @*/
747: PetscErrorCode EPSKrylovSchurSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
748: {

756:   PetscTryMethod(eps,"EPSKrylovSchurSetDimensions_C",(EPS,PetscInt,PetscInt,PetscInt),(eps,nev,ncv,mpd));
757:   return(0);
758: }

760: static PetscErrorCode EPSKrylovSchurGetDimensions_KrylovSchur(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
761: {
762:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

765:   if (nev) *nev = ctx->nev;
766:   if (ncv) *ncv = ctx->ncv;
767:   if (mpd) *mpd = ctx->mpd;
768:   return(0);
769: }

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

775:    Not Collective

777:    Input Parameter:
778: .  eps - the eigenproblem solver context

780:    Output Parameters:
781: +  nev - number of eigenvalues to compute
782: .  ncv - the maximum dimension of the subspace to be used by the subsolve
783: -  mpd - the maximum dimension allowed for the projected problem

785:    Level: advanced

787: .seealso: EPSKrylovSchurSetDimensions()
788: @*/
789: PetscErrorCode EPSKrylovSchurGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
790: {

795:   PetscUseMethod(eps,"EPSKrylovSchurGetDimensions_C",(EPS,PetscInt*,PetscInt*,PetscInt*),(eps,nev,ncv,mpd));
796:   return(0);
797: }

799: static PetscErrorCode EPSKrylovSchurSetSubintervals_KrylovSchur(EPS eps,PetscReal* subint)
800: {
801:   PetscErrorCode  ierr;
802:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
803:   PetscInt        i;

806:   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()");
807:   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");
808:   if (ctx->subintervals) { PetscFree(ctx->subintervals); }
809:   PetscMalloc1(ctx->npart+1,&ctx->subintervals);
810:   for (i=0;i<ctx->npart+1;i++) ctx->subintervals[i] = subint[i];
811:   ctx->subintset = PETSC_TRUE;
812:   eps->state = EPS_STATE_INITIAL;
813:   return(0);
814: }

816: /*@C
817:    EPSKrylovSchurSetSubintervals - Sets the points that delimit the
818:    subintervals to be used in spectrum slicing with several partitions.

820:    Logically Collective on EPS

822:    Input Parameters:
823: +  eps    - the eigenproblem solver context
824: -  subint - array of real values specifying subintervals

826:    Notes:
827:    This function must be called after EPSKrylovSchurSetPartitions(). For npart
828:    partitions, the argument subint must contain npart+1 real values sorted in
829:    ascending order: subint_0, subint_1, ..., subint_npart, where the first
830:    and last values must coincide with the interval endpoints set with
831:    EPSSetInterval().

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

836:    Level: advanced

838: .seealso: EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubintervals(), EPSSetInterval()
839: @*/
840: PetscErrorCode EPSKrylovSchurSetSubintervals(EPS eps,PetscReal *subint)
841: {

846:   PetscTryMethod(eps,"EPSKrylovSchurSetSubintervals_C",(EPS,PetscReal*),(eps,subint));
847:   return(0);
848: }

850: static PetscErrorCode EPSKrylovSchurGetSubintervals_KrylovSchur(EPS eps,PetscReal **subint)
851: {
852:   PetscErrorCode  ierr;
853:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
854:   PetscInt        i;

857:   if (!ctx->subintset) {
858:     if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
859:     if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
860:   }
861:   PetscMalloc1(ctx->npart+1,subint);
862:   for (i=0;i<=ctx->npart;i++) (*subint)[i] = ctx->subintervals[i];
863:   return(0);
864: }

866: /*@C
867:    EPSKrylovSchurGetSubintervals - Returns the points that delimit the
868:    subintervals used in spectrum slicing with several partitions.

870:    Logically Collective on EPS

872:    Input Parameter:
873: .  eps    - the eigenproblem solver context

875:    Output Parameter:
876: .  subint - array of real values specifying subintervals

878:    Notes:
879:    If the user passed values with EPSKrylovSchurSetSubintervals(), then the
880:    same values are returned. Otherwise, the values computed internally are
881:    obtained.

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

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

888:    Fortran Notes:
889:    The calling sequence from Fortran is
890: .vb
891:    EPSKrylovSchurGetSubintervals(eps,subint,ierr)
892:    double precision subint(npart+1) output
893: .ve

895:    Level: advanced

897: .seealso: EPSKrylovSchurSetSubintervals(), EPSKrylovSchurGetPartitions(), EPSSetInterval()
898: @*/
899: PetscErrorCode EPSKrylovSchurGetSubintervals(EPS eps,PetscReal **subint)
900: {

906:   PetscUseMethod(eps,"EPSKrylovSchurGetSubintervals_C",(EPS,PetscReal**),(eps,subint));
907:   return(0);
908: }

910: static PetscErrorCode EPSKrylovSchurGetInertias_KrylovSchur(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
911: {
912:   PetscErrorCode  ierr;
913:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
914:   PetscInt        i,numsh;
915:   EPS_SR          sr = ctx->sr;

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

956: /*@C
957:    EPSKrylovSchurGetInertias - Gets the values of the shifts and their
958:    corresponding inertias in case of doing spectrum slicing for a
959:    computational interval.

961:    Not Collective

963:    Input Parameter:
964: .  eps - the eigenproblem solver context

966:    Output Parameters:
967: +  n        - number of shifts, including the endpoints of the interval
968: .  shifts   - the values of the shifts used internally in the solver
969: -  inertias - the values of the inertia in each shift

971:    Notes:
972:    If called after EPSSolve(), all shifts used internally by the solver are
973:    returned (including both endpoints and any intermediate ones). If called
974:    before EPSSolve() and after EPSSetUp() then only the information of the
975:    endpoints of subintervals is available.

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

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

982:    Fortran Notes:
983:    The calling sequence from Fortran is
984: .vb
985:    EPSKrylovSchurGetInertias(eps,n,shifts,inertias,ierr)
986:    integer n
987:    double precision shifts(*)
988:    integer inertias(*)
989: .ve
990:    The arrays should be at least of length n. The value of n can be determined
991:    by an initial call
992: .vb
993:    EPSKrylovSchurGetInertias(eps,n,PETSC_NULL_REAL,PETSC_NULL_INTEGER,ierr)
994: .ve

996:    Level: advanced

998: .seealso: EPSSetInterval(), EPSKrylovSchurSetSubintervals()
999: @*/
1000: PetscErrorCode EPSKrylovSchurGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
1001: {

1007:   PetscUseMethod(eps,"EPSKrylovSchurGetInertias_C",(EPS,PetscInt*,PetscReal**,PetscInt**),(eps,n,shifts,inertias));
1008:   return(0);
1009: }

1011: static PetscErrorCode EPSKrylovSchurGetSubcommInfo_KrylovSchur(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
1012: {
1013:   PetscErrorCode  ierr;
1014:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1015:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

1018:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1019:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1020:   if (k) *k = (ctx->npart==1)? 0: ctx->subc->color;
1021:   if (n) *n = sr->numEigs;
1022:   if (v) {
1023:     BVCreateVec(sr->V,v);
1024:   }
1025:   return(0);
1026: }

1028: /*@C
1029:    EPSKrylovSchurGetSubcommInfo - Gets information related to the case of
1030:    doing spectrum slicing for a computational interval with multiple
1031:    communicators.

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

1035:    Input Parameter:
1036: .  eps - the eigenproblem solver context

1038:    Output Parameters:
1039: +  k - index of the subinterval for the calling process
1040: .  n - number of eigenvalues found in the k-th subinterval
1041: -  v - a vector owned by processes in the subcommunicator with dimensions
1042:        compatible for locally computed eigenvectors (or NULL)

1044:    Notes:
1045:    This function is only available for spectrum slicing runs.

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

1049:    Level: advanced

1051: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommPairs()
1052: @*/
1053: PetscErrorCode EPSKrylovSchurGetSubcommInfo(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
1054: {

1059:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommInfo_C",(EPS,PetscInt*,PetscInt*,Vec*),(eps,k,n,v));
1060:   return(0);
1061: }

1063: static PetscErrorCode EPSKrylovSchurGetSubcommPairs_KrylovSchur(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1064: {
1065:   PetscErrorCode  ierr;
1066:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1067:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

1070:   EPSCheckSolved(eps,1);
1071:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1072:   if (i<0 || i>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
1073:   if (eig) *eig = sr->eigr[sr->perm[i]];
1074:   if (v) { BVCopyVec(sr->V,sr->perm[i],v); }
1075:   return(0);
1076: }

1078: /*@
1079:    EPSKrylovSchurGetSubcommPairs - Gets the i-th eigenpair stored
1080:    internally in the subcommunicator to which the calling process belongs.

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

1084:    Input Parameter:
1085: +  eps - the eigenproblem solver context
1086: -  i   - index of the solution

1088:    Output Parameters:
1089: +  eig - the eigenvalue
1090: -  v   - the eigenvector

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

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

1100:    Level: advanced

1102: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo(), EPSKrylovSchurGetSubcommMats()
1103: @*/
1104: PetscErrorCode EPSKrylovSchurGetSubcommPairs(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1105: {

1111:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommPairs_C",(EPS,PetscInt,PetscScalar*,Vec),(eps,i,eig,v));
1112:   return(0);
1113: }

1115: static PetscErrorCode EPSKrylovSchurGetSubcommMats_KrylovSchur(EPS eps,Mat *A,Mat *B)
1116: {
1117:   PetscErrorCode  ierr;
1118:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

1121:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1122:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1123:   EPSGetOperators(ctx->eps,A,B);
1124:   return(0);
1125: }

1127: /*@C
1128:    EPSKrylovSchurGetSubcommMats - Gets the eigenproblem matrices stored
1129:    internally in the subcommunicator to which the calling process belongs.

1131:    Collective on the subcommunicator

1133:    Input Parameter:
1134: .  eps - the eigenproblem solver context

1136:    Output Parameters:
1137: +  A  - the matrix associated with the eigensystem
1138: -  B  - the second matrix in the case of generalized eigenproblems

1140:    Notes:
1141:    This is the analog of EPSGetOperators(), but returns the matrices distributed
1142:    differently (in the subcommunicator rather than in the parent communicator).

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

1146:    Level: advanced

1148: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo()
1149: @*/
1150: PetscErrorCode EPSKrylovSchurGetSubcommMats(EPS eps,Mat *A,Mat *B)
1151: {

1156:   PetscTryMethod(eps,"EPSKrylovSchurGetSubcommMats_C",(EPS,Mat*,Mat*),(eps,A,B));
1157:   return(0);
1158: }

1160: static PetscErrorCode EPSKrylovSchurUpdateSubcommMats_KrylovSchur(EPS eps,PetscScalar a,PetscScalar ap,Mat Au,PetscScalar b,PetscScalar bp, Mat Bu,MatStructure str,PetscBool globalup)
1161: {
1162:   PetscErrorCode  ierr;
1163:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data,*subctx;
1164:   Mat             A,B=NULL,Ag,Bg=NULL;
1165:   PetscBool       reuse=PETSC_TRUE;

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

1173:   MatScale(A,a);
1174:   if (Au) {
1175:     MatAXPY(A,ap,Au,str);
1176:   }
1177:   if (B) MatScale(B,b);
1178:   if (Bu) {
1179:     MatAXPY(B,bp,Bu,str);
1180:   }
1181:   EPSSetOperators(ctx->eps,A,B);

1183:   /* Update stored matrix state */
1184:   subctx = (EPS_KRYLOVSCHUR*)ctx->eps->data;
1185:   PetscObjectStateGet((PetscObject)A,&subctx->Astate);
1186:   if (B) { PetscObjectStateGet((PetscObject)B,&subctx->Bstate); }

1188:   /* Update matrices in the parent communicator if requested by user */
1189:   if (globalup) {
1190:     if (ctx->npart>1) {
1191:       if (!ctx->isrow) {
1192:         MatGetOwnershipIS(Ag,&ctx->isrow,&ctx->iscol);
1193:         reuse = PETSC_FALSE;
1194:       }
1195:       if (str==DIFFERENT_NONZERO_PATTERN) reuse = PETSC_FALSE;
1196:       if (ctx->submata && !reuse) {
1197:         MatDestroyMatrices(1,&ctx->submata);
1198:       }
1199:       MatCreateSubMatrices(A,1,&ctx->isrow,&ctx->iscol,(reuse)?MAT_REUSE_MATRIX:MAT_INITIAL_MATRIX,&ctx->submata);
1200:       MatCreateMPIMatConcatenateSeqMat(((PetscObject)Ag)->comm,ctx->submata[0],PETSC_DECIDE,MAT_REUSE_MATRIX,&Ag);
1201:       if (B) {
1202:         if (ctx->submatb && !reuse) {
1203:           MatDestroyMatrices(1,&ctx->submatb);
1204:         }
1205:         MatCreateSubMatrices(B,1,&ctx->isrow,&ctx->iscol,(reuse)?MAT_REUSE_MATRIX:MAT_INITIAL_MATRIX,&ctx->submatb);
1206:         MatCreateMPIMatConcatenateSeqMat(((PetscObject)Bg)->comm,ctx->submatb[0],PETSC_DECIDE,MAT_REUSE_MATRIX,&Bg);
1207:       }
1208:     }
1209:     PetscObjectStateGet((PetscObject)Ag,&ctx->Astate);
1210:     if (Bg) { PetscObjectStateGet((PetscObject)Bg,&ctx->Bstate); }
1211:   }
1212:   EPSSetOperators(eps,Ag,Bg);
1213:   return(0);
1214: }

1216: /*@
1217:    EPSKrylovSchurUpdateSubcommMats - Update the eigenproblem matrices stored
1218:    internally in the subcommunicator to which the calling process belongs.

1220:    Collective on EPS

1222:    Input Parameters:
1223: +  eps - the eigenproblem solver context
1224: .  s   - scalar that multiplies the existing A matrix
1225: .  a   - scalar used in the axpy operation on A
1226: .  Au  - matrix used in the axpy operation on A
1227: .  t   - scalar that multiplies the existing B matrix
1228: .  b   - scalar used in the axpy operation on B
1229: .  Bu  - matrix used in the axpy operation on B
1230: .  str - structure flag
1231: -  globalup - flag indicating if global matrices must be updated

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

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

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

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

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

1250:    Level: advanced

1252: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommMats()
1253: @*/
1254: PetscErrorCode EPSKrylovSchurUpdateSubcommMats(EPS eps,PetscScalar s,PetscScalar a,Mat Au,PetscScalar t,PetscScalar b,Mat Bu,MatStructure str,PetscBool globalup)
1255: {

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

1272: PetscErrorCode EPSSetFromOptions_KrylovSchur(PetscOptionItems *PetscOptionsObject,EPS eps)
1273: {
1274:   PetscErrorCode  ierr;
1275:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1276:   PetscBool       flg,lock,b,f1,f2,f3;
1277:   PetscReal       keep;
1278:   PetscInt        i,j,k;

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

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

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

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

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

1297:     i = 1;
1298:     j = k = PETSC_DECIDE;
1299:     PetscOptionsInt("-eps_krylovschur_nev","Number of eigenvalues to compute in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",40,&i,&f1);
1300:     PetscOptionsInt("-eps_krylovschur_ncv","Number of basis vectors in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&j,&f2);
1301:     PetscOptionsInt("-eps_krylovschur_mpd","Maximum dimension of projected problem in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&k,&f3);
1302:     if (f1 || f2 || f3) { EPSKrylovSchurSetDimensions(eps,i,j,k); }

1304:   PetscOptionsTail();
1305:   return(0);
1306: }

1308: PetscErrorCode EPSView_KrylovSchur(EPS eps,PetscViewer viewer)
1309: {
1310:   PetscErrorCode  ierr;
1311:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1312:   PetscBool       isascii,isfilt;

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

1335: PetscErrorCode EPSDestroy_KrylovSchur(EPS eps)
1336: {

1340:   PetscFree(eps->data);
1341:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",NULL);
1342:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",NULL);
1343:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetLocking_C",NULL);
1344:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetLocking_C",NULL);
1345:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetPartitions_C",NULL);
1346:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetPartitions_C",NULL);
1347:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDetectZeros_C",NULL);
1348:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDetectZeros_C",NULL);
1349:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",NULL);
1350:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",NULL);
1351:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetSubintervals_C",NULL);
1352:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubintervals_C",NULL);
1353:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",NULL);
1354:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommInfo_C",NULL);
1355:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommPairs_C",NULL);
1356:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommMats_C",NULL);
1357:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurUpdateSubcommMats_C",NULL);
1358:   return(0);
1359: }

1361: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
1362: {
1364:   PetscBool      isfilt;

1367:   PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilt);
1368:   if (eps->which==EPS_ALL && !isfilt) {
1369:     EPSReset_KrylovSchur_Slice(eps);
1370:   }
1371:   return(0);
1372: }

1374: PetscErrorCode EPSSetDefaultST_KrylovSchur(EPS eps)
1375: {

1379:   if (eps->which==EPS_ALL) {
1380:     if (!((PetscObject)eps->st)->type_name) {
1381:       STSetType(eps->st,STSINVERT);
1382:     }
1383:   }
1384:   return(0);
1385: }

1387: SLEPC_EXTERN PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
1388: {
1389:   EPS_KRYLOVSCHUR *ctx;
1390:   PetscErrorCode  ierr;

1393:   PetscNewLog(eps,&ctx);
1394:   eps->data   = (void*)ctx;
1395:   ctx->lock   = PETSC_TRUE;
1396:   ctx->nev    = 1;
1397:   ctx->npart  = 1;
1398:   ctx->detect = PETSC_FALSE;
1399:   ctx->global = PETSC_TRUE;

1401:   eps->useds = PETSC_TRUE;
1402:   eps->hasts = PETSC_TRUE;

1404:   /* solve and computevectors determined at setup */
1405:   eps->ops->setup          = EPSSetUp_KrylovSchur;
1406:   eps->ops->setfromoptions = EPSSetFromOptions_KrylovSchur;
1407:   eps->ops->destroy        = EPSDestroy_KrylovSchur;
1408:   eps->ops->reset          = EPSReset_KrylovSchur;
1409:   eps->ops->view           = EPSView_KrylovSchur;
1410:   eps->ops->backtransform  = EPSBackTransform_Default;
1411:   eps->ops->setdefaultst   = EPSSetDefaultST_KrylovSchur;

1413:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",EPSKrylovSchurSetRestart_KrylovSchur);
1414:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",EPSKrylovSchurGetRestart_KrylovSchur);
1415:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetLocking_C",EPSKrylovSchurSetLocking_KrylovSchur);
1416:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetLocking_C",EPSKrylovSchurGetLocking_KrylovSchur);
1417:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetPartitions_C",EPSKrylovSchurSetPartitions_KrylovSchur);
1418:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetPartitions_C",EPSKrylovSchurGetPartitions_KrylovSchur);
1419:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDetectZeros_C",EPSKrylovSchurSetDetectZeros_KrylovSchur);
1420:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDetectZeros_C",EPSKrylovSchurGetDetectZeros_KrylovSchur);
1421:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",EPSKrylovSchurSetDimensions_KrylovSchur);
1422:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",EPSKrylovSchurGetDimensions_KrylovSchur);
1423:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetSubintervals_C",EPSKrylovSchurSetSubintervals_KrylovSchur);
1424:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubintervals_C",EPSKrylovSchurGetSubintervals_KrylovSchur);
1425:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",EPSKrylovSchurGetInertias_KrylovSchur);
1426:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommInfo_C",EPSKrylovSchurGetSubcommInfo_KrylovSchur);
1427:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommPairs_C",EPSKrylovSchurGetSubcommPairs_KrylovSchur);
1428:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommMats_C",EPSKrylovSchurGetSubcommMats_KrylovSchur);
1429:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurUpdateSubcommMats_C",EPSKrylovSchurUpdateSubcommMats_KrylovSchur);
1430:   return(0);
1431: }