Actual source code: ks-slice.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 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: #define InertiaMismatch(h,d) \
 45:   do { \
 46:     SETERRQ1(PetscObjectComm((PetscObject)h),1,"Mismatch between number of values found and information from inertia%s",d?"":", consider using EPSKrylovSchurSetDetectZeros()"); \
 47:   } while (0)

 49: static PetscErrorCode EPSSliceResetSR(EPS eps)
 50: {
 51:   PetscErrorCode  ierr;
 52:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
 53:   EPS_SR          sr=ctx->sr;
 54:   EPS_shift       s;

 57:   if (sr) {
 58:     if (ctx->npart>1) {
 59:       BVDestroy(&sr->V);
 60:       PetscFree4(sr->eigr,sr->eigi,sr->errest,sr->perm);
 61:     }
 62:     /* Reviewing list of shifts to free memory */
 63:     s = sr->s0;
 64:     if (s) {
 65:       while (s->neighb[1]) {
 66:         s = s->neighb[1];
 67:         PetscFree(s->neighb[0]);
 68:       }
 69:       PetscFree(s);
 70:     }
 71:     PetscFree(sr);
 72:   }
 73:   ctx->sr = NULL;
 74:   return(0);
 75: }

 77: PetscErrorCode EPSReset_KrylovSchur_Slice(EPS eps)
 78: {
 79:   PetscErrorCode  ierr;
 80:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

 83:   if (!ctx->global) return(0);
 84:   /* Destroy auxiliary EPS */
 85:   EPSSliceResetSR(ctx->eps);
 86:   EPSDestroy(&ctx->eps);
 87:   if (ctx->npart>1) {
 88:     PetscSubcommDestroy(&ctx->subc);
 89:     if (ctx->commset) {
 90:       MPI_Comm_free(&ctx->commrank);
 91:       ctx->commset = PETSC_FALSE;
 92:     }
 93:   }
 94:   PetscFree(ctx->subintervals);
 95:   PetscFree(ctx->nconv_loc);
 96:   EPSSliceResetSR(eps);
 97:   PetscFree(ctx->inertias);
 98:   PetscFree(ctx->shifts);
 99:   if (ctx->npart>1) {
100:     ISDestroy(&ctx->isrow);
101:     ISDestroy(&ctx->iscol);
102:     MatDestroyMatrices(1,&ctx->submata);
103:     MatDestroyMatrices(1,&ctx->submatb);
104:   }
105:   return(0);
106: }

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

129:   /* allocate space for eigenvalues and friends */
130:   k = PetscMax(1,sr->numEigs);
131:   PetscFree4(sr->eigr,sr->eigi,sr->errest,sr->perm);
132:   PetscMalloc4(k,&sr->eigr,k,&sr->eigi,k,&sr->errest,k,&sr->perm);
133:   cnt = 2*k*sizeof(PetscScalar) + 2*k*sizeof(PetscReal) + k*sizeof(PetscInt);
134:   PetscLogObjectMemory((PetscObject)eps,cnt);

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

158: static PetscErrorCode EPSSliceGetEPS(EPS eps)
159: {
160:   PetscErrorCode     ierr;
161:   EPS_KRYLOVSCHUR    *ctx=(EPS_KRYLOVSCHUR*)eps->data,*ctx_local;
162:   BV                 V;
163:   BVType             type;
164:   PetscReal          eta;
165:   BVOrthogType       orthog_type;
166:   BVOrthogRefineType orthog_ref;
167:   BVOrthogBlockType  ob_type;
168:   Mat                A,B=NULL,Ar=NULL,Br=NULL;
169:   PetscInt           i;
170:   PetscReal          h,a,b,zero;
171:   PetscMPIInt        rank;
172:   EPS_SR             sr=ctx->sr;
173:   PC                 pc;
174:   PCType             pctype;
175:   KSP                ksp;
176:   KSPType            ksptype;
177:   STType             sttype;
178:   PetscObjectState   Astate,Bstate=0;
179:   PetscObjectId      Aid,Bid=0;
180:   MatSolverType      stype;

183:   EPSGetOperators(eps,&A,&B);
184:   if (ctx->npart==1) {
185:     if (!ctx->eps) { EPSCreate(((PetscObject)eps)->comm,&ctx->eps); }
186:     EPSSetType(ctx->eps,((PetscObject)eps)->type_name);
187:     EPSSetOperators(ctx->eps,A,B);
188:     a = eps->inta; b = eps->intb;
189:   } else {
190:     PetscObjectStateGet((PetscObject)A,&Astate);
191:     PetscObjectGetId((PetscObject)A,&Aid);
192:     if (B) {
193:       PetscObjectStateGet((PetscObject)B,&Bstate);
194:       PetscObjectGetId((PetscObject)B,&Bid);
195:     }
196:     if (!ctx->subc) {
197:       /* Create context for subcommunicators */
198:       PetscSubcommCreate(PetscObjectComm((PetscObject)eps),&ctx->subc);
199:       PetscSubcommSetNumber(ctx->subc,ctx->npart);
200:       PetscSubcommSetType(ctx->subc,PETSC_SUBCOMM_CONTIGUOUS);
201:       PetscLogObjectMemory((PetscObject)eps,sizeof(PetscSubcomm));

203:       /* Duplicate matrices */
204:       MatCreateRedundantMatrix(A,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Ar);
205:       ctx->Astate = Astate;
206:       ctx->Aid = Aid;
207:       if (B) {
208:         MatCreateRedundantMatrix(B,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Br);
209:         ctx->Bstate = Bstate;
210:         ctx->Bid = Bid;
211:       }
212:     } else {
213:       if (ctx->Astate != Astate || (B && ctx->Bstate != Bstate) || ctx->Aid != Aid || (B && ctx->Bid != Bid)) {
214:         EPSGetOperators(ctx->eps,&Ar,&Br);
215:         MatCreateRedundantMatrix(A,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Ar);
216:         ctx->Astate = Astate;
217:         ctx->Aid = Aid;
218:         if (B) {
219:           MatCreateRedundantMatrix(B,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Br);
220:           ctx->Bstate = Bstate;
221:           ctx->Bid = Bid;
222:         }
223:         EPSSetOperators(ctx->eps,Ar,Br);
224:         MatDestroy(&Ar);
225:         MatDestroy(&Br);
226:       }
227:     }

229:     /* Determine subintervals */
230:     if (!ctx->subintset) { /* uniform distribution if no set by user */
231:       if (!sr->hasEnd) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Global interval must be bounded for splitting it in uniform subintervals");
232:       h = (eps->intb-eps->inta)/ctx->npart;
233:       a = eps->inta+ctx->subc->color*h;
234:       b = (ctx->subc->color==ctx->npart-1)?eps->intb:eps->inta+(ctx->subc->color+1)*h;
235:       PetscFree(ctx->subintervals);
236:       PetscMalloc1(ctx->npart+1,&ctx->subintervals);
237:       for (i=0;i<ctx->npart;i++) ctx->subintervals[i] = eps->inta+h*i;
238:       ctx->subintervals[ctx->npart] = eps->intb;
239:     } else {
240:       a = ctx->subintervals[ctx->subc->color];
241:       b = ctx->subintervals[ctx->subc->color+1];
242:     }

244:     if (!ctx->eps) {
245:       /* Create auxiliary EPS */
246:       EPSCreate(PetscSubcommChild(ctx->subc),&ctx->eps);
247:       EPSSetOperators(ctx->eps,Ar,Br);
248:       MatDestroy(&Ar);
249:       MatDestroy(&Br);
250:     }

252:     /* Create subcommunicator grouping processes with same rank */
253:     if (ctx->commset) { MPI_Comm_free(&ctx->commrank); }
254:     MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
255:     MPI_Comm_split(((PetscObject)eps)->comm,rank,ctx->subc->color,&ctx->commrank);
256:     ctx->commset = PETSC_TRUE;
257:   }
258:   EPSSetType(ctx->eps,((PetscObject)eps)->type_name);

260:   /* Transfer options for ST, KSP and PC */
261:   STGetType(eps->st,&sttype);
262:   STSetType(ctx->eps->st,sttype);
263:   STGetKSP(eps->st,&ksp);
264:   KSPGetType(ksp,&ksptype);
265:   KSPGetPC(ksp,&pc);
266:   PCGetType(pc,&pctype);
267:   PCFactorGetMatSolverType(pc,&stype);
268:   PCFactorGetZeroPivot(pc,&zero);
269:   STGetKSP(ctx->eps->st,&ksp);
270:   KSPSetType(ksp,ksptype);
271:   KSPGetPC(ksp,&pc);
272:   PCSetType(pc,pctype);
273:   PCFactorSetZeroPivot(pc,zero);
274:   if (stype) { PCFactorSetMatSolverType(pc,stype); }

276:   EPSSetConvergenceTest(ctx->eps,eps->conv);
277:   EPSSetInterval(ctx->eps,a,b);
278:   ctx_local = (EPS_KRYLOVSCHUR*)ctx->eps->data;
279:   ctx_local->npart = ctx->npart;
280:   ctx_local->detect = ctx->detect;
281:   ctx_local->global = PETSC_FALSE;
282:   ctx_local->eps = eps;
283:   ctx_local->subc = ctx->subc;
284:   ctx_local->commrank = ctx->commrank;

286:   EPSSetDimensions(ctx->eps,ctx->nev,ctx->ncv,ctx->mpd);
287:   EPSKrylovSchurSetLocking(ctx->eps,ctx->lock);

289:   /* transfer options from eps->V */
290:   EPSGetBV(ctx->eps,&V);
291:   if (!eps->V) { EPSGetBV(eps,&eps->V); }
292:   if (!((PetscObject)(eps->V))->type_name) {
293:     BVSetType(V,BVSVEC);
294:   } else {
295:     BVGetType(eps->V,&type);
296:     BVSetType(V,type);
297:   }
298:   BVGetOrthogonalization(eps->V,&orthog_type,&orthog_ref,&eta,&ob_type);
299:   BVSetOrthogonalization(V,orthog_type,orthog_ref,eta,ob_type);
300:   ctx->eps->which = eps->which;
301:   ctx->eps->max_it = eps->max_it;
302:   ctx->eps->tol = eps->tol;
303:   ctx->eps->purify = eps->purify;
304:   if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
305:   EPSSetProblemType(ctx->eps,eps->problem_type);
306:   EPSSetUp(ctx->eps);
307:   ctx->eps->nconv = 0;
308:   ctx->eps->its   = 0;
309:   for (i=0;i<ctx->eps->ncv;i++) {
310:     ctx->eps->eigr[i]   = 0.0;
311:     ctx->eps->eigi[i]   = 0.0;
312:     ctx->eps->errest[i] = 0.0;
313:   }
314:   return(0);
315: }

