Actual source code: ks-slice.c

slepc-3.23.1 2025-05-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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 with spectrum slicing for symmetric eigenproblems

 15:    References:

 17:        [1] R.G. Grimes et al., "A shifted block Lanczos algorithm for
 18:            solving sparse symmetric generalized eigenproblems", SIAM J.
 19:            Matrix Anal. Appl. 15(1):228-272, 1994.

 21:        [2] C. Campos and J.E. Roman, "Spectrum slicing strategies based
 22:            on restarted Lanczos methods", Numer. Algor. 60(2):279-295,
 23:            2012.
 24: */

 26: #include <slepc/private/epsimpl.h>
 27: #include "krylovschur.h"

 29: static PetscBool  cited = PETSC_FALSE;
 30: static const char citation[] =
 31:   "@Article{slepc-slice,\n"
 32:   "   author = \"C. Campos and J. E. Roman\",\n"
 33:   "   title = \"Strategies for spectrum slicing based on restarted {Lanczos} methods\",\n"
 34:   "   journal = \"Numer. Algorithms\",\n"
 35:   "   volume = \"60\",\n"
 36:   "   number = \"2\",\n"
 37:   "   pages = \"279--295\",\n"
 38:   "   year = \"2012,\"\n"
 39:   "   doi = \"https://doi.org/10.1007/s11075-012-9564-z\"\n"
 40:   "}\n";

 42: #define SLICE_PTOL PETSC_SQRT_MACHINE_EPSILON

 44: static PetscErrorCode EPSSliceResetSR(EPS eps)
 45: {
 46:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
 47:   EPS_SR          sr=ctx->sr;
 48:   EPS_shift       s;

 50:   PetscFunctionBegin;
 51:   if (sr) {
 52:     if (ctx->npart>1) {
 53:       PetscCall(BVDestroy(&sr->V));
 54:       PetscCall(PetscFree4(sr->eigr,sr->eigi,sr->errest,sr->perm));
 55:     }
 56:     /* Reviewing list of shifts to free memory */
 57:     s = sr->s0;
 58:     if (s) {
 59:       while (s->neighb[1]) {
 60:         s = s->neighb[1];
 61:         PetscCall(PetscFree(s->neighb[0]));
 62:       }
 63:       PetscCall(PetscFree(s));
 64:     }
 65:     PetscCall(PetscFree(sr));
 66:   }
 67:   ctx->sr = NULL;
 68:   PetscFunctionReturn(PETSC_SUCCESS);
 69: }

 71: PetscErrorCode EPSReset_KrylovSchur_Slice(EPS eps)
 72: {
 73:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

 75:   PetscFunctionBegin;
 76:   if (!ctx->global) PetscFunctionReturn(PETSC_SUCCESS);
 77:   /* Reset auxiliary EPS */
 78:   PetscCall(EPSSliceResetSR(ctx->eps));
 79:   PetscCall(EPSReset(ctx->eps));
 80:   PetscCall(EPSSliceResetSR(eps));
 81:   PetscCall(PetscFree(ctx->inertias));
 82:   PetscCall(PetscFree(ctx->shifts));
 83:   PetscFunctionReturn(PETSC_SUCCESS);
 84: }

 86: PetscErrorCode EPSDestroy_KrylovSchur_Slice(EPS eps)
 87: {
 88:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

 90:   PetscFunctionBegin;
 91:   if (!ctx->global) PetscFunctionReturn(PETSC_SUCCESS);
 92:   /* Destroy auxiliary EPS */
 93:   PetscCall(EPSReset_KrylovSchur_Slice(eps));
 94:   PetscCall(EPSDestroy(&ctx->eps));
 95:   if (ctx->npart>1) {
 96:     PetscCall(PetscSubcommDestroy(&ctx->subc));
 97:     if (ctx->commset) {
 98:       PetscCallMPI(MPI_Comm_free(&ctx->commrank));
 99:       ctx->commset = PETSC_FALSE;
100:     }
101:     PetscCall(ISDestroy(&ctx->isrow));
102:     PetscCall(ISDestroy(&ctx->iscol));
103:     PetscCall(MatDestroyMatrices(1,&ctx->submata));
104:     PetscCall(MatDestroyMatrices(1,&ctx->submatb));
105:   }
106:   PetscCall(PetscFree(ctx->subintervals));
107:   PetscCall(PetscFree(ctx->nconv_loc));
108:   PetscFunctionReturn(PETSC_SUCCESS);
109: }

111: /*
112:   EPSSliceAllocateSolution - Allocate memory storage for common variables such
113:   as eigenvalues and eigenvectors. The argument extra is used for methods
114:   that require a working basis slightly larger than ncv.
115: */
116: static PetscErrorCode EPSSliceAllocateSolution(EPS eps,PetscInt extra)
117: {
118:   EPS_KRYLOVSCHUR    *ctx=(EPS_KRYLOVSCHUR*)eps->data;
119:   PetscReal          eta;
120:   PetscInt           k;
121:   BVType             type;
122:   BVOrthogType       orthog_type;
123:   BVOrthogRefineType orthog_ref;
124:   BVOrthogBlockType  ob_type;
125:   Mat                matrix;
126:   Vec                t;
127:   EPS_SR             sr = ctx->sr;

129:   PetscFunctionBegin;
130:   /* allocate space for eigenvalues and friends */
131:   k = PetscMax(1,sr->numEigs);
132:   PetscCall(PetscFree4(sr->eigr,sr->eigi,sr->errest,sr->perm));
133:   PetscCall(PetscMalloc4(k,&sr->eigr,k,&sr->eigi,k,&sr->errest,k,&sr->perm));

135:   /* allocate sr->V and transfer options from eps->V */
136:   PetscCall(BVDestroy(&sr->V));
137:   PetscCall(BVCreate(PetscObjectComm((PetscObject)eps),&sr->V));
138:   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
139:   if (!((PetscObject)eps->V)->type_name) PetscCall(BVSetType(sr->V,BVMAT));
140:   else {
141:     PetscCall(BVGetType(eps->V,&type));
142:     PetscCall(BVSetType(sr->V,type));
143:   }
144:   PetscCall(STMatCreateVecsEmpty(eps->st,&t,NULL));
145:   PetscCall(BVSetSizesFromVec(sr->V,t,k));
146:   PetscCall(VecDestroy(&t));
147:   PetscCall(EPS_SetInnerProduct(eps));
148:   PetscCall(BVGetMatrix(eps->V,&matrix,NULL));
149:   PetscCall(BVSetMatrix(sr->V,matrix,PETSC_FALSE));
150:   PetscCall(BVGetOrthogonalization(eps->V,&orthog_type,&orthog_ref,&eta,&ob_type));
151:   PetscCall(BVSetOrthogonalization(sr->V,orthog_type,orthog_ref,eta,ob_type));
152:   PetscFunctionReturn(PETSC_SUCCESS);
153: }

155: static PetscErrorCode EPSSliceGetEPS(EPS eps)
156: {
157:   EPS_KRYLOVSCHUR    *ctx=(EPS_KRYLOVSCHUR*)eps->data,*ctx_local;
158:   BV                 V;
159:   BVType             type;
160:   PetscReal          eta;
161:   BVOrthogType       orthog_type;
162:   BVOrthogRefineType orthog_ref;
163:   BVOrthogBlockType  ob_type;
164:   PetscInt           i;
165:   PetscReal          h,a,b;
166:   PetscRandom        rand;
167:   EPS_SR             sr=ctx->sr;

169:   PetscFunctionBegin;
170:   if (!ctx->eps) PetscCall(EPSKrylovSchurGetChildEPS(eps,&ctx->eps));

172:   /* Determine subintervals */
173:   if (ctx->npart==1) {
174:     a = eps->inta; b = eps->intb;
175:   } else {
176:     if (!ctx->subintset) { /* uniform distribution if no set by user */
177:       PetscCheck(sr->hasEnd,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Global interval must be bounded for splitting it in uniform subintervals");
178:       h = (eps->intb-eps->inta)/ctx->npart;
179:       a = eps->inta+ctx->subc->color*h;
180:       b = (ctx->subc->color==ctx->npart-1)?eps->intb:eps->inta+(ctx->subc->color+1)*h;
181:       PetscCall(PetscFree(ctx->subintervals));
182:       PetscCall(PetscMalloc1(ctx->npart+1,&ctx->subintervals));
183:       for (i=0;i<ctx->npart;i++) ctx->subintervals[i] = eps->inta+h*i;
184:       ctx->subintervals[ctx->npart] = eps->intb;
185:     } else {
186:       a = ctx->subintervals[ctx->subc->color];
187:       b = ctx->subintervals[ctx->subc->color+1];
188:     }
189:   }
190:   PetscCall(EPSSetInterval(ctx->eps,a,b));
191:   PetscCall(EPSSetConvergenceTest(ctx->eps,eps->conv));
192:   PetscCall(EPSSetDimensions(ctx->eps,ctx->nev,ctx->ncv,ctx->mpd));
193:   PetscCall(EPSKrylovSchurSetLocking(ctx->eps,ctx->lock));

195:   ctx_local = (EPS_KRYLOVSCHUR*)ctx->eps->data;
196:   ctx_local->detect = ctx->detect;

198:   /* transfer options from eps->V */
199:   PetscCall(EPSGetBV(ctx->eps,&V));
200:   PetscCall(BVGetRandomContext(V,&rand));  /* make sure the random context is available when duplicating */
201:   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
202:   if (!((PetscObject)eps->V)->type_name) PetscCall(BVSetType(V,BVMAT));
203:   else {
204:     PetscCall(BVGetType(eps->V,&type));
205:     PetscCall(BVSetType(V,type));
206:   }
207:   PetscCall(BVGetOrthogonalization(eps->V,&orthog_type,&orthog_ref,&eta,&ob_type));
208:   PetscCall(BVSetOrthogonalization(V,orthog_type,orthog_ref,eta,ob_type));

210:   ctx->eps->which = eps->which;
211:   ctx->eps->max_it = eps->max_it;
212:   ctx->eps->tol = eps->tol;
213:   ctx->eps->purify = eps->purify;
214:   if (eps->tol==(PetscReal)PETSC_DETERMINE) eps->tol = SLEPC_DEFAULT_TOL;
215:   PetscCall(EPSSetProblemType(ctx->eps,eps->problem_type));
216:   PetscCall(EPSSetUp(ctx->eps));
217:   ctx->eps->nconv = 0;
218:   ctx->eps->its   = 0;
219:   for (i=0;i<ctx->eps->ncv;i++) {
220:     ctx->eps->eigr[i]   = 0.0;
221:     ctx->eps->eigi[i]   = 0.0;
222:     ctx->eps->errest[i] = 0.0;
223:   }
224:   PetscFunctionReturn(PETSC_SUCCESS);
225: }

227: static PetscErrorCode EPSSliceGetInertia(EPS eps,PetscReal shift,PetscInt *inertia,PetscInt *zeros)
228: {
229:   KSP            ksp,kspr;
230:   PC             pc;
231:   Mat            F;
232:   PetscReal      nzshift=shift;
233:   PetscBool      flg;

235:   PetscFunctionBegin;
236:   if (shift >= PETSC_MAX_REAL) { /* Right-open interval */
237:     if (inertia) *inertia = eps->n;
238:   } else if (shift <= PETSC_MIN_REAL) {
239:     if (inertia) *inertia = 0;
240:     if (zeros) *zeros = 0;
241:   } else {
242:     /* If the shift is zero, perturb it to a very small positive value.
243:        The goal is that the nonzero pattern is the same in all cases and reuse
244:        the symbolic factorizations */
245:     nzshift = (shift==0.0)? 10.0/PETSC_MAX_REAL: shift;
246:     PetscCall(STSetShift(eps->st,nzshift));
247:     PetscCall(STGetKSP(eps->st,&ksp));
248:     PetscCall(KSPGetPC(ksp,&pc));
249:     PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCREDUNDANT,&flg));
250:     if (flg) {
251:       PetscCall(PCRedundantGetKSP(pc,&kspr));
252:       PetscCall(KSPGetPC(kspr,&pc));
253:     }
254:     PetscCall(PCFactorGetMatrix(pc,&F));
255:     PetscCall(MatGetInertia(F,inertia,zeros,NULL));
256:   }
257:   if (inertia) PetscCall(PetscInfo(eps,"Computed inertia at shift %g: %" PetscInt_FMT "\n",(double)nzshift,*inertia));
258:   PetscFunctionReturn(PETSC_SUCCESS);
259: }

261: /*
262:    Dummy backtransform operation
263:  */
264: static PetscErrorCode EPSBackTransform_Skip(EPS eps)
265: {
266:   PetscFunctionBegin;
267:   PetscFunctionReturn(PETSC_SUCCESS);
268: }