317: static PetscErrorCode EPSSliceGetInertia(EPS eps,PetscReal shift,PetscInt *inertia,PetscInt *zeros)
318: {
320:   KSP            ksp;
321:   PC             pc;
322:   Mat            F;
323:   PetscReal      nzshift;

326:   if (shift >= PETSC_MAX_REAL) { /* Right-open interval */
327:     if (inertia) *inertia = eps->n;
328:   } else if (shift <= PETSC_MIN_REAL) {
329:     if (inertia) *inertia = 0;
330:     if (zeros) *zeros = 0;
331:   } else {
332:     /* If the shift is zero, perturb it to a very small positive value.
333:        The goal is that the nonzero pattern is the same in all cases and reuse
334:        the symbolic factorizations */
335:     nzshift = (shift==0.0)? 10.0/PETSC_MAX_REAL: shift;
336:     STSetShift(eps->st,nzshift);
337:     STSetUp(eps->st);
338:     STGetKSP(eps->st,&ksp);
339:     KSPGetPC(ksp,&pc);
340:     PCFactorGetMatrix(pc,&F);
341:     MatGetInertia(F,inertia,zeros,NULL);
342:   }
343:   return(0);
344: }

346: /*
347:    Dummy backtransform operation
348:  */
349: static PetscErrorCode EPSBackTransform_Skip(EPS eps)
350: {
352:   return(0);
353: }

355: PetscErrorCode EPSSetUp_KrylovSchur_Slice(EPS eps)
356: {
357:   PetscErrorCode  ierr;
358:   PetscBool       issinv;
359:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data,*ctx_glob;
360:   EPS_SR          sr,sr_loc,sr_glob;
361:   PetscInt        nEigs,dssz=1,i,zeros=0,off=0,method;
362:   PetscMPIInt     nproc,rank=0,aux;
363:   PetscReal       r;
364:   MPI_Request     req;
365:   Mat             A,B=NULL;
366:   SlepcSC         sc;
367:   PetscInt        flg=0;
368:   DSParallelType  ptype;

371:   if (ctx->global) {
372:     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");
373:     if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
374:     if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs cannot be used with spectrum slicing");
375:     PetscObjectTypeCompareAny((PetscObject)eps->st,&issinv,STSINVERT,STCAYLEY,"");
376:     if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
377:     if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL*1e-2;  /* use tighter tolerance */
378:     if (!eps->max_it) eps->max_it = 100;
379:     if (ctx->nev==1) ctx->nev = PetscMin(40,eps->n);  /* nev not set, use default value */
380:     if (eps->n>10 && ctx->nev<10) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"nev cannot be less than 10 in spectrum slicing runs");
381:   }
382:   eps->ops->backtransform = EPSBackTransform_Skip;

384:   /* create spectrum slicing context and initialize it */
385:   EPSSliceResetSR(eps);
386:   PetscNewLog(eps,&sr);
387:   ctx->sr = sr;
388:   sr->itsKs = 0;
389:   sr->nleap = 0;
390:   sr->nMAXCompl = eps->nev/4;
391:   sr->iterCompl = eps->max_it/4;
392:   sr->sPres = NULL;
393:   sr->nS = 0;