270: PetscErrorCode EPSSetUp_KrylovSchur_Slice(EPS eps)
271: {
272:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data,*ctx_glob;
273:   EPS_SR          sr,sr_loc,sr_glob;
274:   PetscInt        nEigs,dssz=1,i,zeros=0,off=0,method,hiteig=0;
275:   PetscMPIInt     nproc,rank=0,aux;
276:   PetscReal       r;
277:   MPI_Request     req;
278:   Mat             A,B=NULL;
279:   DSParallelType  ptype;
280:   MPI_Comm        child;

282:   PetscFunctionBegin;
283:   if (eps->nev==0) eps->nev = 1;
284:   if (ctx->global) {
285:     EPSCheckHermitianDefiniteCondition(eps,PETSC_TRUE," with spectrum slicing");
286:     EPSCheckSinvertCayleyCondition(eps,PETSC_TRUE," with spectrum slicing");
287:     PetscCheck(eps->inta!=eps->intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues unless you provide a computational interval with EPSSetInterval()");
288:     PetscCheck(eps->intb<PETSC_MAX_REAL || eps->inta>PETSC_MIN_REAL,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The defined computational interval should have at least one of their sides bounded");
289:     PetscCheck(eps->nds==0,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing not supported in combination with deflation space");
290:     EPSCheckUnsupportedCondition(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING,PETSC_TRUE," with spectrum slicing");
291:     EPSCheckIgnoredCondition(eps,EPS_FEATURE_BALANCE,PETSC_TRUE," with spectrum slicing");
292:     if (eps->tol==(PetscReal)PETSC_DETERMINE) {
293: #if defined(PETSC_USE_REAL_SINGLE)
294:       eps->tol = SLEPC_DEFAULT_TOL;
295: #else
296:       /* use tighter tolerance */
297:       eps->tol = SLEPC_DEFAULT_TOL*1e-2;
298: #endif
299:     }
300:     if (eps->max_it==PETSC_DETERMINE) eps->max_it = 100;
301:     if (ctx->nev==1) ctx->nev = PetscMin(40,eps->n);  /* nev not set, use default value */
302:     PetscCheck(eps->n<=10 || ctx->nev>=10,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"nev cannot be less than 10 in spectrum slicing runs");
303:   }
304:   eps->ops->backtransform = EPSBackTransform_Skip;

306:   /* create spectrum slicing context and initialize it */
307:   PetscCall(EPSSliceResetSR(eps));
308:   PetscCall(PetscNew(&sr));
309:   ctx->sr = sr;
310:   sr->itsKs = 0;
311:   sr->nleap = 0;
312:   sr->nMAXCompl = eps->nev/4;
313:   sr->iterCompl = eps->max_it/4;
314:   sr->sPres = NULL;
315:   sr->nS = 0;

317:   if (ctx->npart==1 || ctx->global) {
318:     /* check presence of ends and finding direction */
319:     if ((eps->inta > PETSC_MIN_REAL && !(ctx->subintervals && ctx->subintervals[0]==ctx->subintervals[1])) || eps->intb >= PETSC_MAX_REAL) {
320:       sr->int0 = eps->inta;
321:       sr->int1 = eps->intb;
322:       sr->dir = 1;
323:       if (eps->intb >= PETSC_MAX_REAL) { /* Right-open interval */
324:         sr->hasEnd = PETSC_FALSE;
325:       } else sr->hasEnd = PETSC_TRUE;
326:     } else {
327:       sr->int0 = eps->intb;
328:       sr->int1 = eps->inta;
329:       sr->dir = -1;
330:       sr->hasEnd = PetscNot(eps->inta <= PETSC_MIN_REAL);
331:     }
332:   }
333:   if (ctx->global) {
334:     PetscCall(EPSSetDimensions_Default(eps,&ctx->nev,&ctx->ncv,&ctx->mpd));
335:     /* create subintervals and initialize auxiliary eps for slicing runs */
336:     PetscCall(EPSKrylovSchurGetChildEPS(eps,&ctx->eps));
337:     /* prevent computation of factorization in global eps */
338:     PetscCall(STSetTransform(eps->st,PETSC_FALSE));
339:     PetscCall(EPSSliceGetEPS(eps));
340:     sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
341:     if (ctx->npart>1) {
342:       PetscCall(PetscSubcommGetChild(ctx->subc,&child));
343:       if ((sr->dir>0&&ctx->subc->color==0)||(sr->dir<0&&ctx->subc->color==ctx->npart-1)) sr->inertia0 = sr_loc->inertia0;
344:       PetscCallMPI(MPI_Comm_rank(child,&rank));
345:       if (!rank) {
346:         PetscCall(PetscMPIIntCast((sr->dir>0)?0:ctx->npart-1,&aux));
347:         PetscCallMPI(MPI_Bcast(&sr->inertia0,1,MPIU_INT,aux,ctx->commrank));
348:       }
349:       PetscCallMPI(MPI_Bcast(&sr->inertia0,1,MPIU_INT,0,child));
350:       PetscCall(PetscFree(ctx->nconv_loc));
351:       PetscCall(PetscMalloc1(ctx->npart,&ctx->nconv_loc));
352:       PetscCallMPI(MPI_Comm_size(((PetscObject)eps)->comm,&nproc));
353:       if (sr->dir<0) off = 1;
354:       if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
355:         PetscCall(PetscMPIIntCast(sr_loc->numEigs,&aux));
356:         PetscCallMPI(MPI_Allgather(&aux,1,MPI_INT,ctx->nconv_loc,1,MPI_INT,ctx->commrank));
357:         PetscCallMPI(MPI_Allgather(sr_loc->dir==sr->dir?&sr_loc->int0:&sr_loc->int1,1,MPIU_REAL,ctx->subintervals+off,1,MPIU_REAL,ctx->commrank));
358:       } else {
359:         PetscCallMPI(MPI_Comm_rank(child,&rank));
360:         if (!rank) {
361:           PetscCall(PetscMPIIntCast(sr_loc->numEigs,&aux));
362:           PetscCallMPI(MPI_Allgather(&aux,1,MPI_INT,ctx->nconv_loc,1,MPI_INT,ctx->commrank));
363:           PetscCallMPI(MPI_Allgather(sr_loc->dir==sr->dir?&sr_loc->int0:&sr_loc->int1,1,MPIU_REAL,ctx->subintervals+off,1,MPIU_REAL,ctx->commrank));
364:         }
365:         PetscCall(PetscMPIIntCast(ctx->npart,&aux));
366:         PetscCallMPI(MPI_Bcast(ctx->nconv_loc,aux,MPI_INT,0,child));
367:         PetscCallMPI(MPI_Bcast(ctx->subintervals+off,aux,MPIU_REAL,0,child));
368:       }
369:       nEigs = 0;
370:       for (i=0;i<ctx->npart;i++) nEigs += ctx->nconv_loc[i];
371:     } else {
372:       nEigs = sr_loc->numEigs;
373:       sr->inertia0 = sr_loc->inertia0;
374:       sr->dir = sr_loc->dir;
375:     }
376:     sr->inertia1 = sr->inertia0+sr->dir*nEigs;
377:     sr->numEigs = nEigs;
378:     eps->nev = nEigs;
379:     eps->ncv = nEigs;
380:     eps->mpd = nEigs;
381:   } else {
382:     ctx_glob = (EPS_KRYLOVSCHUR*)ctx->eps->data;
383:     sr_glob = ctx_glob->sr;
384:     if (ctx->npart>1) {
385:       sr->dir = sr_glob->dir;
386:       sr->int0 = (sr->dir==1)?eps->inta:eps->intb;
387:       sr->int1 = (sr->dir==1)?eps->intb:eps->inta;
388:       if ((sr->dir>0&&ctx->subc->color==ctx->npart-1)||(sr->dir<0&&ctx->subc->color==0)) sr->hasEnd = sr_glob->hasEnd;
389:       else sr->hasEnd = PETSC_TRUE;
390:     }
391:     /* sets first shift */
392:     PetscCall(STSetShift(eps->st,(sr->int0==0.0)?10.0/PETSC_MAX_REAL:sr->int0));
393:     PetscCall(STSetUp(eps->st));

395:     /* compute inertia0 */
396:     PetscCall(EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,ctx->detect?&zeros:NULL));
397:     /* undocumented option to control what to do when an eigenvalue is found:
398:        - error out if it's the endpoint of the user-provided interval (or sub-interval)
399:        - if it's an endpoint computed internally:
400:           + if hiteig=0 error out
401:           + else if hiteig=1 the subgroup that hit the eigenvalue does nothing
402:           + otherwise the subgroup that hit the eigenvalue perturbs the shift and recomputes inertia
403:     */
404:     PetscCall(PetscOptionsGetInt(NULL,NULL,"-eps_krylovschur_hiteigenvalue",&hiteig,NULL));
405:     if (zeros) { /* error in factorization */
406:       PetscCheck(sr->int0!=ctx->eps->inta && sr->int0!=ctx->eps->intb,((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in the interval endpoint");
407:       PetscCheck(!ctx_glob->subintset || hiteig,((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
408:       if (hiteig==1) { /* idle subgroup */
409:         sr->inertia0 = -1;
410:       } else { /* perturb shift */
411:         sr->int0 *= (1.0+SLICE_PTOL);
412:         PetscCall(EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,&zeros));
413:         PetscCheck(zeros==0,((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",(double)sr->int1);
414:       }
415:     }
416:     if (ctx->npart>1) {
417:       PetscCall(PetscSubcommGetChild(ctx->subc,&child));
418:       /* inertia1 is received from neighbour */
419:       PetscCallMPI(MPI_Comm_rank(child,&rank));
420:       if (!rank) {
421:         if (sr->inertia0!=-1 && ((sr->dir>0 && ctx->subc->color>0) || (sr->dir<0 && ctx->subc->color<ctx->npart-1))) { /* send inertia0 to neighbour0 */
422:           PetscCall(PetscMPIIntCast(ctx->subc->color-sr->dir,&aux));
423:           PetscCallMPI(MPI_Isend(&sr->inertia0,1,MPIU_INT,aux,0,ctx->commrank,&req));
424:           PetscCallMPI(MPI_Isend(&sr->int0,1,MPIU_REAL,aux,0,ctx->commrank,&req));
425:         }
426:         if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)|| (sr->dir<0 && ctx->subc->color>0)) { /* receive inertia1 from neighbour1 */
427:           PetscCall(PetscMPIIntCast(ctx->subc->color+sr->dir,&aux));
428:           PetscCallMPI(MPI_Recv(&sr->inertia1,1,MPIU_INT,aux,0,ctx->commrank,MPI_STATUS_IGNORE));
429:           PetscCallMPI(MPI_Recv(&sr->int1,1,MPIU_REAL,aux,0,ctx->commrank,MPI_STATUS_IGNORE));
430:         }
431:         if (sr->inertia0==-1 && !(sr->dir>0 && ctx->subc->color==ctx->npart-1) && !(sr->dir<0 && ctx->subc->color==0)) {
432:           sr->inertia0 = sr->inertia1; sr->int0 = sr->int1;
433:           PetscCall(PetscMPIIntCast(ctx->subc->color-sr->dir,&aux));
434:           PetscCallMPI(MPI_Isend(&sr->inertia0,1,MPIU_INT,aux,0,ctx->commrank,&req));
435:           PetscCallMPI(MPI_Isend(&sr->int0,1,MPIU_REAL,aux,0,ctx->commrank,&req));
436:         }
437:       }
438:       if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)||(sr->dir<0 && ctx->subc->color>0)) {
439:         PetscCallMPI(MPI_Bcast(&sr->inertia1,1,MPIU_INT,0,child));
440:         PetscCallMPI(MPI_Bcast(&sr->int1,1,MPIU_REAL,0,child));
441:       } else sr_glob->inertia1 = sr->inertia1;
442:     }

444:     /* last process in eps comm computes inertia1 */
445:     if (ctx->npart==1 || ((sr->dir>0 && ctx->subc->color==ctx->npart-1) || (sr->dir<0 && ctx->subc->color==0))) {
446:       PetscCall(EPSSliceGetInertia(eps,sr->int1,&sr->inertia1,ctx->detect?&zeros:NULL));
447:       PetscCheck(zeros==0,((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
448:       if (!rank && sr->inertia0==-1) {
449:         sr->inertia0 = sr->inertia1; sr->int0 = sr->int1;
450:         PetscCall(PetscMPIIntCast(ctx->subc->color-sr->dir,&aux));
451:         PetscCallMPI(MPI_Isend(&sr->inertia0,1,MPIU_INT,aux,0,ctx->commrank,&req));
452:         PetscCallMPI(MPI_Isend(&sr->int0,1,MPIU_REAL,aux,0,ctx->commrank,&req));
453:       }
454:       if (sr->hasEnd) {
455:         sr->dir = -sr->dir; r = sr->int0; sr->int0 = sr->int1; sr->int1 = r;
456:         i = sr->inertia0; sr->inertia0 = sr->inertia1; sr->inertia1 = i;
457:       }
458:     }

460:     /* number of eigenvalues in interval */
461:     sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
462:     if (ctx->npart>1) {
463:       /* memory allocate for subinterval eigenpairs */
464:       PetscCall(EPSSliceAllocateSolution(eps,1));
465:     }
466:     dssz = eps->ncv+1;
467:     PetscCall(DSGetParallel(ctx->eps->ds,&ptype));
468:     PetscCall(DSSetParallel(eps->ds,ptype));
469:     PetscCall(DSGetMethod(ctx->eps->ds,&method));
470:     PetscCall(DSSetMethod(eps->ds,method));
471:   }
472:   PetscCall(DSSetType(eps->ds,DSHEP));
473:   PetscCall(DSSetCompact(eps->ds,PETSC_TRUE));
474:   PetscCall(DSAllocate(eps->ds,dssz));
475:   /* keep state of subcomm matrices to check that the user does not modify them */
476:   PetscCall(EPSGetOperators(eps,&A,&B));
477:   PetscCall(MatGetState(A,&ctx->Astate));
478:   PetscCall(PetscObjectGetId((PetscObject)A,&ctx->Aid));
479:   if (B) {
480:     PetscCall(MatGetState(B,&ctx->Bstate));
481:     PetscCall(PetscObjectGetId((PetscObject)B,&ctx->Bid));
482:   } else {
483:     ctx->Bstate=0;
484:     ctx->Bid=0;
485:   }
486:   PetscFunctionReturn(PETSC_SUCCESS);
487: }

489: static PetscErrorCode EPSSliceGatherEigenVectors(EPS eps)
490: {
491:   Vec             v,vg,v_loc;
492:   IS              is1,is2;
493:   VecScatter      vec_sc;
494:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
495:   PetscInt        nloc,m0,n0,i,si,idx,*idx1,*idx2,j;
496:   PetscScalar     *array;
497:   EPS_SR          sr_loc;
498:   BV              V_loc;

500:   PetscFunctionBegin;
501:   sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
502:   V_loc = sr_loc->V;

504:   /* Gather parallel eigenvectors */
505:   PetscCall(BVGetColumn(eps->V,0,&v));
506:   PetscCall(VecGetOwnershipRange(v,&n0,&m0));
507:   PetscCall(BVRestoreColumn(eps->V,0,&v));
508:   PetscCall(BVGetColumn(ctx->eps->V,0,&v));
509:   PetscCall(VecGetLocalSize(v,&nloc));
510:   PetscCall(BVRestoreColumn(ctx->eps->V,0,&v));
511:   PetscCall(PetscMalloc2(m0-n0,&idx1,m0-n0,&idx2));
512:   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)eps),nloc,PETSC_DECIDE,&vg));
513:   idx = -1;
514:   for (si=0;si<ctx->npart;si++) {
515:     j = 0;
516:     for (i=n0;i<m0;i++) {
517:       idx1[j]   = i;
518:       idx2[j++] = i+eps->n*si;
519:     }
520:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx1,PETSC_COPY_VALUES,&is1));
521:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx2,PETSC_COPY_VALUES,&is2));
522:     PetscCall(BVGetColumn(eps->V,0,&v));
523:     PetscCall(VecScatterCreate(v,is1,vg,is2,&vec_sc));
524:     PetscCall(BVRestoreColumn(eps->V,0,&v));
525:     PetscCall(ISDestroy(&is1));
526:     PetscCall(ISDestroy(&is2));
527:     for (i=0;i<ctx->nconv_loc[si];i++) {
528:       PetscCall(BVGetColumn(eps->V,++idx,&v));
529:       if (ctx->subc->color==si) {
530:         PetscCall(BVGetColumn(V_loc,i,&v_loc));
531:         PetscCall(VecGetArray(v_loc,&array));
532:         PetscCall(VecPlaceArray(vg,array));
533:       }
534:       PetscCall(VecScatterBegin(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE));
535:       PetscCall(VecScatterEnd(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE));
536:       if (ctx->subc->color==si) {
537:         PetscCall(VecResetArray(vg));
538:         PetscCall(VecRestoreArray(v_loc,&array));
539:         PetscCall(BVRestoreColumn(V_loc,i,&v_loc));
540:       }
541:       PetscCall(BVRestoreColumn(eps->V,idx,&v));
542:     }
543:     PetscCall(VecScatterDestroy(&vec_sc));
544:   }
545:   PetscCall(PetscFree2(idx1,idx2));
546:   PetscCall(VecDestroy(&vg));
547:   PetscFunctionReturn(PETSC_SUCCESS);
548: }

550: /*
551:   EPSComputeVectors_Slice - Recover Eigenvectors from subcomunicators
552:  */
553: PetscErrorCode EPSComputeVectors_Slice(EPS eps)
554: {
555:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

557:   PetscFunctionBegin;
558:   if (ctx->global && ctx->npart>1) {
559:     PetscCall(EPSComputeVectors(ctx->eps));
560:     PetscCall(EPSSliceGatherEigenVectors(eps));
561:   }
562:   PetscFunctionReturn(PETSC_SUCCESS);
563: }

565: static PetscErrorCode EPSSliceGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
566: {
567:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
568:   PetscInt        i=0,j,tmpi;
569:   PetscReal       v,tmpr;
570:   EPS_shift       s;

572:   PetscFunctionBegin;
573:   PetscCheck(eps->state,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
574:   PetscCheck(ctx->sr,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
575:   if (!ctx->sr->s0) {  /* EPSSolve not called yet */
576:     *n = 2;
577:   } else {
578:     *n = 1;
579:     s = ctx->sr->s0;
580:     while (s) {
581:       (*n)++;
582:       s = s->neighb[1];
583:     }
584:   }
585:   PetscCall(PetscMalloc1(*n,shifts));
586:   PetscCall(PetscMalloc1(*n,inertias));
587:   if (!ctx->sr->s0) {  /* EPSSolve not called yet */
588:     (*shifts)[0]   = ctx->sr->int0;
589:     (*shifts)[1]   = ctx->sr->int1;
590:     (*inertias)[0] = ctx->sr->inertia0;
591:     (*inertias)[1] = ctx->sr->inertia1;
592:   } else {
593:     s = ctx->sr->s0;
594:     while (s) {
595:       (*shifts)[i]     = s->value;
596:       (*inertias)[i++] = s->inertia;
597:       s = s->neighb[1];
598:     }
599:     (*shifts)[i]   = ctx->sr->int1;
600:     (*inertias)[i] = ctx->sr->inertia1;
601:   }
602:   /* remove possible duplicate in last position */
603:   if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
604:   /* sort result */
605:   for (i=0;i<*n;i++) {
606:     v = (*shifts)[i];
607:     for (j=i+1;j<*n;j++) {
608:       if (v > (*shifts)[j]) {
609:         SlepcSwap((*shifts)[i],(*shifts)[j],tmpr);
610:         SlepcSwap((*inertias)[i],(*inertias)[j],tmpi);
611:         v = (*shifts)[i];
612:       }
613:     }
614:   }
615:   PetscFunctionReturn(PETSC_SUCCESS);
616: }

618: static PetscErrorCode EPSSliceGatherSolution(EPS eps)
619: {
620:   PetscMPIInt     rank,nproc,*disp,*ns_loc,aux;
621:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
622:   PetscInt        i,idx,j,*perm_loc,off=0,*inertias_loc,ns;
623:   PetscScalar     *eigr_loc;
624:   EPS_SR          sr_loc;
625:   PetscReal       *shifts_loc;
626:   MPI_Comm        child;

628:   PetscFunctionBegin;
629:   eps->nconv = 0;
630:   for (i=0;i<ctx->npart;i++) eps->nconv += ctx->nconv_loc[i];
631:   sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

633:   /* Gather the shifts used and the inertias computed */
634:   PetscCall(EPSSliceGetInertias(ctx->eps,&ns,&shifts_loc,&inertias_loc));
635:   if (ctx->sr->dir>0 && shifts_loc[ns-1]==sr_loc->int1 && ctx->subc->color<ctx->npart-1) ns--;
636:   if (ctx->sr->dir<0 && shifts_loc[ns-1]==sr_loc->int0 && ctx->subc->color>0) {
637:     ns--;
638:     for (i=0;i<ns;i++) {
639:       inertias_loc[i] = inertias_loc[i+1];
640:       shifts_loc[i] = shifts_loc[i+1];
641:     }
642:   }
643:   PetscCall(PetscMalloc1(ctx->npart,&ns_loc));
644:   PetscCall(PetscSubcommGetChild(ctx->subc,&child));
645:   PetscCallMPI(MPI_Comm_rank(child,&rank));
646:   PetscCall(PetscMPIIntCast(ns,&aux));
647:   if (!rank) PetscCallMPI(MPI_Allgather(&aux,1,MPI_INT,ns_loc,1,MPI_INT,ctx->commrank));
648:   PetscCall(PetscMPIIntCast(ctx->npart,&aux));
649:   PetscCallMPI(MPI_Bcast(ns_loc,aux,MPI_INT,0,child));
650:   ctx->nshifts = 0;
651:   for (i=0;i<ctx->npart;i++) ctx->nshifts += ns_loc[i];
652:   PetscCall(PetscFree(ctx->inertias));
653:   PetscCall(PetscFree(ctx->shifts));
654:   PetscCall(PetscMalloc1(ctx->nshifts,&ctx->inertias));
655:   PetscCall(PetscMalloc1(ctx->nshifts,&ctx->shifts));

657:   /* Gather eigenvalues (same ranks have fully set of eigenvalues)*/
658:   eigr_loc = sr_loc->eigr;
659:   perm_loc = sr_loc->perm;
660:   PetscCallMPI(MPI_Comm_size(((PetscObject)eps)->comm,&nproc));
661:   PetscCall(PetscMalloc1(ctx->npart,&disp));
662:   disp[0] = 0;
663:   for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ctx->nconv_loc[i-1];
664:   if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
665:     PetscCall(PetscMPIIntCast(sr_loc->numEigs,&aux));
666:     PetscCallMPI(MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank)); /* eigenvalues */
667:     PetscCallMPI(MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank)); /* perm */
668:     for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
669:     PetscCall(PetscMPIIntCast(ns,&aux));
670:     PetscCallMPI(MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank)); /* shifts */
671:     PetscCallMPI(MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank)); /* inertias */
672:     PetscCallMPI(MPIU_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank));
673:   } else { /* subcommunicators with different size */
674:     if (!rank) {
675:       PetscCall(PetscMPIIntCast(sr_loc->numEigs,&aux));
676:       PetscCallMPI(MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank)); /* eigenvalues */
677:       PetscCallMPI(MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank)); /* perm */
678:       for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
679:       PetscCall(PetscMPIIntCast(ns,&aux));
680:       PetscCallMPI(MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank)); /* shifts */
681:       PetscCallMPI(MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank)); /* inertias */
682:       PetscCallMPI(MPIU_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank));
683:     }
684:     PetscCall(PetscMPIIntCast(eps->nconv,&aux));
685:     PetscCallMPI(MPI_Bcast(eps->eigr,aux,MPIU_SCALAR,0,child));
686:     PetscCallMPI(MPI_Bcast(eps->perm,aux,MPIU_INT,0,child));
687:     PetscCall(PetscMPIIntCast(ctx->nshifts,&aux));
688:     PetscCallMPI(MPI_Bcast(ctx->shifts,aux,MPIU_REAL,0,child));
689:     PetscCallMPI(MPI_Bcast(ctx->inertias,aux,MPIU_INT,0,child));
690:     PetscCallMPI(MPI_Bcast(&eps->its,1,MPIU_INT,0,child));
691:   }
692:   /* Update global array eps->perm */
693:   idx = ctx->nconv_loc[0];
694:   for (i=1;i<ctx->npart;i++) {
695:     off += ctx->nconv_loc[i-1];
696:     for (j=0;j<ctx->nconv_loc[i];j++) eps->perm[idx++] += off;
697:   }