395:   if (ctx->npart==1 || ctx->global) {
396:     /* check presence of ends and finding direction */
397:     if ((eps->inta > PETSC_MIN_REAL && !(ctx->subintervals && ctx->subintervals[0]==ctx->subintervals[1])) || eps->intb >= PETSC_MAX_REAL) {
398:       sr->int0 = eps->inta;
399:       sr->int1 = eps->intb;
400:       sr->dir = 1;
401:       if (eps->intb >= PETSC_MAX_REAL) { /* Right-open interval */
402:         sr->hasEnd = PETSC_FALSE;
403:       } else sr->hasEnd = PETSC_TRUE;
404:     } else {
405:       sr->int0 = eps->intb;
406:       sr->int1 = eps->inta;
407:       sr->dir = -1;
408:       sr->hasEnd = PetscNot(eps->inta <= PETSC_MIN_REAL);
409:     }
410:   }
411:   if (ctx->global) {
412:     /* prevent computation of factorization in global eps */
413:     STSetTransform(eps->st,PETSC_FALSE);
414:     EPSSetDimensions_Default(eps,ctx->nev,&ctx->ncv,&ctx->mpd);
415:     /* create subintervals and initialize auxiliary eps for slicing runs */
416:     EPSSliceGetEPS(eps);
417:     sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
418:     if (ctx->npart>1) {
419:       if ((sr->dir>0&&ctx->subc->color==0)||(sr->dir<0&&ctx->subc->color==ctx->npart-1)) sr->inertia0 = sr_loc->inertia0;
420:       MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
421:       if (!rank) {
422:         MPI_Bcast(&sr->inertia0,1,MPIU_INT,(sr->dir>0)?0:ctx->npart-1,ctx->commrank);
423:       }
424:       MPI_Bcast(&sr->inertia0,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
425:       PetscFree(ctx->nconv_loc);
426:       PetscMalloc1(ctx->npart,&ctx->nconv_loc);
427:       MPI_Comm_size(((PetscObject)eps)->comm,&nproc);
428:       if (sr->dir<0) off = 1;
429:       if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
430:         PetscMPIIntCast(sr_loc->numEigs,&aux);
431:         MPI_Allgather(&aux,1,MPI_INT,ctx->nconv_loc,1,MPI_INT,ctx->commrank);
432:         MPI_Allgather(sr_loc->dir==sr->dir?&sr_loc->int0:&sr_loc->int1,1,MPIU_REAL,ctx->subintervals+off,1,MPIU_REAL,ctx->commrank);
433:       } else {
434:         MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
435:         if (!rank) {
436:           PetscMPIIntCast(sr_loc->numEigs,&aux);
437:           MPI_Allgather(&aux,1,MPI_INT,ctx->nconv_loc,1,MPI_INT,ctx->commrank);
438:           MPI_Allgather(sr_loc->dir==sr->dir?&sr_loc->int0:&sr_loc->int1,1,MPIU_REAL,ctx->subintervals+off,1,MPIU_REAL,ctx->commrank);
439:         }
440:         PetscMPIIntCast(ctx->npart,&aux);
441:         MPI_Bcast(ctx->nconv_loc,aux,MPI_INT,0,PetscSubcommChild(ctx->subc));
442:         MPI_Bcast(ctx->subintervals+off,aux,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
443:       }
444:       nEigs = 0;
445:       for (i=0;i<ctx->npart;i++) nEigs += ctx->nconv_loc[i];
446:     } else {
447:       nEigs = sr_loc->numEigs;
448:       sr->inertia0 = sr_loc->inertia0;
449:     }
450:     sr->inertia1 = sr->inertia0+sr->dir*nEigs;
451:     sr->numEigs = nEigs;
452:     eps->nev = nEigs;
453:     eps->ncv = nEigs;
454:     eps->mpd = nEigs;
455:   } else {
456:     ctx_glob = (EPS_KRYLOVSCHUR*)ctx->eps->data;
457:     sr_glob = ctx_glob->sr;
458:     if (ctx->npart>1) {
459:       sr->dir = sr_glob->dir;
460:       sr->int0 = (sr->dir==1)?eps->inta:eps->intb;
461:       sr->int1 = (sr->dir==1)?eps->intb:eps->inta;
462:       if ((sr->dir>0&&ctx->subc->color==ctx->npart-1)||(sr->dir<0&&ctx->subc->color==0)) sr->hasEnd = sr_glob->hasEnd;
463:       else sr->hasEnd = PETSC_TRUE;
464:     }

466:     /* compute inertia0 */
467:     EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,ctx->detect?&zeros:NULL);
468:     PetscOptionsGetInt(NULL,NULL,"-eps_krylovschur_hiteigenvalue",&flg,NULL);
469:     if (zeros) { /* error in factorization */
470:       if (sr->int0==ctx->eps->inta || sr->int0==ctx->eps->intb) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in the interval endpoint");
471:       else if(ctx_glob->subintset && !flg) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
472:       else {
473:         if (flg==1) { /* idle subgroup */
474:           sr->inertia0 = -1;
475:         } else { /* perturb shift */
476:           sr->int0 *= (1.0+SLICE_PTOL);
477:           EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,&zeros);
478:           if (zeros) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",sr->int1);
479:         }
480:       }
481:     }
482:     if (ctx->npart>1) {
483:       /* inertia1 is received from neighbour */
484:       MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
485:       if (!rank) {
486:         if ( sr->inertia0!=-1 && ((sr->dir>0 && ctx->subc->color>0) || (sr->dir<0 && ctx->subc->color<ctx->npart-1)) ) { /* send inertia0 to neighbour0 */
487:           MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
488:           MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
489:         }
490:         if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)|| (sr->dir<0 && ctx->subc->color>0)) { /* receive inertia1 from neighbour1 */
491:           MPI_Recv(&(sr->inertia1),1,MPIU_INT,ctx->subc->color+sr->dir,0,ctx->commrank,MPI_STATUS_IGNORE);
492:           MPI_Recv(&(sr->int1),1,MPIU_REAL,ctx->subc->color+sr->dir,0,ctx->commrank,MPI_STATUS_IGNORE);
493:         }
494:         if ( sr->inertia0==-1 && !(sr->dir>0 && ctx->subc->color==ctx->npart-1) && !(sr->dir<0 && ctx->subc->color==0)) {
495:           sr->inertia0 = sr->inertia1; sr->int0 = sr->int1;
496:           MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
497:           MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
498:         }
499:       }
500:       if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)||(sr->dir<0 && ctx->subc->color>0)) {
501:         MPI_Bcast(&sr->inertia1,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
502:         MPI_Bcast(&sr->int1,1,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
503:       } else sr_glob->inertia1 = sr->inertia1;
504:     }

506:     /* last process in eps comm computes inertia1 */
507:     if (ctx->npart==1 || ((sr->dir>0 && ctx->subc->color==ctx->npart-1) || (sr->dir<0 && ctx->subc->color==0))) {
508:       EPSSliceGetInertia(eps,sr->int1,&sr->inertia1,ctx->detect?&zeros:NULL);
509:       if (zeros) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
510:       if (!rank && sr->inertia0==-1) {
511:         sr->inertia0 = sr->inertia1; sr->int0 = sr->int1;
512:         MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
513:         MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
514:       }
515:       if (sr->hasEnd) {
516:         sr->dir = -sr->dir; r = sr->int0; sr->int0 = sr->int1; sr->int1 = r;
517:         i = sr->inertia0; sr->inertia0 = sr->inertia1; sr->inertia1 = i;
518:       }
519:     }

521:     /* number of eigenvalues in interval */
522:     sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
523:     if (ctx->npart>1) {
524:       /* memory allocate for subinterval eigenpairs */
525:       EPSSliceAllocateSolution(eps,1);
526:     }
527:     dssz = eps->ncv+1;
528:     if (sr->numEigs>0) {
529:       DSGetSlepcSC(eps->ds,&sc);
530:       sc->rg            = NULL;
531:       sc->comparison    = SlepcCompareLargestMagnitude;
532:       sc->comparisonctx = NULL;
533:       sc->map           = NULL;
534:       sc->mapobj        = NULL;
535:     }
536:     DSGetParallel(ctx->eps->ds,&ptype);
537:     DSSetParallel(eps->ds,ptype);
538:     DSGetMethod(ctx->eps->ds,&method);
539:     DSSetMethod(eps->ds,method);
540:   }
541:   DSSetType(eps->ds,DSHEP);
542:   DSSetCompact(eps->ds,PETSC_TRUE);
543:   DSAllocate(eps->ds,dssz);
544:   /* keep state of subcomm matrices to check that the user does not modify them */
545:   EPSGetOperators(eps,&A,&B);
546:   PetscObjectStateGet((PetscObject)A,&ctx->Astate);
547:   PetscObjectGetId((PetscObject)A,&ctx->Aid);
548:   if (B) {
549:     PetscObjectStateGet((PetscObject)B,&ctx->Bstate);
550:     PetscObjectGetId((PetscObject)B,&ctx->Bid);
551:   } else {
552:     ctx->Bstate=0;
553:     ctx->Bid=0;
554:   }
555:   return(0);
556: }