699:   /* Gather parallel eigenvectors */
700:   PetscCall(PetscFree(ns_loc));
701:   PetscCall(PetscFree(disp));
702:   PetscCall(PetscFree(shifts_loc));
703:   PetscCall(PetscFree(inertias_loc));
704:   PetscFunctionReturn(PETSC_SUCCESS);
705: }

707: /*
708:    Fills the fields of a shift structure
709: */
710: static PetscErrorCode EPSCreateShift(EPS eps,PetscReal val,EPS_shift neighb0,EPS_shift neighb1)
711: {
712:   EPS_shift       s,*pending2;
713:   PetscInt        i;
714:   EPS_SR          sr;
715:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

717:   PetscFunctionBegin;
718:   sr = ctx->sr;
719:   if ((neighb0 && val==neighb0->value) || (neighb1 && val==neighb1->value)) {
720:     sr->nPend++;
721:     PetscFunctionReturn(PETSC_SUCCESS);
722:   }
723:   PetscCall(PetscNew(&s));
724:   s->value = val;
725:   s->neighb[0] = neighb0;
726:   if (neighb0) neighb0->neighb[1] = s;
727:   s->neighb[1] = neighb1;
728:   if (neighb1) neighb1->neighb[0] = s;
729:   s->comp[0] = PETSC_FALSE;
730:   s->comp[1] = PETSC_FALSE;
731:   s->index = -1;
732:   s->neigs = 0;
733:   s->nconv[0] = s->nconv[1] = 0;
734:   s->nsch[0] = s->nsch[1]=0;
735:   /* Inserts in the stack of pending shifts */
736:   /* If needed, the array is resized */
737:   if (sr->nPend >= sr->maxPend) {
738:     sr->maxPend *= 2;
739:     PetscCall(PetscMalloc1(sr->maxPend,&pending2));
740:     for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
741:     PetscCall(PetscFree(sr->pending));
742:     sr->pending = pending2;
743:   }
744:   sr->pending[sr->nPend++]=s;
745:   PetscFunctionReturn(PETSC_SUCCESS);
746: }

748: /* Prepare for Rational Krylov update */
749: static PetscErrorCode EPSPrepareRational(EPS eps)
750: {
751:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
752:   PetscInt        dir,i,k,ld,nv;
753:   PetscScalar     *A;
754:   EPS_SR          sr = ctx->sr;
755:   Vec             v;

757:   PetscFunctionBegin;
758:   PetscCall(DSGetLeadingDimension(eps->ds,&ld));
759:   dir = (sr->sPres->neighb[0] == sr->sPrev)?1:-1;
760:   dir*=sr->dir;
761:   k = 0;
762:   for (i=0;i<sr->nS;i++) {
763:     if (dir*PetscRealPart(sr->S[i])>0.0) {
764:       sr->S[k] = sr->S[i];
765:       sr->S[sr->nS+k] = sr->S[sr->nS+i];
766:       PetscCall(BVGetColumn(sr->Vnext,k,&v));
767:       PetscCall(BVCopyVec(eps->V,eps->nconv+i,v));
768:       PetscCall(BVRestoreColumn(sr->Vnext,k,&v));
769:       k++;
770:       if (k>=sr->nS/2) break;
771:     }
772:   }
773:   /* Copy to DS */
774:   PetscCall(DSSetCompact(eps->ds,PETSC_FALSE));  /* make sure DS_MAT_A is allocated */
775:   PetscCall(DSGetArray(eps->ds,DS_MAT_A,&A));
776:   PetscCall(PetscArrayzero(A,ld*ld));
777:   for (i=0;i<k;i++) {
778:     A[i*(1+ld)] = sr->S[i];
779:     A[k+i*ld] = sr->S[sr->nS+i];
780:   }
781:   sr->nS = k;
782:   PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&A));
783:   PetscCall(DSGetDimensions(eps->ds,&nv,NULL,NULL,NULL));
784:   PetscCall(DSSetDimensions(eps->ds,nv,0,k));
785:   /* Append u to V */
786:   PetscCall(BVGetColumn(sr->Vnext,sr->nS,&v));
787:   PetscCall(BVCopyVec(eps->V,sr->nv,v));
788:   PetscCall(BVRestoreColumn(sr->Vnext,sr->nS,&v));
789:   PetscFunctionReturn(PETSC_SUCCESS);
790: }

792: /* Provides next shift to be computed */
793: static PetscErrorCode EPSExtractShift(EPS eps)
794: {
795:   PetscInt        iner,zeros=0;
796:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
797:   EPS_SR          sr;
798:   PetscReal       newShift,diam,ptol;
799:   EPS_shift       sPres;

801:   PetscFunctionBegin;
802:   sr = ctx->sr;
803:   if (sr->nPend > 0) {
804:     if (sr->sPres==sr->pending[sr->nPend-1]) {
805:       eps->reason = EPS_CONVERGED_ITERATING;
806:       eps->its = 0;
807:       sr->nPend--;
808:       sr->sPres->rep = PETSC_TRUE;
809:       PetscFunctionReturn(PETSC_SUCCESS);
810:     }
811:     sr->sPrev = sr->sPres;
812:     sr->sPres = sr->pending[--sr->nPend];
813:     sPres = sr->sPres;
814:     PetscCall(EPSSliceGetInertia(eps,sPres->value,&iner,ctx->detect?&zeros:NULL));
815:     if (zeros) {
816:       diam = PetscMin(PetscAbsReal(sPres->neighb[0]->value-sPres->value),PetscAbsReal(sPres->neighb[1]->value-sPres->value));
817:       ptol = PetscMin(SLICE_PTOL,diam/2);
818:       newShift = sPres->value*(1.0+ptol);
819:       if (sr->dir*(sPres->neighb[0] && newShift-sPres->neighb[0]->value) < 0) newShift = (sPres->value+sPres->neighb[0]->value)/2;
820:       else if (sPres->neighb[1] && sr->dir*(sPres->neighb[1]->value-newShift) < 0) newShift = (sPres->value+sPres->neighb[1]->value)/2;
821:       PetscCall(EPSSliceGetInertia(eps,newShift,&iner,&zeros));
822:       PetscCheck(zeros==0,((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",(double)newShift);
823:       sPres->value = newShift;
824:     }
825:     sr->sPres->inertia = iner;
826:     eps->target = sr->sPres->value;
827:     eps->reason = EPS_CONVERGED_ITERATING;
828:     eps->its = 0;
829:   } else sr->sPres = NULL;
830:   PetscFunctionReturn(PETSC_SUCCESS);
831: }

833: /*
834:    Symmetric KrylovSchur adapted to spectrum slicing:
835:    Allows searching an specific amount of eigenvalues in the subintervals left and right.
836:    Returns whether the search has succeeded
837: */
838: static PetscErrorCode EPSKrylovSchur_Slice(EPS eps)
839: {
840:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
841:   PetscInt        i,k,l,ld,nv,*iwork,j,count0,count1,iterCompl=0,n0,n1;
842:   Mat             U,Op,T;
843:   PetscScalar     *Q,*A;
844:   PetscReal       *a,*b,beta,lambda;
845:   EPS_shift       sPres;
846:   PetscBool       breakdown,complIterating,sch0,sch1;
847:   EPS_SR          sr = ctx->sr;

849:   PetscFunctionBegin;
850:   /* Spectrum slicing data */
851:   sPres = sr->sPres;
852:   complIterating =PETSC_FALSE;
853:   sch1 = sch0 = PETSC_TRUE;
854:   PetscCall(DSGetLeadingDimension(eps->ds,&ld));
855:   PetscCall(PetscMalloc1(2*ld,&iwork));
856:   count0=0;count1=0; /* Found on both sides */
857:   if (!sPres->rep && sr->nS > 0 && (sPres->neighb[0] == sr->sPrev || sPres->neighb[1] == sr->sPrev)) {
858:     /* Rational Krylov */
859:     PetscCall(DSTranslateRKS(eps->ds,sr->sPrev->value-sPres->value));
860:     PetscCall(DSGetDimensions(eps->ds,NULL,NULL,&l,NULL));
861:     PetscCall(DSSetDimensions(eps->ds,l+1,0,0));
862:     PetscCall(BVSetActiveColumns(eps->V,0,l+1));
863:     PetscCall(DSGetMat(eps->ds,DS_MAT_Q,&U));
864:     PetscCall(BVMultInPlace(eps->V,U,0,l+1));
865:     PetscCall(DSRestoreMat(eps->ds,DS_MAT_Q,&U));
866:   } else {
867:     /* Get the starting Lanczos vector */
868:     PetscCall(EPSGetStartVector(eps,0,NULL));
869:     l = 0;
870:   }
871:   /* Restart loop */
872:   while (eps->reason == EPS_CONVERGED_ITERATING) {
873:     eps->its++; sr->itsKs++;
874:     /* Compute an nv-step Lanczos factorization */
875:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
876:     PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l));
877:     PetscCall(DSGetMat(eps->ds,DS_MAT_T,&T));
878:     PetscCall(STGetOperator(eps->st,&Op));
879:     PetscCall(BVMatLanczos(eps->V,Op,T,eps->nconv+l,&nv,&beta,&breakdown));
880:     PetscCall(STRestoreOperator(eps->st,&Op));
881:     sr->nv = nv;
882:     PetscCall(DSRestoreMat(eps->ds,DS_MAT_T,&T));
883:     PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l));
884:     if (l==0) PetscCall(DSSetState(eps->ds,DS_STATE_INTERMEDIATE));
885:     else PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
886:     PetscCall(BVSetActiveColumns(eps->V,eps->nconv,nv));

888:     /* Solve projected problem and compute residual norm estimates */
889:     if (eps->its == 1 && l > 0) { /* After rational update, DS_MAT_A is available */
890:       PetscCall(DSGetArray(eps->ds,DS_MAT_A,&A));
891:       PetscCall(DSGetArrayReal(eps->ds,DS_MAT_T,&a));
892:       b = a + ld;
893:       k = eps->nconv+l;
894:       A[k*ld+k-1] = A[(k-1)*ld+k];
895:       A[k*ld+k] = a[k];
896:       for (j=k+1; j< nv; j++) {
897:         A[j*ld+j] = a[j];
898:         A[j*ld+j-1] = b[j-1] ;
899:         A[(j-1)*ld+j] = b[j-1];
900:       }
901:       PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&A));
902:       PetscCall(DSRestoreArrayReal(eps->ds,DS_MAT_T,&a));
903:       PetscCall(DSSolve(eps->ds,eps->eigr,NULL));
904:       PetscCall(DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL));
905:       PetscCall(DSSetCompact(eps->ds,PETSC_TRUE));
906:     } else { /* Restart */
907:       PetscCall(DSSolve(eps->ds,eps->eigr,NULL));
908:       PetscCall(DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL));
909:     }
910:     PetscCall(DSSynchronize(eps->ds,eps->eigr,NULL));

912:     /* Residual */
913:     PetscCall(EPSKrylovConvergence(eps,PETSC_TRUE,eps->nconv,nv-eps->nconv,beta,0.0,1.0,&k));
914:     /* Checking values obtained for completing */
915:     for (i=0;i<k;i++) {
916:       sr->back[i]=eps->eigr[i];
917:     }
918:     PetscCall(STBackTransform(eps->st,k,sr->back,eps->eigi));
919:     count0=count1=0;
920:     for (i=0;i<k;i++) {
921:       lambda = PetscRealPart(sr->back[i]);
922:       if ((sr->dir*(sPres->value - lambda) > 0) && (sr->dir*(lambda - sPres->ext[0]) > 0)) count0++;
923:       if ((sr->dir*(lambda - sPres->value) > 0) && (sr->dir*(sPres->ext[1] - lambda) > 0)) count1++;
924:     }
925:     if (k>eps->nev && eps->ncv-k<5) eps->reason = EPS_CONVERGED_TOL;
926:     else {
927:       /* Checks completion */
928:       if ((!sch0||count0 >= sPres->nsch[0]) && (!sch1 ||count1 >= sPres->nsch[1])) {
929:         eps->reason = EPS_CONVERGED_TOL;
930:       } else {
931:         if (!complIterating && eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
932:         if (complIterating) {
933:           if (--iterCompl <= 0) eps->reason = EPS_DIVERGED_ITS;
934:         } else if (k >= eps->nev) {
935:           n0 = sPres->nsch[0]-count0;
936:           n1 = sPres->nsch[1]-count1;
937:           if (sr->iterCompl>0 && ((n0>0 && n0<= sr->nMAXCompl)||(n1>0&&n1<=sr->nMAXCompl))) {
938:             /* Iterating for completion*/
939:             complIterating = PETSC_TRUE;
940:             if (n0 >sr->nMAXCompl)sch0 = PETSC_FALSE;
941:             if (n1 >sr->nMAXCompl)sch1 = PETSC_FALSE;
942:             iterCompl = sr->iterCompl;
943:           } else eps->reason = EPS_CONVERGED_TOL;
944:         }
945:       }
946:     }
947:     /* Update l */
948:     if (eps->reason == EPS_CONVERGED_ITERATING) l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
949:     else l = nv-k;
950:     if (breakdown) l=0;
951:     if (!ctx->lock && l>0 && eps->reason == EPS_CONVERGED_ITERATING) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */

953:     if (eps->reason == EPS_CONVERGED_ITERATING) {
954:       if (breakdown) {
955:         /* Start a new Lanczos factorization */
956:         PetscCall(PetscInfo(eps,"Breakdown in Krylov-Schur method (it=%" PetscInt_FMT " norm=%g)\n",eps->its,(double)beta));
957:         PetscCall(EPSGetStartVector(eps,k,&breakdown));
958:         if (breakdown) {
959:           eps->reason = EPS_DIVERGED_BREAKDOWN;
960:           PetscCall(PetscInfo(eps,"Unable to generate more start vectors\n"));
961:         }
962:       } else {
963:         /* Prepare the Rayleigh quotient for restart */
964:         PetscCall(DSGetArrayReal(eps->ds,DS_MAT_T,&a));
965:         PetscCall(DSGetArray(eps->ds,DS_MAT_Q,&Q));
966:         b = a + ld;
967:         for (i=k;i<k+l;i++) {
968:           a[i] = PetscRealPart(eps->eigr[i]);
969:           b[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
970:         }
971:         PetscCall(DSRestoreArrayReal(eps->ds,DS_MAT_T,&a));
972:         PetscCall(DSRestoreArray(eps->ds,DS_MAT_Q,&Q));
973:       }
974:     }
975:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
976:     PetscCall(DSGetMat(eps->ds,DS_MAT_Q,&U));
977:     PetscCall(BVMultInPlace(eps->V,U,eps->nconv,k+l));
978:     PetscCall(DSRestoreMat(eps->ds,DS_MAT_Q,&U));

980:     /* Normalize u and append it to V */
981:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) PetscCall(BVCopyColumn(eps->V,nv,k+l));
982:     eps->nconv = k;
983:     if (eps->reason != EPS_CONVERGED_ITERATING) {
984:       /* Store approximated values for next shift */
985:       PetscCall(DSGetArray(eps->ds,DS_MAT_Q,&Q));
986:       sr->nS = l;
987:       for (i=0;i<l;i++) {
988:         sr->S[i] = eps->eigr[i+k];/* Diagonal elements */
989:         sr->S[i+l] = Q[nv-1+(i+k)*ld]*beta; /* Out of diagonal elements */
990:       }
991:       PetscCall(DSRestoreArray(eps->ds,DS_MAT_Q,&Q));
992:     }
993:   }
994:   /* Check for completion */
995:   for (i=0;i< eps->nconv; i++) {
996:     if (sr->dir*PetscRealPart(eps->eigr[i])>0) sPres->nconv[1]++;
997:     else sPres->nconv[0]++;
998:   }
999:   sPres->comp[0] = PetscNot(count0 < sPres->nsch[0]);
1000:   sPres->comp[1] = PetscNot(count1 < sPres->nsch[1]);
1001:   PetscCall(PetscInfo(eps,"Lanczos: %" PetscInt_FMT " evals in [%g,%g]%s and %" PetscInt_FMT " evals in [%g,%g]%s\n",count0,(double)(sr->dir==1?sPres->ext[0]:sPres->value),(double)(sr->dir==1?sPres->value:sPres->ext[0]),sPres->comp[0]?"*":"",count1,(double)(sr->dir==1?sPres->value:sPres->ext[1]),(double)(sr->dir==1?sPres->ext[1]:sPres->value),sPres->comp[1]?"*":""));
1002:   PetscCheck(count0<=sPres->nsch[0] && count1<=sPres->nsch[1],PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Mismatch between number of values found and information from inertia%s",ctx->detect?"":", consider using EPSKrylovSchurSetDetectZeros()");
1003:   PetscCall(PetscFree(iwork));
1004:   PetscFunctionReturn(PETSC_SUCCESS);
1005: }