558: static PetscErrorCode EPSSliceGatherEigenVectors(EPS eps)
559: {
560:   PetscErrorCode  ierr;
561:   Vec             v,vg,v_loc;
562:   IS              is1,is2;
563:   VecScatter      vec_sc;
564:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
565:   PetscInt        nloc,m0,n0,i,si,idx,*idx1,*idx2,j;
566:   PetscScalar     *array;
567:   EPS_SR          sr_loc;
568:   BV              V_loc;

571:   sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
572:   V_loc = sr_loc->V;

574:   /* Gather parallel eigenvectors */
575:   BVGetColumn(eps->V,0,&v);
576:   VecGetOwnershipRange(v,&n0,&m0);
577:   BVRestoreColumn(eps->V,0,&v);
578:   BVGetColumn(ctx->eps->V,0,&v);
579:   VecGetLocalSize(v,&nloc);
580:   BVRestoreColumn(ctx->eps->V,0,&v);
581:   PetscMalloc2(m0-n0,&idx1,m0-n0,&idx2);
582:   VecCreateMPI(PetscObjectComm((PetscObject)eps),nloc,PETSC_DECIDE,&vg);
583:   idx = -1;
584:   for (si=0;si<ctx->npart;si++) {
585:     j = 0;
586:     for (i=n0;i<m0;i++) {
587:       idx1[j]   = i;
588:       idx2[j++] = i+eps->n*si;
589:     }
590:     ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx1,PETSC_COPY_VALUES,&is1);
591:     ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx2,PETSC_COPY_VALUES,&is2);
592:     BVGetColumn(eps->V,0,&v);
593:     VecScatterCreate(v,is1,vg,is2,&vec_sc);
594:     BVRestoreColumn(eps->V,0,&v);
595:     ISDestroy(&is1);
596:     ISDestroy(&is2);
597:     for (i=0;i<ctx->nconv_loc[si];i++) {
598:       BVGetColumn(eps->V,++idx,&v);
599:       if (ctx->subc->color==si) {
600:         BVGetColumn(V_loc,i,&v_loc);
601:         VecGetArray(v_loc,&array);
602:         VecPlaceArray(vg,array);
603:       }
604:       VecScatterBegin(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE);
605:       VecScatterEnd(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE);
606:       if (ctx->subc->color==si) {
607:         VecResetArray(vg);
608:         VecRestoreArray(v_loc,&array);
609:         BVRestoreColumn(V_loc,i,&v_loc);
610:       }
611:       BVRestoreColumn(eps->V,idx,&v);
612:     }
613:     VecScatterDestroy(&vec_sc);
614:   }
615:   PetscFree2(idx1,idx2);
616:   VecDestroy(&vg);
617:   return(0);
618: }

620: /*
621:   EPSComputeVectors_Slice - Recover Eigenvectors from subcomunicators
622:  */
623: PetscErrorCode EPSComputeVectors_Slice(EPS eps)
624: {
625:   PetscErrorCode  ierr;
626:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

629:   if (ctx->global && ctx->npart>1) {
630:     EPSComputeVectors(ctx->eps);
631:     EPSSliceGatherEigenVectors(eps);
632:   }
633:   return(0);
634: }

636: #define SWAP(a,b,t) {t=a;a=b;b=t;}

638: static PetscErrorCode EPSSliceGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
639: {
640:   PetscErrorCode  ierr;
641:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
642:   PetscInt        i=0,j,tmpi;
643:   PetscReal       v,tmpr;
644:   EPS_shift       s;

647:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
648:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
649:   if (!ctx->sr->s0) {  /* EPSSolve not called yet */
650:     *n = 2;
651:   } else {
652:     *n = 1;
653:     s = ctx->sr->s0;
654:     while (s) {
655:       (*n)++;
656:       s = s->neighb[1];
657:     }
658:   }
659:   PetscMalloc1(*n,shifts);
660:   PetscMalloc1(*n,inertias);
661:   if (!ctx->sr->s0) {  /* EPSSolve not called yet */
662:     (*shifts)[0]   = ctx->sr->int0;
663:     (*shifts)[1]   = ctx->sr->int1;
664:     (*inertias)[0] = ctx->sr->inertia0;
665:     (*inertias)[1] = ctx->sr->inertia1;
666:   } else {
667:     s = ctx->sr->s0;
668:     while (s) {
669:       (*shifts)[i]     = s->value;
670:       (*inertias)[i++] = s->inertia;
671:       s = s->neighb[1];
672:     }
673:     (*shifts)[i]   = ctx->sr->int1;
674:     (*inertias)[i] = ctx->sr->inertia1;
675:   }
676:   /* remove possible duplicate in last position */
677:   if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
678:   /* sort result */
679:   for (i=0;i<*n;i++) {
680:     v = (*shifts)[i];
681:     for (j=i+1;j<*n;j++) {
682:       if (v > (*shifts)[j]) {
683:         SWAP((*shifts)[i],(*shifts)[j],tmpr);
684:         SWAP((*inertias)[i],(*inertias)[j],tmpi);
685:         v = (*shifts)[i];
686:       }
687:     }
688:   }
689:   return(0);
690: }

692: static PetscErrorCode EPSSliceGatherSolution(EPS eps)
693: {
694:   PetscErrorCode  ierr;
695:   PetscMPIInt     rank,nproc;
696:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
697:   PetscInt        i,idx,j;
698:   PetscInt        *perm_loc,off=0,*inertias_loc,ns;
699:   PetscScalar     *eigr_loc;
700:   EPS_SR          sr_loc;
701:   PetscReal       *shifts_loc;
702:   PetscMPIInt     *disp,*ns_loc,aux;

705:   eps->nconv = 0;
706:   for (i=0;i<ctx->npart;i++) eps->nconv += ctx->nconv_loc[i];
707:   sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

709:   /* Gather the shifts used and the inertias computed */
710:   EPSSliceGetInertias(ctx->eps,&ns,&shifts_loc,&inertias_loc);
711:   if (ctx->sr->dir>0 && shifts_loc[ns-1]==sr_loc->int1 && ctx->subc->color<ctx->npart-1) ns--;
712:   if (ctx->sr->dir<0 && shifts_loc[ns-1]==sr_loc->int0 && ctx->subc->color>0) {
713:     ns--;
714:     for (i=0;i<ns;i++) {
715:       inertias_loc[i] = inertias_loc[i+1];
716:       shifts_loc[i] = shifts_loc[i+1];
717:     }
718:   }
719:   PetscMalloc1(ctx->npart,&ns_loc);
720:   MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
721:   PetscMPIIntCast(ns,&aux);
722:   if (!rank) { MPI_Allgather(&aux,1,MPI_INT,ns_loc,1,MPI_INT,ctx->commrank); }
723:   PetscMPIIntCast(ctx->npart,&aux);
724:   MPI_Bcast(ns_loc,aux,MPI_INT,0,PetscSubcommChild(ctx->subc));
725:   ctx->nshifts = 0;
726:   for (i=0;i<ctx->npart;i++) ctx->nshifts += ns_loc[i];
727:   PetscFree(ctx->inertias);
728:   PetscFree(ctx->shifts);
729:   PetscMalloc1(ctx->nshifts,&ctx->inertias);
730:   PetscMalloc1(ctx->nshifts,&ctx->shifts);

732:   /* Gather eigenvalues (same ranks have fully set of eigenvalues)*/
733:   eigr_loc = sr_loc->eigr;
734:   perm_loc = sr_loc->perm;
735:   MPI_Comm_size(((PetscObject)eps)->comm,&nproc);
736:   PetscMalloc1(ctx->npart,&disp);
737:   disp[0] = 0;
738:   for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ctx->nconv_loc[i-1];
739:   if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
740:     PetscMPIIntCast(sr_loc->numEigs,&aux);
741:     MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank); /* eigenvalues */
742:     MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank); /* perm */
743:     for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
744:     PetscMPIIntCast(ns,&aux);
745:     MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank); /* shifts */
746:     MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank); /* inertias */
747:     MPI_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank);
748:   } else { /* subcommunicators with different size */
749:     MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
750:     if (!rank) {
751:       PetscMPIIntCast(sr_loc->numEigs,&aux);
752:       MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank); /* eigenvalues */
753:       MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank); /* perm */
754:       for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
755:       PetscMPIIntCast(ns,&aux);
756:       MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank); /* shifts */
757:       MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank); /* inertias */
758:       MPI_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank);
759:     }
760:     PetscMPIIntCast(eps->nconv,&aux);
761:     MPI_Bcast(eps->eigr,aux,MPIU_SCALAR,0,PetscSubcommChild(ctx->subc));
762:     MPI_Bcast(eps->perm,aux,MPIU_INT,0,PetscSubcommChild(ctx->subc));
763:     MPI_Bcast(ctx->shifts,ctx->nshifts,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
764:     PetscMPIIntCast(ctx->nshifts,&aux);
765:     MPI_Bcast(ctx->inertias,aux,MPIU_INT,0,PetscSubcommChild(ctx->subc));
766:     MPI_Bcast(&eps->its,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
767:   }
768:   /* Update global array eps->perm */
769:   idx = ctx->nconv_loc[0];
770:   for (i=1;i<ctx->npart;i++) {
771:     off += ctx->nconv_loc[i-1];
772:     for (j=0;j<ctx->nconv_loc[i];j++) eps->perm[idx++] += off;
773:   }

775:   /* Gather parallel eigenvectors */
776:   PetscFree(ns_loc);
777:   PetscFree(disp);
778:   PetscFree(shifts_loc);
779:   PetscFree(inertias_loc);
780:   return(0);
781: }

783: /*
784:    Fills the fields of a shift structure
785: */
786: static PetscErrorCode EPSCreateShift(EPS eps,PetscReal val,EPS_shift neighb0,EPS_shift neighb1)
787: {
788:   PetscErrorCode  ierr;
789:   EPS_shift       s,*pending2;
790:   PetscInt        i;
791:   EPS_SR          sr;
792:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

795:   sr = ctx->sr;
796:   PetscNewLog(eps,&s);
797:   s->value = val;
798:   s->neighb[0] = neighb0;
799:   if (neighb0) neighb0->neighb[1] = s;
800:   s->neighb[1] = neighb1;
801:   if (neighb1) neighb1->neighb[0] = s;
802:   s->comp[0] = PETSC_FALSE;
803:   s->comp[1] = PETSC_FALSE;
804:   s->index = -1;
805:   s->neigs = 0;
806:   s->nconv[0] = s->nconv[1] = 0;
807:   s->nsch[0] = s->nsch[1]=0;
808:   /* Inserts in the stack of pending shifts */
809:   /* If needed, the array is resized */
810:   if (sr->nPend >= sr->maxPend) {
811:     sr->maxPend *= 2;
812:     PetscMalloc1(sr->maxPend,&pending2);
813:     PetscLogObjectMemory((PetscObject)eps,sizeof(EPS_shift));
814:     for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
815:     PetscFree(sr->pending);
816:     sr->pending = pending2;
817:   }
818:   sr->pending[sr->nPend++]=s;
819:   return(0);
820: }

822: /* Prepare for Rational Krylov update */
823: static PetscErrorCode EPSPrepareRational(EPS eps)
824: {
825:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
826:   PetscErrorCode  ierr;
827:   PetscInt        dir,i,k,ld,nv;
828:   PetscScalar     *A;
829:   EPS_SR          sr = ctx->sr;
830:   Vec             v;

833:   DSGetLeadingDimension(eps->ds,&ld);
834:   dir = (sr->sPres->neighb[0] == sr->sPrev)?1:-1;
835:   dir*=sr->dir;
836:   k = 0;
837:   for (i=0;i<sr->nS;i++) {
838:     if (dir*PetscRealPart(sr->S[i])>0.0) {
839:       sr->S[k] = sr->S[i];
840:       sr->S[sr->nS+k] = sr->S[sr->nS+i];
841:       BVGetColumn(sr->Vnext,k,&v);
842:       BVCopyVec(eps->V,eps->nconv+i,v);
843:       BVRestoreColumn(sr->Vnext,k,&v);
844:       k++;
845:       if (k>=sr->nS/2)break;
846:     }
847:   }
848:   /* Copy to DS */
849:   DSGetArray(eps->ds,DS_MAT_A,&A);
850:   PetscMemzero(A,ld*ld*sizeof(PetscScalar));
851:   for (i=0;i<k;i++) {
852:     A[i*(1+ld)] = sr->S[i];
853:     A[k+i*ld] = sr->S[sr->nS+i];
854:   }
855:   sr->nS = k;
856:   DSRestoreArray(eps->ds,DS_MAT_A,&A);
857:   DSGetDimensions(eps->ds,&nv,NULL,NULL,NULL,NULL);
858:   DSSetDimensions(eps->ds,nv,0,0,k);
859:   /* Append u to V */
860:   BVGetColumn(sr->Vnext,sr->nS,&v);
861:   BVCopyVec(eps->V,sr->nv,v);
862:   BVRestoreColumn(sr->Vnext,sr->nS,&v);
863:   return(0);
864: }

866: /* Provides next shift to be computed */
867: static PetscErrorCode EPSExtractShift(EPS eps)
868: {
869:   PetscErrorCode  ierr;
870:   PetscInt        iner,zeros=0;
871:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
872:   EPS_SR          sr;
873:   PetscReal       newShift;
874:   EPS_shift       sPres;

877:   sr = ctx->sr;
878:   if (sr->nPend > 0) {
879:     sr->sPrev = sr->sPres;
880:     sr->sPres = sr->pending[--sr->nPend];
881:     sPres = sr->sPres;
882:     EPSSliceGetInertia(eps,sPres->value,&iner,ctx->detect?&zeros:NULL);
883:     if (zeros) {
884:       newShift = sPres->value*(1.0+SLICE_PTOL);
885:       if (sr->dir*(sPres->neighb[0] && newShift-sPres->neighb[0]->value) < 0) newShift = (sPres->value+sPres->neighb[0]->value)/2;
886:       else if (sPres->neighb[1] && sr->dir*(sPres->neighb[1]->value-newShift) < 0) newShift = (sPres->value+sPres->neighb[1]->value)/2;
887:       EPSSliceGetInertia(eps,newShift,&iner,&zeros);
888:       if (zeros) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",newShift);
889:       sPres->value = newShift;
890:     }
891:     sr->sPres->inertia = iner;
892:     eps->target = sr->sPres->value;
893:     eps->reason = EPS_CONVERGED_ITERATING;
894:     eps->its = 0;
895:   } else sr->sPres = NULL;
896:   return(0);
897: }