1007: /*
1008:   Obtains value of subsequent shift
1009: */
1010: static PetscErrorCode EPSGetNewShiftValue(EPS eps,PetscInt side,PetscReal *newS)
1011: {
1012:   PetscReal       lambda,d_prev;
1013:   PetscInt        i,idxP;
1014:   EPS_SR          sr;
1015:   EPS_shift       sPres,s;
1016:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

1018:   PetscFunctionBegin;
1019:   sr = ctx->sr;
1020:   sPres = sr->sPres;
1021:   if (sPres->neighb[side]) {
1022:     /* Completing a previous interval */
1023:     *newS = (sPres->value + sPres->neighb[side]->value)/2;
1024:     if (PetscAbsReal(sPres->value - *newS)/PetscAbsReal(sPres->value)<=100*PETSC_SQRT_MACHINE_EPSILON) *newS = sPres->value;
1025:   } else { /* (Only for side=1). Creating a new interval. */
1026:     if (sPres->neigs==0) {/* No value has been accepted*/
1027:       if (sPres->neighb[0]) {
1028:         /* Multiplying by 10 the previous distance */
1029:         *newS = sPres->value + 10*sr->dir*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
1030:         sr->nleap++;
1031:         /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
1032:         PetscCheck(sr->hasEnd || sr->nleap<=5,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Unable to compute the wanted eigenvalues with open interval");
1033:       } else { /* First shift */
1034:         PetscCheck(eps->nconv!=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"First shift renders no information");
1035:         /* Unaccepted values give information for next shift */
1036:         idxP=0;/* Number of values left from shift */
1037:         for (i=0;i<eps->nconv;i++) {
1038:           lambda = PetscRealPart(eps->eigr[i]);
1039:           if (sr->dir*(lambda - sPres->value) <0) idxP++;
1040:           else break;
1041:         }
1042:         /* Avoiding subtraction of eigenvalues (might be the same).*/
1043:         if (idxP>0) {
1044:           d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[0]))/(idxP+0.3);
1045:         } else {
1046:           d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[eps->nconv-1]))/(eps->nconv+0.3);
1047:         }
1048:         *newS = sPres->value + (sr->dir*d_prev*eps->nev)/2;
1049:       }
1050:     } else { /* Accepted values found */
1051:       sr->nleap = 0;
1052:       /* Average distance of values in previous subinterval */
1053:       s = sPres->neighb[0];
1054:       while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
1055:         s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
1056:       }
1057:       if (s) {
1058:         d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
1059:       } else { /* First shift. Average distance obtained with values in this shift */
1060:         /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
1061:         if (sr->dir*(PetscRealPart(sr->eigr[0])-sPres->value)>0 && PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0]))/PetscRealPart(sr->eigr[0])) > PetscSqrtReal(eps->tol)) {
1062:           d_prev =  PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0])))/(sPres->neigs+0.3);
1063:         } else {
1064:           d_prev = PetscAbsReal(PetscRealPart(sr->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
1065:         }
1066:       }
1067:       /* Average distance is used for next shift by adding it to value on the right or to shift */
1068:       if (sr->dir*(PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
1069:         *newS = PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1])+ (sr->dir*d_prev*eps->nev)/2;
1070:       } else { /* Last accepted value is on the left of shift. Adding to shift */
1071:         *newS = sPres->value + (sr->dir*d_prev*eps->nev)/2;
1072:       }
1073:     }
1074:     /* End of interval can not be surpassed */
1075:     if (sr->dir*(sr->int1 - *newS) < 0) *newS = sr->int1;
1076:   }/* of neighb[side]==null */
1077:   PetscFunctionReturn(PETSC_SUCCESS);
1078: }

1080: /*
1081:   Function for sorting an array of real values
1082: */
1083: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
1084: {
1085:   PetscReal re;
1086:   PetscInt  i,j,tmp;

1088:   PetscFunctionBegin;
1089:   if (!prev) for (i=0;i<nr;i++) perm[i] = i;
1090:   /* Insertion sort */
1091:   for (i=1;i<nr;i++) {
1092:     re = PetscRealPart(r[perm[i]]);
1093:     j = i-1;
1094:     while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
1095:       tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
1096:     }
1097:   }
1098:   PetscFunctionReturn(PETSC_SUCCESS);
1099: }

1101: /* Stores the pairs obtained since the last shift in the global arrays */
1102: static PetscErrorCode EPSStoreEigenpairs(EPS eps)
1103: {
1104:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1105:   PetscReal       lambda,err,norm;
1106:   PetscInt        i,count;
1107:   PetscBool       iscayley;
1108:   EPS_SR          sr = ctx->sr;
1109:   EPS_shift       sPres;
1110:   Vec             v,w;

1112:   PetscFunctionBegin;
1113:   sPres = sr->sPres;
1114:   sPres->index = sr->indexEig;
1115:   count = sr->indexEig;
1116:   /* Back-transform */
1117:   PetscCall(STBackTransform(eps->st,eps->nconv,eps->eigr,eps->eigi));
1118:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley));
1119:   /* Sort eigenvalues */
1120:   PetscCall(sortRealEigenvalues(eps->eigr,eps->perm,eps->nconv,PETSC_FALSE,sr->dir));
1121:   /* Values stored in global array */
1122:   for (i=0;i<eps->nconv;i++) {
1123:     lambda = PetscRealPart(eps->eigr[eps->perm[i]]);
1124:     err = eps->errest[eps->perm[i]];

1126:     if (sr->dir*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
1127:       PetscCheck(count<sr->numEigs,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Unexpected error in Spectrum Slicing");
1128:       sr->eigr[count] = lambda;
1129:       sr->errest[count] = err;
1130:       /* Explicit purification */
1131:       PetscCall(BVGetColumn(eps->V,eps->perm[i],&w));
1132:       if (eps->purify) {
1133:         PetscCall(BVGetColumn(sr->V,count,&v));
1134:         PetscCall(STApply(eps->st,w,v));
1135:         PetscCall(BVRestoreColumn(sr->V,count,&v));
1136:       } else PetscCall(BVInsertVec(sr->V,count,w));
1137:       PetscCall(BVRestoreColumn(eps->V,eps->perm[i],&w));
1138:       PetscCall(BVNormColumn(sr->V,count,NORM_2,&norm));
1139:       PetscCall(BVScaleColumn(sr->V,count,1.0/norm));
1140:       count++;
1141:     }
1142:   }
1143:   sPres->neigs = count - sr->indexEig;
1144:   sr->indexEig = count;
1145:   /* Global ordering array updating */
1146:   PetscCall(sortRealEigenvalues(sr->eigr,sr->perm,count,PETSC_TRUE,sr->dir));
1147:   PetscFunctionReturn(PETSC_SUCCESS);
1148: }