899: /*
900:    Symmetric KrylovSchur adapted to spectrum slicing:
901:    Allows searching an specific amount of eigenvalues in the subintervals left and right.
902:    Returns whether the search has succeeded
903: */
904: static PetscErrorCode EPSKrylovSchur_Slice(EPS eps)
905: {
906:   PetscErrorCode  ierr;
907:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
908:   PetscInt        i,k,l,ld,nv,*iwork,j;
909:   Mat             U;
910:   PetscScalar     *Q,*A;
911:   PetscReal       *a,*b,beta;
912:   PetscBool       breakdown;
913:   PetscInt        count0,count1;
914:   PetscReal       lambda;
915:   EPS_shift       sPres;
916:   PetscBool       complIterating;
917:   PetscBool       sch0,sch1;
918:   PetscInt        iterCompl=0,n0,n1;
919:   EPS_SR          sr = ctx->sr;

922:   /* Spectrum slicing data */
923:   sPres = sr->sPres;
924:   complIterating =PETSC_FALSE;
925:   sch1 = sch0 = PETSC_TRUE;
926:   DSGetLeadingDimension(eps->ds,&ld);
927:   PetscMalloc1(2*ld,&iwork);
928:   count0=0;count1=0; /* Found on both sides */
929:   if (sr->nS > 0 && (sPres->neighb[0] == sr->sPrev || sPres->neighb[1] == sr->sPrev)) {
930:     /* Rational Krylov */
931:     DSTranslateRKS(eps->ds,sr->sPrev->value-sPres->value);
932:     DSGetDimensions(eps->ds,NULL,NULL,NULL,&l,NULL);
933:     DSSetDimensions(eps->ds,l+1,0,0,0);
934:     BVSetActiveColumns(eps->V,0,l+1);
935:     DSGetMat(eps->ds,DS_MAT_Q,&U);
936:     BVMultInPlace(eps->V,U,0,l+1);
937:     MatDestroy(&U);
938:   } else {
939:     /* Get the starting Lanczos vector */
940:     EPSGetStartVector(eps,0,NULL);
941:     l = 0;
942:   }
943:   /* Restart loop */
944:   while (eps->reason == EPS_CONVERGED_ITERATING) {
945:     eps->its++; sr->itsKs++;
946:     /* Compute an nv-step Lanczos factorization */
947:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
948:     DSGetArrayReal(eps->ds,DS_MAT_T,&a);
949:     b = a + ld;
950:     EPSFullLanczos(eps,a,b,eps->nconv+l,&nv,&breakdown);
951:     sr->nv = nv;
952:     beta = b[nv-1];
953:     DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
954:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
955:     if (l==0) {
956:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
957:     } else {
958:       DSSetState(eps->ds,DS_STATE_RAW);
959:     }
960:     BVSetActiveColumns(eps->V,eps->nconv,nv);

962:     /* Solve projected problem and compute residual norm estimates */
963:     if (eps->its == 1 && l > 0) {/* After rational update */
964:       DSGetArray(eps->ds,DS_MAT_A,&A);
965:       DSGetArrayReal(eps->ds,DS_MAT_T,&a);
966:       b = a + ld;
967:       k = eps->nconv+l;
968:       A[k*ld+k-1] = A[(k-1)*ld+k];
969:       A[k*ld+k] = a[k];
970:       for (j=k+1; j< nv; j++) {
971:         A[j*ld+j] = a[j];
972:         A[j*ld+j-1] = b[j-1] ;
973:         A[(j-1)*ld+j] = b[j-1];
974:       }
975:       DSRestoreArray(eps->ds,DS_MAT_A,&A);
976:       DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
977:       DSSolve(eps->ds,eps->eigr,NULL);
978:       DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
979:       DSSetCompact(eps->ds,PETSC_TRUE);
980:     } else { /* Restart */
981:       DSSolve(eps->ds,eps->eigr,NULL);
982:       DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
983:     }
984:     DSSynchronize(eps->ds,eps->eigr,NULL);

986:     /* Residual */
987:     EPSKrylovConvergence(eps,PETSC_TRUE,eps->nconv,nv-eps->nconv,beta,0.0,1.0,&k);
988:     /* Checking values obtained for completing */
989:     for (i=0;i<k;i++) {
990:       sr->back[i]=eps->eigr[i];
991:     }
992:     STBackTransform(eps->st,k,sr->back,eps->eigi);
993:     count0=count1=0;
994:     for (i=0;i<k;i++) {
995:       lambda = PetscRealPart(sr->back[i]);
996:       if (((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0)) count0++;
997:       if (((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0)) count1++;
998:     }
999:     if (k>eps->nev && eps->ncv-k<5) eps->reason = EPS_CONVERGED_TOL;
1000:     else {
1001:       /* Checks completion */
1002:       if ((!sch0||count0 >= sPres->nsch[0]) && (!sch1 ||count1 >= sPres->nsch[1])) {
1003:         eps->reason = EPS_CONVERGED_TOL;
1004:       } else {
1005:         if (!complIterating && eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
1006:         if (complIterating) {
1007:           if (--iterCompl <= 0) eps->reason = EPS_DIVERGED_ITS;
1008:         } else if (k >= eps->nev) {
1009:           n0 = sPres->nsch[0]-count0;
1010:           n1 = sPres->nsch[1]-count1;
1011:           if (sr->iterCompl>0 && ((n0>0 && n0<= sr->nMAXCompl)||(n1>0&&n1<=sr->nMAXCompl))) {
1012:             /* Iterating for completion*/
1013:             complIterating = PETSC_TRUE;
1014:             if (n0 >sr->nMAXCompl)sch0 = PETSC_FALSE;
1015:             if (n1 >sr->nMAXCompl)sch1 = PETSC_FALSE;
1016:             iterCompl = sr->iterCompl;
1017:           } else eps->reason = EPS_CONVERGED_TOL;
1018:         }
1019:       }
1020:     }
1021:     /* Update l */
1022:     if (eps->reason == EPS_CONVERGED_ITERATING) l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
1023:     else l = nv-k;
1024:     if (breakdown) l=0;
1025:     if (!ctx->lock && l>0 && eps->reason == EPS_CONVERGED_ITERATING) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */

1027:     if (eps->reason == EPS_CONVERGED_ITERATING) {
1028:       if (breakdown) {
1029:         /* Start a new Lanczos factorization */
1030:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
1031:         EPSGetStartVector(eps,k,&breakdown);
1032:         if (breakdown) {
1033:           eps->reason = EPS_DIVERGED_BREAKDOWN;
1034:           PetscInfo(eps,"Unable to generate more start vectors\n");
1035:         }
1036:       } else {
1037:         /* Prepare the Rayleigh quotient for restart */
1038:         DSGetArrayReal(eps->ds,DS_MAT_T,&a);
1039:         DSGetArray(eps->ds,DS_MAT_Q,&Q);
1040:         b = a + ld;
1041:         for (i=k;i<k+l;i++) {
1042:           a[i] = PetscRealPart(eps->eigr[i]);
1043:           b[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
1044:         }
1045:         DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
1046:         DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
1047:       }
1048:     }
1049:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1050:     DSGetMat(eps->ds,DS_MAT_Q,&U);
1051:     BVMultInPlace(eps->V,U,eps->nconv,k+l);
1052:     MatDestroy(&U);

1054:     /* Normalize u and append it to V */
1055:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
1056:       BVCopyColumn(eps->V,nv,k+l);
1057:     }
1058:     eps->nconv = k;
1059:     if (eps->reason != EPS_CONVERGED_ITERATING) {
1060:       /* Store approximated values for next shift */
1061:       DSGetArray(eps->ds,DS_MAT_Q,&Q);
1062:       sr->nS = l;
1063:       for (i=0;i<l;i++) {
1064:         sr->S[i] = eps->eigr[i+k];/* Diagonal elements */
1065:         sr->S[i+l] = Q[nv-1+(i+k)*ld]*beta; /* Out of diagonal elements */
1066:       }
1067:       DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
1068:     }
1069:   }
1070:   /* Check for completion */
1071:   for (i=0;i< eps->nconv; i++) {
1072:     if ((sr->dir)*PetscRealPart(eps->eigr[i])>0) sPres->nconv[1]++;
1073:     else sPres->nconv[0]++;
1074:   }
1075:   sPres->comp[0] = PetscNot(count0 < sPres->nsch[0]);
1076:   sPres->comp[1] = PetscNot(count1 < sPres->nsch[1]);
1077:   if (count0 > sPres->nsch[0] || count1 > sPres->nsch[1]) InertiaMismatch(eps,ctx->detect);
1078:   PetscFree(iwork);
1079:   return(0);
1080: }

1082: /*
1083:   Obtains value of subsequent shift
1084: */
1085: static PetscErrorCode EPSGetNewShiftValue(EPS eps,PetscInt side,PetscReal *newS)
1086: {
1087:   PetscReal       lambda,d_prev;
1088:   PetscInt        i,idxP;
1089:   EPS_SR          sr;
1090:   EPS_shift       sPres,s;
1091:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

1094:   sr = ctx->sr;
1095:   sPres = sr->sPres;
1096:   if (sPres->neighb[side]) {
1097:   /* Completing a previous interval */
1098:     if (!sPres->neighb[side]->neighb[side] && sPres->neighb[side]->nconv[side]==0) { /* One of the ends might be too far from eigenvalues */
1099:       if (side) *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[sr->indexEig-1]]))/2;
1100:       else *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[0]]))/2;
1101:     } else *newS=(sPres->value + sPres->neighb[side]->value)/2;
1102:   } else { /* (Only for side=1). Creating a new interval. */
1103:     if (sPres->neigs==0) {/* No value has been accepted*/
1104:       if (sPres->neighb[0]) {
1105:         /* Multiplying by 10 the previous distance */
1106:         *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
1107:         sr->nleap++;
1108:         /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
1109:         if (!sr->hasEnd && sr->nleap > 5) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to compute the wanted eigenvalues with open interval");
1110:       } else { /* First shift */
1111:         if (eps->nconv != 0) {
1112:           /* Unaccepted values give information for next shift */
1113:           idxP=0;/* Number of values left from shift */
1114:           for (i=0;i<eps->nconv;i++) {
1115:             lambda = PetscRealPart(eps->eigr[i]);
1116:             if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
1117:             else break;
1118:           }
1119:           /* Avoiding subtraction of eigenvalues (might be the same).*/
1120:           if (idxP>0) {
1121:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[0]))/(idxP+0.3);
1122:           } else {
1123:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[eps->nconv-1]))/(eps->nconv+0.3);
1124:           }
1125:           *newS = sPres->value + ((sr->dir)*d_prev*eps->nev)/2;
1126:         } else { /* No values found, no information for next shift */
1127:           SETERRQ(PetscObjectComm((PetscObject)eps),1,"First shift renders no information");
1128:         }
1129:       }
1130:     } else { /* Accepted values found */
1131:       sr->nleap = 0;
1132:       /* Average distance of values in previous subinterval */
1133:       s = sPres->neighb[0];
1134:       while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
1135:         s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
1136:       }
1137:       if (s) {
1138:         d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
1139:       } else { /* First shift. Average distance obtained with values in this shift */
1140:         /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
1141:         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)) {
1142:           d_prev =  PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0])))/(sPres->neigs+0.3);
1143:         } else {
1144:           d_prev = PetscAbsReal(PetscRealPart(sr->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
1145:         }
1146:       }
1147:       /* Average distance is used for next shift by adding it to value on the right or to shift */
1148:       if ((sr->dir)*(PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
1149:         *newS = PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(eps->nev))/2;
1150:       } else { /* Last accepted value is on the left of shift. Adding to shift */
1151:         *newS = sPres->value + ((sr->dir)*d_prev*(eps->nev))/2;
1152:       }
1153:     }
1154:     /* End of interval can not be surpassed */
1155:     if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
1156:   }/* of neighb[side]==null */
1157:   return(0);
1158: }

1160: /*
1161:   Function for sorting an array of real values
1162: */
1163: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
1164: {
1165:   PetscReal re;
1166:   PetscInt  i,j,tmp;

1169:   if (!prev) for (i=0;i<nr;i++) perm[i] = i;
1170:   /* Insertion sort */
1171:   for (i=1;i<nr;i++) {
1172:     re = PetscRealPart(r[perm[i]]);
1173:     j = i-1;
1174:     while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
1175:       tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
1176:     }
1177:   }
1178:   return(0);
1179: }

1181: /* Stores the pairs obtained since the last shift in the global arrays */
1182: static PetscErrorCode EPSStoreEigenpairs(EPS eps)
1183: {
1184:   PetscErrorCode  ierr;
1185:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1186:   PetscReal       lambda,err,norm;
1187:   PetscInt        i,count;
1188:   PetscBool       iscayley;
1189:   EPS_SR          sr = ctx->sr;
1190:   EPS_shift       sPres;
1191:   Vec             v,w;

1194:   sPres = sr->sPres;
1195:   sPres->index = sr->indexEig;
1196:   count = sr->indexEig;
1197:   /* Back-transform */
1198:   STBackTransform(eps->st,eps->nconv,eps->eigr,eps->eigi);
1199:   PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
1200:   /* Sort eigenvalues */
1201:   sortRealEigenvalues(eps->eigr,eps->perm,eps->nconv,PETSC_FALSE,sr->dir);
1202:   /* Values stored in global array */
1203:   for (i=0;i<eps->nconv;i++) {
1204:     lambda = PetscRealPart(eps->eigr[eps->perm[i]]);
1205:     err = eps->errest[eps->perm[i]];

1207:     if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
1208:       if (count>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error in Spectrum Slicing");
1209:       sr->eigr[count] = lambda;
1210:       sr->errest[count] = err;
1211:       /* Explicit purification */
1212:       BVGetColumn(eps->V,eps->perm[i],&w);
1213:       if (eps->purify) {
1214:         BVGetColumn(sr->V,count,&v);
1215:         STApply(eps->st,w,v);
1216:         BVRestoreColumn(sr->V,count,&v);
1217:       } else {
1218:         BVInsertVec(sr->V,count,w);
1219:       }
1220:       BVRestoreColumn(eps->V,eps->perm[i],&w);
1221:       BVNormColumn(sr->V,count,NORM_2,&norm);
1222:       BVScaleColumn(sr->V,count,1.0/norm);
1223:       count++;
1224:     }
1225:   }
1226:   sPres->neigs = count - sr->indexEig;
1227:   sr->indexEig = count;
1228:   /* Global ordering array updating */
1229:   sortRealEigenvalues(sr->eigr,sr->perm,count,PETSC_TRUE,sr->dir);
1230:   return(0);
1231: }

1233: static PetscErrorCode EPSLookForDeflation(EPS eps)
1234: {
1235:   PetscErrorCode  ierr;
1236:   PetscReal       val;
1237:   PetscInt        i,count0=0,count1=0;
1238:   EPS_shift       sPres;
1239:   PetscInt        ini,fin,k,idx0,idx1;
1240:   EPS_SR          sr;
1241:   Vec             v;
1242:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

1245:   sr = ctx->sr;
1246:   sPres = sr->sPres;

1248:   if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
1249:   else ini = 0;
1250:   fin = sr->indexEig;
1251:   /* Selection of ends for searching new values */
1252:   if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
1253:   else sPres->ext[0] = sPres->neighb[0]->value;
1254:   if (!sPres->neighb[1]) {
1255:     if (sr->hasEnd) sPres->ext[1] = sr->int1;
1256:     else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
1257:   } else sPres->ext[1] = sPres->neighb[1]->value;
1258:   /* Selection of values between right and left ends */
1259:   for (i=ini;i<fin;i++) {
1260:     val=PetscRealPart(sr->eigr[sr->perm[i]]);
1261:     /* Values to the right of left shift */
1262:     if ((sr->dir)*(val - sPres->ext[1]) < 0) {
1263:       if ((sr->dir)*(val - sPres->value) < 0) count0++;
1264:       else count1++;
1265:     } else break;
1266:   }
1267:   /* The number of values on each side are found */
1268:   if (sPres->neighb[0]) {
1269:     sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
1270:     if (sPres->nsch[0]<0) InertiaMismatch(eps,ctx->detect);
1271:   } else sPres->nsch[0] = 0;

1273:   if (sPres->neighb[1]) {
1274:     sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
1275:     if (sPres->nsch[1]<0) InertiaMismatch(eps,ctx->detect);
1276:   } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);