1150: static PetscErrorCode EPSLookForDeflation(EPS eps)
1151: {
1152:   PetscReal       val;
1153:   PetscInt        i,count0=0,count1=0;
1154:   EPS_shift       sPres;
1155:   PetscInt        ini,fin,k,idx0,idx1;
1156:   EPS_SR          sr;
1157:   Vec             v;
1158:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

1160:   PetscFunctionBegin;
1161:   sr = ctx->sr;
1162:   sPres = sr->sPres;

1164:   if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
1165:   else ini = 0;
1166:   fin = sr->indexEig;
1167:   /* Selection of ends for searching new values */
1168:   if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
1169:   else sPres->ext[0] = sPres->neighb[0]->value;
1170:   if (!sPres->neighb[1]) {
1171:     if (sr->hasEnd) sPres->ext[1] = sr->int1;
1172:     else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
1173:   } else sPres->ext[1] = sPres->neighb[1]->value;
1174:   /* Selection of values between right and left ends */
1175:   for (i=ini;i<fin;i++) {
1176:     val=PetscRealPart(sr->eigr[sr->perm[i]]);
1177:     /* Values to the right of left shift */
1178:     if (sr->dir*(val - sPres->ext[1]) < 0) {
1179:       if (sr->dir*(val - sPres->value) < 0) count0++;
1180:       else count1++;
1181:     } else break;
1182:   }
1183:   /* The number of values on each side are found */
1184:   if (sPres->neighb[0]) {
1185:     sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
1186:     PetscCheck(sPres->nsch[0]>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Mismatch between number of values found and information from inertia%s",ctx->detect?"":", consider using EPSKrylovSchurSetDetectZeros()");
1187:   } else sPres->nsch[0] = 0;

1189:   if (sPres->neighb[1]) {
1190:     sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
1191:     PetscCheck(sPres->nsch[1]>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Mismatch between number of values found and information from inertia%s",ctx->detect?"":", consider using EPSKrylovSchurSetDetectZeros()");
1192:   } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);

1194:   /* Completing vector of indexes for deflation */
1195:   idx0 = ini;
1196:   idx1 = ini+count0+count1;
1197:   k=0;
1198:   for (i=idx0;i<idx1;i++) sr->idxDef[k++]=sr->perm[i];
1199:   PetscCall(BVDuplicateResize(eps->V,k+eps->ncv+1,&sr->Vnext));
1200:   PetscCall(BVSetNumConstraints(sr->Vnext,k));
1201:   for (i=0;i<k;i++) {
1202:     PetscCall(BVGetColumn(sr->Vnext,-i-1,&v));
1203:     PetscCall(BVCopyVec(sr->V,sr->idxDef[i],v));
1204:     PetscCall(BVRestoreColumn(sr->Vnext,-i-1,&v));
1205:   }

1207:   /* For rational Krylov */
1208:   if (!sr->sPres->rep && sr->nS>0 && (sr->sPrev == sr->sPres->neighb[0] || sr->sPrev == sr->sPres->neighb[1])) PetscCall(EPSPrepareRational(eps));
1209:   eps->nconv = 0;
1210:   /* Get rid of temporary Vnext */
1211:   PetscCall(BVDestroy(&eps->V));
1212:   eps->V = sr->Vnext;
1213:   sr->Vnext = NULL;
1214:   PetscFunctionReturn(PETSC_SUCCESS);
1215: }

1217: PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS eps)
1218: {
1219:   PetscInt         i,lds,ti;
1220:   PetscReal        newS;
1221:   EPS_KRYLOVSCHUR  *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1222:   EPS_SR           sr=ctx->sr;
1223:   Mat              A,B=NULL;
1224:   PetscObjectState Astate,Bstate=0;
1225:   PetscObjectId    Aid,Bid=0;

1227:   PetscFunctionBegin;
1228:   PetscCall(PetscCitationsRegister(citation,&cited));
1229:   if (ctx->global) {
1230:     PetscCall(EPSSolve_KrylovSchur_Slice(ctx->eps));
1231:     ctx->eps->state = EPS_STATE_SOLVED;
1232:     eps->reason = EPS_CONVERGED_TOL;
1233:     if (ctx->npart>1) {
1234:       /* Gather solution from subsolvers */
1235:       PetscCall(EPSSliceGatherSolution(eps));
1236:     } else {
1237:       eps->nconv = sr->numEigs;
1238:       eps->its   = ctx->eps->its;
1239:       PetscCall(PetscFree(ctx->inertias));
1240:       PetscCall(PetscFree(ctx->shifts));
1241:       PetscCall(EPSSliceGetInertias(ctx->eps,&ctx->nshifts,&ctx->shifts,&ctx->inertias));
1242:     }
1243:   } else {
1244:     if (ctx->npart==1) {
1245:       sr->eigr   = ctx->eps->eigr;
1246:       sr->eigi   = ctx->eps->eigi;
1247:       sr->perm   = ctx->eps->perm;
1248:       sr->errest = ctx->eps->errest;
1249:       sr->V      = ctx->eps->V;
1250:     }
1251:     /* Check that the user did not modify subcomm matrices */
1252:     PetscCall(EPSGetOperators(eps,&A,&B));
1253:     PetscCall(MatGetState(A,&Astate));
1254:     PetscCall(PetscObjectGetId((PetscObject)A,&Aid));
1255:     if (B) {
1256:       PetscCall(MatGetState(B,&Bstate));
1257:       PetscCall(PetscObjectGetId((PetscObject)B,&Bid));
1258:     }
1259:     PetscCheck(Astate==ctx->Astate && (!B || Bstate==ctx->Bstate) && Aid==ctx->Aid && (!B || Bid==ctx->Bid),PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Subcomm matrices have been modified by user");
1260:     /* Only with eigenvalues present in the interval ...*/
1261:     if (sr->numEigs==0) {
1262:       eps->reason = EPS_CONVERGED_TOL;
1263:       PetscFunctionReturn(PETSC_SUCCESS);
1264:     }
1265:     /* Array of pending shifts */
1266:     sr->maxPend = 100; /* Initial size */
1267:     sr->nPend = 0;
1268:     PetscCall(PetscMalloc1(sr->maxPend,&sr->pending));
1269:     PetscCall(EPSCreateShift(eps,sr->int0,NULL,NULL));
1270:     /* extract first shift */
1271:     sr->sPrev = NULL;
1272:     sr->sPres = sr->pending[--sr->nPend];
1273:     sr->sPres->inertia = sr->inertia0;
1274:     eps->target = sr->sPres->value;
1275:     sr->s0 = sr->sPres;
1276:     sr->indexEig = 0;
1277:     /* Memory reservation for auxiliary variables */
1278:     lds = PetscMin(eps->mpd,eps->ncv);
1279:     PetscCall(PetscCalloc1(lds*lds,&sr->S));
1280:     PetscCall(PetscMalloc1(eps->ncv,&sr->back));
1281:     for (i=0;i<sr->numEigs;i++) {
1282:       sr->eigr[i]   = 0.0;
1283:       sr->eigi[i]   = 0.0;
1284:       sr->errest[i] = 0.0;
1285:       sr->perm[i]   = i;
1286:     }
1287:     /* Vectors for deflation */
1288:     PetscCall(PetscMalloc1(sr->numEigs,&sr->idxDef));
1289:     sr->indexEig = 0;
1290:     /* Main loop */
1291:     while (sr->sPres) {
1292:       /* Search for deflation */
1293:       PetscCall(EPSLookForDeflation(eps));
1294:       /* KrylovSchur */
1295:       PetscCall(EPSKrylovSchur_Slice(eps));

1297:       PetscCall(EPSStoreEigenpairs(eps));
1298:       /* Select new shift */
1299:       if (!sr->sPres->comp[1]) {
1300:         PetscCall(EPSGetNewShiftValue(eps,1,&newS));
1301:         PetscCall(EPSCreateShift(eps,newS,sr->sPres,sr->sPres->neighb[1]));
1302:       }
1303:       if (!sr->sPres->comp[0]) {
1304:         /* Completing earlier interval */
1305:         PetscCall(EPSGetNewShiftValue(eps,0,&newS));
1306:         PetscCall(EPSCreateShift(eps,newS,sr->sPres->neighb[0],sr->sPres));
1307:       }
1308:       /* Preparing for a new search of values */
1309:       PetscCall(EPSExtractShift(eps));
1310:     }

1312:     /* Updating eps values prior to exit */
1313:     PetscCall(PetscFree(sr->S));
1314:     PetscCall(PetscFree(sr->idxDef));
1315:     PetscCall(PetscFree(sr->pending));
1316:     PetscCall(PetscFree(sr->back));
1317:     PetscCall(BVDuplicateResize(eps->V,eps->ncv+1,&sr->Vnext));
1318:     PetscCall(BVSetNumConstraints(sr->Vnext,0));
1319:     PetscCall(BVDestroy(&eps->V));
1320:     eps->V      = sr->Vnext;
1321:     eps->nconv  = sr->indexEig;
1322:     eps->reason = EPS_CONVERGED_TOL;
1323:     eps->its    = sr->itsKs;
1324:     eps->nds    = 0;
1325:     if (sr->dir<0) {
1326:       for (i=0;i<eps->nconv/2;i++) {
1327:         ti = sr->perm[i]; sr->perm[i] = sr->perm[eps->nconv-1-i]; sr->perm[eps->nconv-1-i] = ti;
1328:       }
1329:     }
1330:   }
1331:   PetscFunctionReturn(PETSC_SUCCESS);
1332: }