1278:   /* Completing vector of indexes for deflation */
1279:   idx0 = ini;
1280:   idx1 = ini+count0+count1;
1281:   k=0;
1282:   for (i=idx0;i<idx1;i++) sr->idxDef[k++]=sr->perm[i];
1283:   BVDuplicateResize(eps->V,k+eps->ncv+1,&sr->Vnext);
1284:   BVSetNumConstraints(sr->Vnext,k);
1285:   for (i=0;i<k;i++) {
1286:     BVGetColumn(sr->Vnext,-i-1,&v);
1287:     BVCopyVec(sr->V,sr->idxDef[i],v);
1288:     BVRestoreColumn(sr->Vnext,-i-1,&v);
1289:   }

1291:   /* For rational Krylov */
1292:   if (sr->nS>0 && (sr->sPrev == sr->sPres->neighb[0] || sr->sPrev == sr->sPres->neighb[1])) {
1293:     EPSPrepareRational(eps);
1294:   }
1295:   eps->nconv = 0;
1296:   /* Get rid of temporary Vnext */
1297:   BVDestroy(&eps->V);
1298:   eps->V = sr->Vnext;
1299:   sr->Vnext = NULL;
1300:   return(0);
1301: }

1303: PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS eps)
1304: {
1305:   PetscErrorCode   ierr;
1306:   PetscInt         i,lds,ti;
1307:   PetscReal        newS;
1308:   EPS_KRYLOVSCHUR  *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1309:   EPS_SR           sr=ctx->sr;
1310:   Mat              A,B=NULL;
1311:   PetscObjectState Astate,Bstate=0;
1312:   PetscObjectId    Aid,Bid=0;

1315:   PetscCitationsRegister(citation,&cited);
1316:   if (ctx->global) {
1317:     EPSSolve_KrylovSchur_Slice(ctx->eps);
1318:     ctx->eps->state = EPS_STATE_SOLVED;
1319:     eps->reason = EPS_CONVERGED_TOL;
1320:     if (ctx->npart>1) {
1321:       /* Gather solution from subsolvers */
1322:       EPSSliceGatherSolution(eps);
1323:     } else {
1324:       eps->nconv = sr->numEigs;
1325:       eps->its   = ctx->eps->its;
1326:       PetscFree(ctx->inertias);
1327:       PetscFree(ctx->shifts);
1328:       EPSSliceGetInertias(ctx->eps,&ctx->nshifts,&ctx->shifts,&ctx->inertias);
1329:     }
1330:   } else {
1331:     if (ctx->npart==1) {
1332:       sr->eigr   = ctx->eps->eigr;
1333:       sr->eigi   = ctx->eps->eigi;
1334:       sr->perm   = ctx->eps->perm;
1335:       sr->errest = ctx->eps->errest;
1336:       sr->V      = ctx->eps->V;
1337:     }
1338:     /* Check that the user did not modify subcomm matrices */
1339:     EPSGetOperators(eps,&A,&B);
1340:     PetscObjectStateGet((PetscObject)A,&Astate);
1341:     PetscObjectGetId((PetscObject)A,&Aid);
1342:     if (B) {
1343:       PetscObjectStateGet((PetscObject)B,&Bstate);
1344:       PetscObjectGetId((PetscObject)B,&Bid);
1345:     }
1346:     if (Astate!=ctx->Astate || (B && Bstate!=ctx->Bstate) || Aid!=ctx->Aid || (B && Bid!=ctx->Bid)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Subcomm matrices have been modified by user");
1347:     /* Only with eigenvalues present in the interval ...*/
1348:     if (sr->numEigs==0) {
1349:       eps->reason = EPS_CONVERGED_TOL;
1350:       return(0);
1351:     }
1352:     /* Array of pending shifts */
1353:     sr->maxPend = 100; /* Initial size */
1354:     sr->nPend = 0;
1355:     PetscMalloc1(sr->maxPend,&sr->pending);
1356:     PetscLogObjectMemory((PetscObject)eps,(sr->maxPend)*sizeof(EPS_shift));
1357:     EPSCreateShift(eps,sr->int0,NULL,NULL);
1358:     /* extract first shift */
1359:     sr->sPrev = NULL;
1360:     sr->sPres = sr->pending[--sr->nPend];
1361:     sr->sPres->inertia = sr->inertia0;
1362:     eps->target = sr->sPres->value;
1363:     sr->s0 = sr->sPres;
1364:     sr->indexEig = 0;
1365:     /* Memory reservation for auxiliary variables */
1366:     lds = PetscMin(eps->mpd,eps->ncv);
1367:     PetscCalloc1(lds*lds,&sr->S);
1368:     PetscMalloc1(eps->ncv,&sr->back);
1369:     PetscLogObjectMemory((PetscObject)eps,(sr->numEigs+2*eps->ncv)*sizeof(PetscScalar));
1370:     for (i=0;i<sr->numEigs;i++) {
1371:       sr->eigr[i]   = 0.0;
1372:       sr->eigi[i]   = 0.0;
1373:       sr->errest[i] = 0.0;
1374:       sr->perm[i]   = i;
1375:     }
1376:     /* Vectors for deflation */
1377:     PetscMalloc1(sr->numEigs,&sr->idxDef);
1378:     PetscLogObjectMemory((PetscObject)eps,sr->numEigs*sizeof(PetscInt));
1379:     sr->indexEig = 0;
1380:     /* Main loop */
1381:     while (sr->sPres) {
1382:       /* Search for deflation */
1383:       EPSLookForDeflation(eps);
1384:       /* KrylovSchur */
1385:       EPSKrylovSchur_Slice(eps);

1387:       EPSStoreEigenpairs(eps);
1388:       /* Select new shift */
1389:       if (!sr->sPres->comp[1]) {
1390:         EPSGetNewShiftValue(eps,1,&newS);
1391:         EPSCreateShift(eps,newS,sr->sPres,sr->sPres->neighb[1]);
1392:       }
1393:       if (!sr->sPres->comp[0]) {
1394:         /* Completing earlier interval */
1395:         EPSGetNewShiftValue(eps,0,&newS);
1396:         EPSCreateShift(eps,newS,sr->sPres->neighb[0],sr->sPres);
1397:       }
1398:       /* Preparing for a new search of values */
1399:       EPSExtractShift(eps);
1400:     }

1402:     /* Updating eps values prior to exit */
1403:     PetscFree(sr->S);
1404:     PetscFree(sr->idxDef);
1405:     PetscFree(sr->pending);
1406:     PetscFree(sr->back);
1407:     BVDuplicateResize(eps->V,eps->ncv+1,&sr->Vnext);
1408:     BVSetNumConstraints(sr->Vnext,0);
1409:     BVDestroy(&eps->V);
1410:     eps->V      = sr->Vnext;
1411:     eps->nconv  = sr->indexEig;
1412:     eps->reason = EPS_CONVERGED_TOL;
1413:     eps->its    = sr->itsKs;
1414:     eps->nds    = 0;
1415:     if (sr->dir<0) {
1416:       for (i=0;i<eps->nconv/2;i++) {
1417:         ti = sr->perm[i]; sr->perm[i] = sr->perm[eps->nconv-1-i]; sr->perm[eps->nconv-1-i] = ti;
1418:       }
1419:     }
1420:   }
1421:   return(0);
1422: }