Actual source code: ks-slice.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*

  3:    SLEPc eigensolver: "krylovschur"

  5:    Method: Krylov-Schur with spectrum slicing for symmetric eigenproblems

  7:    References:

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

 13:        [2] C. Campos and J.E. Roman, "Spectrum slicing strategies based
 14:            on restarted Lanczos methods", Numer. Algor. 60(2):279-295,
 15:            2012.

 17:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 18:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 19:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

 21:    This file is part of SLEPc.

 23:    SLEPc is free software: you can redistribute it and/or modify it under  the
 24:    terms of version 3 of the GNU Lesser General Public License as published by
 25:    the Free Software Foundation.

 27:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 28:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 29:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 30:    more details.

 32:    You  should have received a copy of the GNU Lesser General  Public  License
 33:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 34:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 35: */

 37: #include <slepc/private/epsimpl.h>
 38:  #include krylovschur.h

 40: static PetscBool  cited = PETSC_FALSE;
 41: static const char citation[] =
 42:   "@Article{slepc-slice,\n"
 43:   "   author = \"C. Campos and J. E. Roman\",\n"
 44:   "   title = \"Strategies for spectrum slicing based on restarted {Lanczos} methods\",\n"
 45:   "   journal = \"Numer. Algorithms\",\n"
 46:   "   volume = \"60\",\n"
 47:   "   number = \"2\",\n"
 48:   "   pages = \"279--295\",\n"
 49:   "   year = \"2012,\"\n"
 50:   "   doi = \"http://dx.doi.org/10.1007/s11075-012-9564-z\"\n"
 51:   "}\n";

 53: #define SLICE_PTOL PETSC_SQRT_MACHINE_EPSILON

 57: static PetscErrorCode EPSSliceResetSR(EPS eps) {
 58:   PetscErrorCode  ierr;
 59:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
 60:   EPS_SR          sr=ctx->sr;
 61:   EPS_shift       s;

 64:   if (sr) {
 65:     if (ctx->npart>1) {
 66:       BVDestroy(&sr->V);
 67:       PetscFree4(sr->eigr,sr->eigi,sr->errest,sr->perm);
 68:     }
 69:     /* Reviewing list of shifts to free memory */
 70:     s = sr->s0;
 71:     if (s) {
 72:       while (s->neighb[1]) {
 73:         s = s->neighb[1];
 74:         PetscFree(s->neighb[0]);
 75:       }
 76:       PetscFree(s);
 77:     }
 78:     PetscFree(sr);
 79:   }
 80:   ctx->sr = NULL;
 81:   return(0);
 82: }

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

 92:   if (!ctx->global) return(0);
 93:   /* Destroy auxiliary EPS */
 94:   EPSSliceResetSR(ctx->eps);
 95:   EPSDestroy(&ctx->eps);
 96:   if (ctx->npart>1) {
 97:     PetscSubcommDestroy(&ctx->subc);
 98:     if (ctx->commset) {
 99:       MPI_Comm_free(&ctx->commrank);
100:       ctx->commset = PETSC_FALSE;
101:     }
102:   }
103:   PetscFree(ctx->subintervals);
104:   PetscFree(ctx->nconv_loc);
105:   EPSSliceResetSR(eps);
106:   PetscFree(ctx->inertias);
107:   PetscFree(ctx->shifts);
108:   if (ctx->npart>1) {
109:     ISDestroy(&ctx->isrow);
110:     ISDestroy(&ctx->iscol);
111:     MatDestroyMatrices(1,&ctx->submata);
112:     MatDestroyMatrices(1,&ctx->submatb);
113:   }
114:   return(0);
115: }

119: /*
120:   EPSSliceAllocateSolution - Allocate memory storage for common variables such
121:   as eigenvalues and eigenvectors. The argument extra is used for methods
122:   that require a working basis slightly larger than ncv.
123: */
124: static PetscErrorCode EPSSliceAllocateSolution(EPS eps,PetscInt extra)
125: {
126:   PetscErrorCode     ierr;
127:   EPS_KRYLOVSCHUR    *ctx=(EPS_KRYLOVSCHUR*)eps->data;
128:   PetscReal          eta;
129:   PetscInt           k;
130:   PetscLogDouble     cnt;
131:   BVType             type;
132:   BVOrthogType       orthog_type;
133:   BVOrthogRefineType orthog_ref;
134:   BVOrthogBlockType  ob_type;
135:   Mat                matrix;
136:   Vec                t;
137:   EPS_SR             sr = ctx->sr;

140:   /* allocate space for eigenvalues and friends */
141:   k = PetscMax(1,sr->numEigs);
142:   PetscFree4(sr->eigr,sr->eigi,sr->errest,sr->perm);
143:   PetscMalloc4(k,&sr->eigr,k,&sr->eigi,k,&sr->errest,k,&sr->perm);
144:   cnt = 2*k*sizeof(PetscScalar) + 2*k*sizeof(PetscReal) + k*sizeof(PetscInt);
145:   PetscLogObjectMemory((PetscObject)eps,cnt);

147:   /* allocate sr->V and transfer options from eps->V */
148:   BVDestroy(&sr->V);
149:   BVCreate(PetscObjectComm((PetscObject)eps),&sr->V);
150:   PetscLogObjectParent((PetscObject)eps,(PetscObject)sr->V);
151:   if (!eps->V) { EPSGetBV(eps,&eps->V); }
152:   if (!((PetscObject)(eps->V))->type_name) {
153:     BVSetType(sr->V,BVSVEC);
154:   } else {
155:     BVGetType(eps->V,&type);
156:     BVSetType(sr->V,type);
157:   }
158:   STMatCreateVecs(eps->st,&t,NULL);
159:   BVSetSizesFromVec(sr->V,t,k);
160:   VecDestroy(&t);
161:   EPS_SetInnerProduct(eps);
162:   BVGetMatrix(eps->V,&matrix,NULL);
163:   BVSetMatrix(sr->V,matrix,PETSC_FALSE);
164:   BVGetOrthogonalization(eps->V,&orthog_type,&orthog_ref,&eta,&ob_type);
165:   BVSetOrthogonalization(sr->V,orthog_type,orthog_ref,eta,ob_type);
166:   return(0);
167: }

171: static PetscErrorCode EPSSliceGetEPS(EPS eps)
172: {
173:   PetscErrorCode     ierr;
174:   EPS_KRYLOVSCHUR    *ctx=(EPS_KRYLOVSCHUR*)eps->data,*ctx_local;
175:   BV                 V;
176:   BVType             type;
177:   PetscReal          eta;
178:   BVOrthogType       orthog_type;
179:   BVOrthogRefineType orthog_ref;
180:   BVOrthogBlockType  ob_type;
181:   Mat                A,B=NULL,Ar,Br=NULL;
182:   PetscInt           i;
183:   PetscReal          h,a,b;
184:   PetscMPIInt        rank;
185:   EPS_SR             sr=ctx->sr;
186:   PC                 pc;
187:   PCType             pctype;
188:   KSP                ksp;
189:   KSPType            ksptype;
190:   STType             sttype;
191:   PetscObjectState   Astate,Bstate=0;
192:   PetscObjectId      Aid,Bid=0;
193:   const MatSolverPackage stype;

196:   EPSGetOperators(eps,&A,&B);
197:   if (ctx->npart==1) {
198:     if (!ctx->eps) { EPSCreate(((PetscObject)eps)->comm,&ctx->eps); }
199:     EPSSetType(ctx->eps,((PetscObject)eps)->type_name);
200:     EPSSetOperators(ctx->eps,A,B);
201:     a = eps->inta; b = eps->intb;
202:   } else {
203:     PetscObjectStateGet((PetscObject)A,&Astate);
204:     PetscObjectGetId((PetscObject)A,&Aid);
205:     if (B) {
206:       PetscObjectStateGet((PetscObject)B,&Bstate);
207:       PetscObjectGetId((PetscObject)B,&Bid);
208:     }
209:     if (!ctx->subc) {
210:       /* Create context for subcommunicators */
211:       PetscSubcommCreate(PetscObjectComm((PetscObject)eps),&ctx->subc);
212:       PetscSubcommSetNumber(ctx->subc,ctx->npart);
213:       PetscSubcommSetType(ctx->subc,PETSC_SUBCOMM_CONTIGUOUS);
214:       PetscLogObjectMemory((PetscObject)eps,sizeof(PetscSubcomm));

216:       /* Duplicate matrices */
217:       MatCreateRedundantMatrix(A,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Ar);
218:       ctx->Astate = Astate;
219:       ctx->Aid = Aid;
220:       if (B) {
221:         MatCreateRedundantMatrix(B,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Br);
222:         ctx->Bstate = Bstate;
223:         ctx->Bid = Bid;
224:       }
225:     } else {
226:       if (ctx->Astate != Astate || (B && ctx->Bstate != Bstate) || ctx->Aid != Aid || (B && ctx->Bid != Bid)) {
227:         EPSGetOperators(ctx->eps,&Ar,&Br);
228:         MatCreateRedundantMatrix(A,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Ar);
229:         ctx->Astate = Astate;
230:         ctx->Aid = Aid;
231:         if (B) {
232:           MatCreateRedundantMatrix(B,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Br);
233:           ctx->Bstate = Bstate;
234:           ctx->Bid = Bid;
235:         }
236:         EPSSetOperators(ctx->eps,Ar,Br);
237:         MatDestroy(&Ar);
238:         MatDestroy(&Br);
239:       }
240:     }

242:     /* Determine subintervals */
243:     if (!ctx->subintset) { /* uniform distribution if no set by user */
244:       if (!sr->hasEnd) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Global interval must be bounded for splitting it in uniform subintervals");
245:       h = (eps->intb-eps->inta)/ctx->npart;
246:       a = eps->inta+ctx->subc->color*h;
247:       b = (ctx->subc->color==ctx->npart-1)?eps->intb:eps->inta+(ctx->subc->color+1)*h;
248:       PetscFree(ctx->subintervals);
249:       PetscMalloc1(ctx->npart+1,&ctx->subintervals);
250:       for (i=0;i<ctx->npart;i++) ctx->subintervals[i] = eps->inta+h*i;
251:       ctx->subintervals[ctx->npart] = eps->intb;
252:     } else {
253:       a = ctx->subintervals[ctx->subc->color];
254:       b = ctx->subintervals[ctx->subc->color+1];
255:     }

257:     if (!ctx->eps) {
258:       /* Create auxiliary EPS */
259:       EPSCreate(PetscSubcommChild(ctx->subc),&ctx->eps);
260:       EPSSetOperators(ctx->eps,Ar,Br);
261:       MatDestroy(&Ar);
262:       MatDestroy(&Br);
263:     }

265:     /* Create subcommunicator grouping processes with same rank */
266:     if (ctx->commset) { MPI_Comm_free(&ctx->commrank); }
267:     MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
268:     MPI_Comm_split(((PetscObject)eps)->comm,rank,ctx->subc->color,&ctx->commrank);
269:     ctx->commset = PETSC_TRUE;
270:   }
271:   EPSSetType(ctx->eps,((PetscObject)eps)->type_name);

273:   /* Transfer options for ST, KSP and PC */
274:   STGetType(eps->st,&sttype);
275:   STSetType(ctx->eps->st,sttype);
276:   STGetKSP(eps->st,&ksp);
277:   KSPGetType(ksp,&ksptype);
278:   KSPGetPC(ksp,&pc);
279:   PCGetType(pc,&pctype);
280:   PCFactorGetMatSolverPackage(pc,&stype);
281:   STGetKSP(ctx->eps->st,&ksp);
282:   KSPSetType(ksp,ksptype);
283:   KSPGetPC(ksp,&pc);
284:   PCSetType(pc,pctype);
285:   if (stype) { PCFactorSetMatSolverPackage(pc,stype); }

287:   EPSSetConvergenceTest(ctx->eps,eps->conv);
288:   EPSSetInterval(ctx->eps,a,b);
289:   ctx_local = (EPS_KRYLOVSCHUR*)ctx->eps->data;
290:   ctx_local->npart = ctx->npart;
291:   ctx_local->detect = ctx->detect;
292:   ctx_local->global = PETSC_FALSE;
293:   ctx_local->eps = eps;
294:   ctx_local->subc = ctx->subc;
295:   ctx_local->commrank = ctx->commrank;

297:   EPSSetDimensions(ctx->eps,ctx->nev,ctx->ncv,ctx->mpd);
298:   EPSKrylovSchurSetLocking(ctx->eps,ctx->lock);

300:   /* transfer options from eps->V */
301:   EPSGetBV(ctx->eps,&V);
302:   if (!eps->V) { EPSGetBV(eps,&eps->V); }
303:   if (!((PetscObject)(eps->V))->type_name) {
304:     BVSetType(V,BVSVEC);
305:   } else {
306:     BVGetType(eps->V,&type);
307:     BVSetType(V,type);
308:   }
309:   BVGetOrthogonalization(eps->V,&orthog_type,&orthog_ref,&eta,&ob_type);
310:   BVSetOrthogonalization(V,orthog_type,orthog_ref,eta,ob_type);
311:   ctx->eps->which = eps->which;
312:   ctx->eps->max_it = eps->max_it;
313:   ctx->eps->tol = eps->tol;
314:   ctx->eps->purify = eps->purify;
315:   if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
316:   EPSSetProblemType(ctx->eps,eps->problem_type);
317:   EPSSetUp(ctx->eps);
318:   ctx->eps->nconv = 0;
319:   ctx->eps->its   = 0;
320:   for (i=0;i<ctx->eps->ncv;i++) {
321:     ctx->eps->eigr[i]   = 0.0;
322:     ctx->eps->eigi[i]   = 0.0;
323:     ctx->eps->errest[i] = 0.0;
324:   }
325:   return(0);
326: }

330: static PetscErrorCode EPSSliceGetInertia(EPS eps,PetscReal shift,PetscInt *inertia,PetscInt *zeros)
331: {
333:   KSP            ksp;
334:   PC             pc;
335:   Mat            F;
336:   PetscReal      nzshift;

339:   if (shift >= PETSC_MAX_REAL) { /* Right-open interval */
340:     if (inertia) *inertia = eps->n;
341:   } else if (shift <= PETSC_MIN_REAL) {
342:     if (inertia) *inertia = 0;
343:     if (zeros) *zeros = 0;
344:   } else {
345:     /* If the shift is zero, perturb it to a very small positive value.
346:        The goal is that the nonzero pattern is the same in all cases and reuse
347:        the symbolic factorizations */
348:     nzshift = (shift==0.0)? 10.0/PETSC_MAX_REAL: shift;
349:     STSetShift(eps->st,nzshift);
350:     STSetUp(eps->st);
351:     STGetKSP(eps->st,&ksp);
352:     KSPGetPC(ksp,&pc);
353:     PCFactorGetMatrix(pc,&F);
354:     MatGetInertia(F,inertia,zeros,NULL);
355:   }
356:   return(0);
357: }

361: PetscErrorCode EPSSetUp_KrylovSchur_Slice(EPS eps)
362: {
363:   PetscErrorCode  ierr;
364:   PetscBool       issinv;
365:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data,*ctx_glob;
366:   EPS_SR          sr,sr_loc,sr_glob;
367:   PetscInt        nEigs,dssz=1,i,zeros=0,off=0;
368:   PetscMPIInt     nproc,rank,aux;
369:   MPI_Request     req;
370:   Mat             A,B=NULL;

373:   if (ctx->global) {
374:     if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Must define a computational interval when using EPS_ALL");
375:     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");
376:     if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
377:     if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs cannot be used with spectrum slicing");
378:     if (!((PetscObject)(eps->st))->type_name) { /* default to shift-and-invert */
379:       STSetType(eps->st,STSINVERT);
380:     }
381:     PetscObjectTypeCompareAny((PetscObject)eps->st,&issinv,STSINVERT,STCAYLEY,"");
382:     if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
383:     if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL*1e-2;  /* use tighter tolerance */
384:     if (!eps->max_it) eps->max_it = 100;
385:     if (ctx->nev==1) ctx->nev = PetscMin(40,eps->n);  /* nev not set, use default value */
386:     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");
387:   }
388:   eps->ops->backtransform = NULL;

390:   /* create spectrum slicing context and initialize it */
391:   EPSSliceResetSR(eps);
392:   PetscNewLog(eps,&sr);
393:   ctx->sr = sr;
394:   sr->itsKs = 0;
395:   sr->nleap = 0;
396:   sr->nMAXCompl = eps->nev/4;
397:   sr->iterCompl = eps->max_it/4;
398:   sr->sPres = NULL;
399:   sr->nS = 0;

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

472:     /* compute inertia0 */
473:     EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,ctx->detect?&zeros:NULL);
474:     if (zeros) { /* error in factorization */
475:       if (ctx->npart==1 || ctx_glob->subintset || ((sr->dir>0 && ctx->subc->color==0) || (sr->dir<0 && ctx->subc->color==ctx->npart-1))) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
476:       else { /* perturb shift */
477:         sr->int0 *= (1.0+SLICE_PTOL);
478:         EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,&zeros);
479:         if (zeros) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",sr->int1);
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->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:       }
495:       if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)||(sr->dir<0 && ctx->subc->color>0)) {
496:         MPI_Bcast(&sr->inertia1,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
497:         MPI_Bcast(&sr->int1,1,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
498:       } else sr_glob->inertia1 = sr->inertia1;
499:     }

501:     /* last process in eps comm computes inertia1 */
502:     if (ctx->npart==1 || ((sr->dir>0 && ctx->subc->color==ctx->npart-1) || (sr->dir<0 && ctx->subc->color==0))) {
503:       EPSSliceGetInertia(eps,sr->int1,&sr->inertia1,ctx->detect?&zeros:NULL);
504:       if (zeros) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
505:     }

507:     /* number of eigenvalues in interval */
508:     sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
509:     if (ctx->npart>1) {
510:       /* memory allocate for subinterval eigenpairs */
511:       EPSSliceAllocateSolution(eps,1);
512:     }
513:     dssz = eps->ncv+1;
514:   }
515:   DSSetType(eps->ds,DSHEP);
516:   DSSetCompact(eps->ds,PETSC_TRUE);
517:   DSAllocate(eps->ds,dssz);
518:   /* keep state of subcomm matrices to check that the user does not modify them */
519:   EPSGetOperators(eps,&A,&B);
520:   PetscObjectStateGet((PetscObject)A,&ctx->Astate);
521:   PetscObjectGetId((PetscObject)A,&ctx->Aid);
522:   if (B) { 
523:     PetscObjectStateGet((PetscObject)B,&ctx->Bstate);
524:     PetscObjectGetId((PetscObject)B,&ctx->Bid);
525:   } else {
526:     ctx->Bstate=0;
527:     ctx->Bid=0;
528:   }
529:   return(0);
530: }

534: static PetscErrorCode EPSSliceGatherEigenVectors(EPS eps)
535: {
536:   PetscErrorCode  ierr;
537:   Vec             v,vg,v_loc;
538:   IS              is1,is2;
539:   VecScatter      vec_sc;
540:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
541:   PetscInt        nloc,m0,n0,i,si,idx,*idx1,*idx2,j;
542:   PetscScalar     *array;
543:   EPS_SR          sr_loc;
544:   BV              V_loc;

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

550:   /* Gather parallel eigenvectors */
551:   BVGetColumn(eps->V,0,&v);
552:   VecGetOwnershipRange(v,&n0,&m0);
553:   BVRestoreColumn(eps->V,0,&v);
554:   BVGetColumn(ctx->eps->V,0,&v);
555:   VecGetLocalSize(v,&nloc);
556:   BVRestoreColumn(ctx->eps->V,0,&v);
557:   PetscMalloc2(m0-n0,&idx1,m0-n0,&idx2);
558:   VecCreateMPI(PetscObjectComm((PetscObject)eps),nloc,PETSC_DECIDE,&vg);
559:   idx = -1;
560:   for (si=0;si<ctx->npart;si++) {
561:     j = 0;
562:     for (i=n0;i<m0;i++) {
563:       idx1[j]   = i;
564:       idx2[j++] = i+eps->n*si;
565:     }
566:     ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx1,PETSC_COPY_VALUES,&is1);
567:     ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx2,PETSC_COPY_VALUES,&is2);
568:     BVGetColumn(eps->V,0,&v);
569:     VecScatterCreate(v,is1,vg,is2,&vec_sc);
570:     BVRestoreColumn(eps->V,0,&v);
571:     ISDestroy(&is1);
572:     ISDestroy(&is2);
573:     for (i=0;i<ctx->nconv_loc[si];i++) {
574:       BVGetColumn(eps->V,++idx,&v);
575:       if (ctx->subc->color==si) {
576:         BVGetColumn(V_loc,i,&v_loc);
577:         VecGetArray(v_loc,&array);
578:         VecPlaceArray(vg,array);
579:       }
580:       VecScatterBegin(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE);
581:       VecScatterEnd(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE);
582:       if (ctx->subc->color==si) {
583:         VecResetArray(vg);
584:         VecRestoreArray(v_loc,&array);
585:         BVRestoreColumn(V_loc,i,&v_loc);
586:       }
587:       BVRestoreColumn(eps->V,idx,&v);
588:     }
589:     VecScatterDestroy(&vec_sc);
590:   }
591:   PetscFree2(idx1,idx2);
592:   VecDestroy(&vg);
593:   return(0);
594: }

598: /*
599:   EPSComputeVectors_Slice - Recover Eigenvectors from subcomunicators
600:  */
601: PetscErrorCode EPSComputeVectors_Slice(EPS eps)
602: {
603:   PetscErrorCode  ierr;
604:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

607:   if (ctx->global && ctx->npart>1) {
608:     EPSComputeVectors(ctx->eps);
609:     EPSSliceGatherEigenVectors(eps);
610:   }
611:   return(0);
612: }

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

618: static PetscErrorCode EPSSliceGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
619: {
620:   PetscErrorCode  ierr;
621:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
622:   PetscInt        i=0,j,tmpi;
623:   PetscReal       v,tmpr;
624:   EPS_shift       s;

627:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
628:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
629:   if (!ctx->sr->s0) {  /* EPSSolve not called yet */
630:     *n = 2;
631:   } else {
632:     *n = 1;
633:     s = ctx->sr->s0;
634:     while (s) {
635:       (*n)++;
636:       s = s->neighb[1];
637:     }
638:   }
639:   PetscMalloc1(*n,shifts);
640:   PetscMalloc1(*n,inertias);
641:   if (!ctx->sr->s0) {  /* EPSSolve not called yet */
642:     (*shifts)[0]   = ctx->sr->int0;
643:     (*shifts)[1]   = ctx->sr->int1;
644:     (*inertias)[0] = ctx->sr->inertia0;
645:     (*inertias)[1] = ctx->sr->inertia1;
646:   } else {
647:     s = ctx->sr->s0;
648:     while (s) {
649:       (*shifts)[i]     = s->value;
650:       (*inertias)[i++] = s->inertia;
651:       s = s->neighb[1];
652:     }
653:     (*shifts)[i]   = ctx->sr->int1;
654:     (*inertias)[i] = ctx->sr->inertia1;
655:   }
656:   /* remove possible duplicate in last position */
657:   if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
658:   /* sort result */
659:   for (i=0;i<*n;i++) {
660:     v = (*shifts)[i];
661:     for (j=i+1;j<*n;j++) {
662:       if (v > (*shifts)[j]) {
663:         SWAP((*shifts)[i],(*shifts)[j],tmpr);
664:         SWAP((*inertias)[i],(*inertias)[j],tmpi);
665:         v = (*shifts)[i];
666:       }
667:     }
668:   }
669:   return(0);
670: }

674: static PetscErrorCode EPSSliceGatherSolution(EPS eps)
675: {
676:   PetscErrorCode  ierr;
677:   PetscMPIInt     rank,nproc;
678:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
679:   PetscInt        i,idx,j;
680:   PetscInt        *perm_loc,off=0,*inertias_loc,ns;
681:   PetscScalar     *eigr_loc;
682:   EPS_SR          sr_loc;
683:   PetscReal       *shifts_loc;
684:   PetscMPIInt     *disp,*ns_loc,aux;

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

691:   /* Gather the shifts used and the inertias computed */
692:   EPSSliceGetInertias(ctx->eps,&ns,&shifts_loc,&inertias_loc);
693:   if (ctx->sr->dir>0 && shifts_loc[ns-1]==sr_loc->int1 && ctx->subc->color<ctx->npart-1) ns--;
694:   if (ctx->sr->dir<0 && shifts_loc[ns-1]==sr_loc->int0 && ctx->subc->color>0) {
695:     ns--;
696:     for (i=0;i<ns;i++) {
697:       inertias_loc[i] = inertias_loc[i+1];
698:       shifts_loc[i] = shifts_loc[i+1];
699:     }
700:   }
701:   PetscMalloc1(ctx->npart,&ns_loc);
702:   MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
703:   PetscMPIIntCast(ns,&aux);
704:   if (rank==0) { MPI_Allgather(&aux,1,MPI_INT,ns_loc,1,MPI_INT,ctx->commrank); }
705:   PetscMPIIntCast(ctx->npart,&aux);
706:   MPI_Bcast(ns_loc,aux,MPI_INT,0,PetscSubcommChild(ctx->subc));
707:   ctx->nshifts = 0;
708:   for (i=0;i<ctx->npart;i++) ctx->nshifts += ns_loc[i];
709:   PetscFree(ctx->inertias);
710:   PetscFree(ctx->shifts);
711:   PetscMalloc1(ctx->nshifts,&ctx->inertias);
712:   PetscMalloc1(ctx->nshifts,&ctx->shifts);

714:   /* Gather eigenvalues (same ranks have fully set of eigenvalues)*/
715:   eigr_loc = sr_loc->eigr;
716:   perm_loc = sr_loc->perm;
717:   MPI_Comm_size(((PetscObject)eps)->comm,&nproc);
718:   PetscMalloc1(ctx->npart,&disp);
719:   disp[0] = 0;
720:   for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ctx->nconv_loc[i-1];
721:   if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
722:     PetscMPIIntCast(sr_loc->numEigs,&aux);
723:     MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank); /* eigenvalues */
724:     MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank); /* perm */
725:     for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
726:     PetscMPIIntCast(ns,&aux);
727:     MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank); /* shifts */
728:     MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank); /* inertias */
729:     MPI_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank);
730:   } else { /* subcommunicators with different size */
731:     MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
732:     if (rank==0) {
733:       PetscMPIIntCast(sr_loc->numEigs,&aux);
734:       MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank); /* eigenvalues */
735:       MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank); /* perm */
736:       for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
737:       PetscMPIIntCast(ns,&aux);
738:       MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank); /* shifts */
739:       MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank); /* inertias */
740:       MPI_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank);
741:     }
742:     PetscMPIIntCast(eps->nconv,&aux);
743:     MPI_Bcast(eps->eigr,aux,MPIU_SCALAR,0,PetscSubcommChild(ctx->subc));
744:     MPI_Bcast(eps->perm,aux,MPIU_INT,0,PetscSubcommChild(ctx->subc));
745:     MPI_Bcast(ctx->shifts,ctx->nshifts,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
746:     PetscMPIIntCast(ctx->nshifts,&aux);
747:     MPI_Bcast(ctx->inertias,aux,MPIU_INT,0,PetscSubcommChild(ctx->subc));
748:     MPI_Bcast(&eps->its,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
749:   }
750:   /* Update global array eps->perm */
751:   idx = ctx->nconv_loc[0];
752:   for (i=1;i<ctx->npart;i++) {
753:     off += ctx->nconv_loc[i-1];
754:     for (j=0;j<ctx->nconv_loc[i];j++) eps->perm[idx++] += off;
755:   }

757:   /* Gather parallel eigenvectors */
758:   PetscFree(ns_loc);
759:   PetscFree(disp);
760:   PetscFree(shifts_loc);
761:   PetscFree(inertias_loc);
762:   return(0);
763: }

765: /*
766:    Fills the fields of a shift structure
767: */
770: static PetscErrorCode EPSCreateShift(EPS eps,PetscReal val,EPS_shift neighb0,EPS_shift neighb1)
771: {
772:   PetscErrorCode  ierr;
773:   EPS_shift       s,*pending2;
774:   PetscInt        i;
775:   EPS_SR          sr;
776:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

779:   sr = ctx->sr;
780:   PetscNewLog(eps,&s);
781:   s->value = val;
782:   s->neighb[0] = neighb0;
783:   if (neighb0) neighb0->neighb[1] = s;
784:   s->neighb[1] = neighb1;
785:   if (neighb1) neighb1->neighb[0] = s;
786:   s->comp[0] = PETSC_FALSE;
787:   s->comp[1] = PETSC_FALSE;
788:   s->index = -1;
789:   s->neigs = 0;
790:   s->nconv[0] = s->nconv[1] = 0;
791:   s->nsch[0] = s->nsch[1]=0;
792:   /* Inserts in the stack of pending shifts */
793:   /* If needed, the array is resized */
794:   if (sr->nPend >= sr->maxPend) {
795:     sr->maxPend *= 2;
796:     PetscMalloc1(sr->maxPend,&pending2);
797:     PetscLogObjectMemory((PetscObject)eps,sizeof(EPS_shift));
798:     for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
799:     PetscFree(sr->pending);
800:     sr->pending = pending2;
801:   }
802:   sr->pending[sr->nPend++]=s;
803:   return(0);
804: }

806: /* Prepare for Rational Krylov update */
809: static PetscErrorCode EPSPrepareRational(EPS eps)
810: {
811:   EPS_KRYLOVSCHUR  *ctx=(EPS_KRYLOVSCHUR*)eps->data;
812:   PetscErrorCode   ierr;
813:   PetscInt         dir,i,k,ld,nv;
814:   PetscScalar      *A;
815:   EPS_SR           sr = ctx->sr;
816:   Vec              v;

819:   DSGetLeadingDimension(eps->ds,&ld);
820:   dir = (sr->sPres->neighb[0] == sr->sPrev)?1:-1;
821:   dir*=sr->dir;
822:   k = 0;
823:   for (i=0;i<sr->nS;i++) {
824:     if (dir*PetscRealPart(sr->S[i])>0.0) {
825:       sr->S[k] = sr->S[i];
826:       sr->S[sr->nS+k] = sr->S[sr->nS+i];
827:       BVGetColumn(sr->Vnext,k,&v);
828:       BVCopyVec(eps->V,eps->nconv+i,v);
829:       BVRestoreColumn(sr->Vnext,k,&v);
830:       k++;
831:       if (k>=sr->nS/2)break;
832:     }
833:   }
834:   /* Copy to DS */
835:   DSGetArray(eps->ds,DS_MAT_A,&A);
836:   PetscMemzero(A,ld*ld*sizeof(PetscScalar));
837:   for (i=0;i<k;i++) {
838:     A[i*(1+ld)] = sr->S[i];
839:     A[k+i*ld] = sr->S[sr->nS+i];
840:   }
841:   sr->nS = k;
842:   DSRestoreArray(eps->ds,DS_MAT_A,&A);
843:   DSGetDimensions(eps->ds,&nv,NULL,NULL,NULL,NULL);
844:   DSSetDimensions(eps->ds,nv,0,0,k);
845:   /* Append u to V */
846:   BVGetColumn(sr->Vnext,sr->nS,&v);
847:   BVCopyVec(eps->V,sr->nv,v);
848:   BVRestoreColumn(sr->Vnext,sr->nS,&v);
849:   return(0);
850: }

852: /* Provides next shift to be computed */
855: static PetscErrorCode EPSExtractShift(EPS eps)
856: {
857:   PetscErrorCode   ierr;
858:   PetscInt         iner,zeros=0;
859:   EPS_KRYLOVSCHUR  *ctx=(EPS_KRYLOVSCHUR*)eps->data;
860:   EPS_SR           sr;
861:   PetscReal        newShift;
862:   EPS_shift        sPres;

865:   sr = ctx->sr;
866:   if (sr->nPend > 0) {
867:     sr->sPrev = sr->sPres;
868:     sr->sPres = sr->pending[--sr->nPend];
869:     sPres = sr->sPres;
870:     EPSSliceGetInertia(eps,sPres->value,&iner,ctx->detect?&zeros:NULL);
871:     if (zeros) {
872:       newShift = sPres->value*(1.0+SLICE_PTOL);
873:       if (sr->dir*(sPres->neighb[0] && newShift-sPres->neighb[0]->value) < 0) newShift = (sPres->value+sPres->neighb[0]->value)/2;
874:       else if (sPres->neighb[1] && sr->dir*(sPres->neighb[1]->value-newShift) < 0) newShift = (sPres->value+sPres->neighb[1]->value)/2;
875:       EPSSliceGetInertia(eps,newShift,&iner,&zeros);
876:       if (zeros) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",newShift);
877:       sPres->value = newShift;
878:     }
879:     sr->sPres->inertia = iner;
880:     eps->target = sr->sPres->value;
881:     eps->reason = EPS_CONVERGED_ITERATING;
882:     eps->its = 0;
883:   } else sr->sPres = NULL;
884:   return(0);
885: }

887: /*
888:    Symmetric KrylovSchur adapted to spectrum slicing:
889:    Allows searching an specific amount of eigenvalues in the subintervals left and right.
890:    Returns whether the search has succeeded
891: */
894: static PetscErrorCode EPSKrylovSchur_Slice(EPS eps)
895: {
896:   PetscErrorCode  ierr;
897:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
898:   PetscInt        i,conv,k,l,ld,nv,*iwork,j,p;
899:   Mat             U;
900:   PetscScalar     *Q,*A,rtmp;
901:   PetscReal       *a,*b,beta;
902:   PetscBool       breakdown;
903:   PetscInt        count0,count1;
904:   PetscReal       lambda;
905:   EPS_shift       sPres;
906:   PetscBool       complIterating;
907:   PetscBool       sch0,sch1;
908:   PetscInt        iterCompl=0,n0,n1;
909:   EPS_SR          sr = ctx->sr;

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

952:     /* Solve projected problem and compute residual norm estimates */
953:     if (eps->its == 1 && l > 0) {/* After rational update */
954:       DSGetArray(eps->ds,DS_MAT_A,&A);
955:       DSGetArrayReal(eps->ds,DS_MAT_T,&a);
956:       b = a + ld;
957:       k = eps->nconv+l;
958:       A[k*ld+k-1] = A[(k-1)*ld+k];
959:       A[k*ld+k] = a[k];
960:       for (j=k+1; j< nv; j++) {
961:         A[j*ld+j] = a[j];
962:         A[j*ld+j-1] = b[j-1] ;
963:         A[(j-1)*ld+j] = b[j-1];
964:       }
965:       DSRestoreArray(eps->ds,DS_MAT_A,&A);
966:       DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
967:       DSSolve(eps->ds,eps->eigr,NULL);
968:       DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
969:       DSSetCompact(eps->ds,PETSC_TRUE);
970:     } else { /* Restart */
971:       DSSolve(eps->ds,eps->eigr,NULL);
972:       DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
973:     }
974:     /* Residual */
975:     EPSKrylovConvergence(eps,PETSC_TRUE,eps->nconv,nv-eps->nconv,beta,1.0,&k);

977:     if (ctx->lock) {
978:       /* Check convergence */
979:       DSGetArrayReal(eps->ds,DS_MAT_T,&a);
980:       b = a + ld;
981:       conv = 0;
982:       j = k = eps->nconv;
983:       for (i=eps->nconv;i<nv;i++) if (eps->errest[i] < eps->tol) conv++;
984:       for (i=eps->nconv;i<nv;i++) {
985:         if (eps->errest[i] < eps->tol) {
986:           iwork[j++]=i;
987:         } else iwork[conv+k++]=i;
988:       }
989:       for (i=eps->nconv;i<nv;i++) {
990:         a[i]=PetscRealPart(eps->eigr[i]);
991:         b[i]=eps->errest[i];
992:       }
993:       for (i=eps->nconv;i<nv;i++) {
994:         eps->eigr[i] = a[iwork[i]];
995:         eps->errest[i] = b[iwork[i]];
996:       }
997:       for (i=eps->nconv;i<nv;i++) {
998:         a[i]=PetscRealPart(eps->eigr[i]);
999:         b[i]=eps->errest[i];
1000:       }
1001:       DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
1002:       DSGetArray(eps->ds,DS_MAT_Q,&Q);
1003:       for (i=eps->nconv;i<nv;i++) {
1004:         p=iwork[i];
1005:         if (p!=i) {
1006:           j=i+1;
1007:           while (iwork[j]!=i) j++;
1008:           iwork[j]=p;iwork[i]=i;
1009:           for (k=0;k<nv;k++) {
1010:             rtmp=Q[k+p*ld];Q[k+p*ld]=Q[k+i*ld];Q[k+i*ld]=rtmp;
1011:           }
1012:         }
1013:       }
1014:       DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
1015:       k=eps->nconv+conv;
1016:     }

1018:     /* Checking values obtained for completing */
1019:     for (i=0;i<k;i++) {
1020:       sr->back[i]=eps->eigr[i];
1021:     }
1022:     STBackTransform(eps->st,k,sr->back,eps->eigi);
1023:     count0=count1=0;
1024:     for (i=0;i<k;i++) {
1025:       lambda = PetscRealPart(sr->back[i]);
1026:       if (((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0)) count0++;
1027:       if (((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0)) count1++;
1028:     }
1029:     if (k>eps->nev && eps->ncv-k<5) eps->reason = EPS_CONVERGED_TOL;
1030:     else {
1031:       /* Checks completion */
1032:       if ((!sch0||count0 >= sPres->nsch[0]) && (!sch1 ||count1 >= sPres->nsch[1])) {
1033:         eps->reason = EPS_CONVERGED_TOL;
1034:       } else {
1035:         if (!complIterating && eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
1036:         if (complIterating) {
1037:           if (--iterCompl <= 0) eps->reason = EPS_DIVERGED_ITS;
1038:         } else if (k >= eps->nev) {
1039:           n0 = sPres->nsch[0]-count0;
1040:           n1 = sPres->nsch[1]-count1;
1041:           if (sr->iterCompl>0 && ((n0>0 && n0<= sr->nMAXCompl)||(n1>0&&n1<=sr->nMAXCompl))) {
1042:             /* Iterating for completion*/
1043:             complIterating = PETSC_TRUE;
1044:             if (n0 >sr->nMAXCompl)sch0 = PETSC_FALSE;
1045:             if (n1 >sr->nMAXCompl)sch1 = PETSC_FALSE;
1046:             iterCompl = sr->iterCompl;
1047:           } else eps->reason = EPS_CONVERGED_TOL;
1048:         }
1049:       }
1050:     }
1051:     /* Update l */
1052:     if (eps->reason == EPS_CONVERGED_ITERATING) l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
1053:     else l = 0;
1054:     if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
1055:     if (breakdown) l=0;

1057:     if (eps->reason == EPS_CONVERGED_ITERATING) {
1058:       if (breakdown) {
1059:         /* Start a new Lanczos factorization */
1060:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
1061:         EPSGetStartVector(eps,k,&breakdown);
1062:         if (breakdown) {
1063:           eps->reason = EPS_DIVERGED_BREAKDOWN;
1064:           PetscInfo(eps,"Unable to generate more start vectors\n");
1065:         }
1066:       } else {
1067:         /* Prepare the Rayleigh quotient for restart */
1068:         DSGetArrayReal(eps->ds,DS_MAT_T,&a);
1069:         DSGetArray(eps->ds,DS_MAT_Q,&Q);
1070:         b = a + ld;
1071:         for (i=k;i<k+l;i++) {
1072:           a[i] = PetscRealPart(eps->eigr[i]);
1073:           b[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
1074:         }
1075:         DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
1076:         DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
1077:       }
1078:     }
1079:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1080:     DSGetMat(eps->ds,DS_MAT_Q,&U);
1081:     BVMultInPlace(eps->V,U,eps->nconv,k+l);
1082:     MatDestroy(&U);

1084:     /* Normalize u and append it to V */
1085:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
1086:       BVCopyColumn(eps->V,nv,k+l);
1087:     }
1088:     eps->nconv = k;
1089:     if (eps->reason != EPS_CONVERGED_ITERATING) {
1090:       /* Store approximated values for next shift */
1091:       DSGetArray(eps->ds,DS_MAT_Q,&Q);
1092:       sr->nS = l;
1093:       for (i=0;i<l;i++) {
1094:         sr->S[i] = eps->eigr[i+k];/* Diagonal elements */
1095:         sr->S[i+l] = Q[nv-1+(i+k)*ld]*beta; /* Out of diagonal elements */
1096:       }
1097:       DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
1098:     }
1099:   }
1100:   /* Check for completion */
1101:   for (i=0;i< eps->nconv; i++) {
1102:     if ((sr->dir)*PetscRealPart(eps->eigr[i])>0) sPres->nconv[1]++;
1103:     else sPres->nconv[0]++;
1104:   }
1105:   sPres->comp[0] = PetscNot(count0 < sPres->nsch[0]);
1106:   sPres->comp[1] = PetscNot(count1 < sPres->nsch[1]);
1107:   if (count0 > sPres->nsch[0] || count1 > sPres->nsch[1])SETERRQ(PetscObjectComm((PetscObject)eps),1,"Mismatch between number of values found and information from inertia, consider using EPSKrylovSchurSetDetectZeros()");
1108:   PetscFree(iwork);
1109:   return(0);
1110: }

1112: /*
1113:   Obtains value of subsequent shift
1114: */
1117: static PetscErrorCode EPSGetNewShiftValue(EPS eps,PetscInt side,PetscReal *newS)
1118: {
1119:   PetscReal       lambda,d_prev;
1120:   PetscInt        i,idxP;
1121:   EPS_SR          sr;
1122:   EPS_shift       sPres,s;
1123:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

1126:   sr = ctx->sr;
1127:   sPres = sr->sPres;
1128:   if (sPres->neighb[side]) {
1129:   /* Completing a previous interval */
1130:     if (!sPres->neighb[side]->neighb[side] && sPres->neighb[side]->nconv[side]==0) { /* One of the ends might be too far from eigenvalues */
1131:       if (side) *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[sr->indexEig-1]]))/2;
1132:       else *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[0]]))/2;
1133:     } else *newS=(sPres->value + sPres->neighb[side]->value)/2;
1134:   } else { /* (Only for side=1). Creating a new interval. */
1135:     if (sPres->neigs==0) {/* No value has been accepted*/
1136:       if (sPres->neighb[0]) {
1137:         /* Multiplying by 10 the previous distance */
1138:         *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
1139:         sr->nleap++;
1140:         /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
1141:         if (!sr->hasEnd && sr->nleap > 5) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to compute the wanted eigenvalues with open interval");
1142:       } else { /* First shift */
1143:         if (eps->nconv != 0) {
1144:           /* Unaccepted values give information for next shift */
1145:           idxP=0;/* Number of values left from shift */
1146:           for (i=0;i<eps->nconv;i++) {
1147:             lambda = PetscRealPart(sr->eigr[i]);
1148:             if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
1149:             else break;
1150:           }
1151:           /* Avoiding subtraction of eigenvalues (might be the same).*/
1152:           if (idxP>0) {
1153:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(sr->eigr[0]))/(idxP+0.3);
1154:           } else {
1155:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(sr->eigr[eps->nconv-1]))/(eps->nconv+0.3);
1156:           }
1157:           *newS = sPres->value + ((sr->dir)*d_prev*eps->nev)/2;
1158:         } else { /* No values found, no information for next shift */
1159:           SETERRQ(PetscObjectComm((PetscObject)eps),1,"First shift renders no information");
1160:         }
1161:       }
1162:     } else { /* Accepted values found */
1163:       sr->nleap = 0;
1164:       /* Average distance of values in previous subinterval */
1165:       s = sPres->neighb[0];
1166:       while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
1167:         s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
1168:       }
1169:       if (s) {
1170:         d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
1171:       } else { /* First shift. Average distance obtained with values in this shift */
1172:         /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
1173:         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)) {
1174:           d_prev =  PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0])))/(sPres->neigs+0.3);
1175:         } else {
1176:           d_prev = PetscAbsReal(PetscRealPart(sr->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
1177:         }
1178:       }
1179:       /* Average distance is used for next shift by adding it to value on the right or to shift */
1180:       if ((sr->dir)*(PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
1181:         *newS = PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(eps->nev))/2;
1182:       } else { /* Last accepted value is on the left of shift. Adding to shift */
1183:         *newS = sPres->value + ((sr->dir)*d_prev*(eps->nev))/2;
1184:       }
1185:     }
1186:     /* End of interval can not be surpassed */
1187:     if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
1188:   }/* of neighb[side]==null */
1189:   return(0);
1190: }

1192: /*
1193:   Function for sorting an array of real values
1194: */
1197: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
1198: {
1199:   PetscReal      re;
1200:   PetscInt       i,j,tmp;

1203:   if (!prev) for (i=0;i<nr;i++) perm[i] = i;
1204:   /* Insertion sort */
1205:   for (i=1;i<nr;i++) {
1206:     re = PetscRealPart(r[perm[i]]);
1207:     j = i-1;
1208:     while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
1209:       tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
1210:     }
1211:   }
1212:   return(0);
1213: }

1215: /* Stores the pairs obtained since the last shift in the global arrays */
1218: static PetscErrorCode EPSStoreEigenpairs(EPS eps)
1219: {
1220:   PetscErrorCode  ierr;
1221:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1222:   PetscReal       lambda,err,norm;
1223:   PetscInt        i,count;
1224:   PetscBool       iscayley;
1225:   EPS_SR          sr = ctx->sr;
1226:   EPS_shift       sPres;
1227:   Vec             v,w;

1230:   sPres = sr->sPres;
1231:   sPres->index = sr->indexEig;
1232:   count = sr->indexEig;
1233:   /* Back-transform */
1234:   STBackTransform(eps->st,eps->nconv,eps->eigr,eps->eigi);
1235:   PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
1236:   /* Sort eigenvalues */
1237:   sortRealEigenvalues(eps->eigr,eps->perm,eps->nconv,PETSC_FALSE,sr->dir);
1238:   /* Values stored in global array */
1239:   for (i=0;i<eps->nconv;i++) {
1240:     lambda = PetscRealPart(eps->eigr[eps->perm[i]]);
1241:     err = eps->errest[eps->perm[i]];

1243:     if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
1244:       if (count>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error in Spectrum Slicing");
1245:       sr->eigr[count] = lambda;
1246:       sr->errest[count] = err;
1247:       /* Explicit purification */
1248:       if (eps->purify) {
1249:         BVGetColumn(sr->V,count,&v);
1250:         BVGetColumn(eps->V,eps->perm[i],&w);
1251:         STApply(eps->st,w,v);
1252:         BVRestoreColumn(sr->V,count,&v);
1253:         BVRestoreColumn(eps->V,eps->perm[i],&w);
1254:         BVNormColumn(sr->V,count,NORM_2,&norm);
1255:         BVScaleColumn(sr->V,count,1.0/norm);
1256:       } else {
1257:         BVGetColumn(eps->V,eps->perm[i],&w);
1258:         BVInsertVec(sr->V,count,w);
1259:         BVRestoreColumn(eps->V,eps->perm[i],&w);
1260:         BVNormColumn(sr->V,count,NORM_2,&norm);
1261:         BVScaleColumn(sr->V,count,1.0/norm);
1262:       }
1263:       count++;
1264:     }
1265:   }
1266:   sPres->neigs = count - sr->indexEig;
1267:   sr->indexEig = count;
1268:   /* Global ordering array updating */
1269:   sortRealEigenvalues(sr->eigr,sr->perm,count,PETSC_TRUE,sr->dir);
1270:   return(0);
1271: }

1275: static PetscErrorCode EPSLookForDeflation(EPS eps)
1276: {
1277:   PetscErrorCode  ierr;
1278:   PetscReal       val;
1279:   PetscInt        i,count0=0,count1=0;
1280:   EPS_shift       sPres;
1281:   PetscInt        ini,fin,k,idx0,idx1;
1282:   EPS_SR          sr;
1283:   Vec             v;
1284:   EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;

1287:   sr = ctx->sr;
1288:   sPres = sr->sPres;

1290:   if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
1291:   else ini = 0;
1292:   fin = sr->indexEig;
1293:   /* Selection of ends for searching new values */
1294:   if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
1295:   else sPres->ext[0] = sPres->neighb[0]->value;
1296:   if (!sPres->neighb[1]) {
1297:     if (sr->hasEnd) sPres->ext[1] = sr->int1;
1298:     else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
1299:   } else sPres->ext[1] = sPres->neighb[1]->value;
1300:   /* Selection of values between right and left ends */
1301:   for (i=ini;i<fin;i++) {
1302:     val=PetscRealPart(sr->eigr[sr->perm[i]]);
1303:     /* Values to the right of left shift */
1304:     if ((sr->dir)*(val - sPres->ext[1]) < 0) {
1305:       if ((sr->dir)*(val - sPres->value) < 0) count0++;
1306:       else count1++;
1307:     } else break;
1308:   }
1309:   /* The number of values on each side are found */
1310:   if (sPres->neighb[0]) {
1311:     sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
1312:     if (sPres->nsch[0]<0)SETERRQ(PetscObjectComm((PetscObject)eps),1,"Mismatch between number of values found and information from inertia, consider using EPSKrylovSchurSetDetectZeros()");
1313:   } else sPres->nsch[0] = 0;

1315:   if (sPres->neighb[1]) {
1316:     sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
1317:     if (sPres->nsch[1]<0)SETERRQ(PetscObjectComm((PetscObject)eps),1,"Mismatch between number of values found and information from inertia, consider using EPSKrylovSchurSetDetectZeros()");
1318:   } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);

1320:   /* Completing vector of indexes for deflation */
1321:   idx0 = ini;
1322:   idx1 = ini+count0+count1;
1323:   k=0;
1324:   for (i=idx0;i<idx1;i++) sr->idxDef[k++]=sr->perm[i];
1325:   BVDuplicateResize(eps->V,k+eps->ncv+1,&sr->Vnext);
1326:   BVSetNumConstraints(sr->Vnext,k);
1327:   for (i=0;i<k;i++) {
1328:     BVGetColumn(sr->Vnext,-i-1,&v);
1329:     BVCopyVec(sr->V,sr->idxDef[i],v);
1330:     BVRestoreColumn(sr->Vnext,-i-1,&v);
1331:   }

1333:   /* For rational Krylov */
1334:   if (sr->nS>0 && (sr->sPrev == sr->sPres->neighb[0] || sr->sPrev == sr->sPres->neighb[1])) {
1335:     EPSPrepareRational(eps);
1336:   }
1337:   eps->nconv = 0;
1338:   /* Get rid of temporary Vnext */
1339:   BVDestroy(&eps->V);
1340:   eps->V = sr->Vnext;
1341:   sr->Vnext = NULL;
1342:   return(0);
1343: }

1347: PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS eps)
1348: {
1349:   PetscErrorCode   ierr;
1350:   PetscInt         i,lds;
1351:   PetscReal        newS;
1352:   EPS_KRYLOVSCHUR  *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1353:   EPS_SR           sr=ctx->sr;
1354:   Mat              A,B=NULL;
1355:   PetscObjectState Astate,Bstate=0;
1356:   PetscObjectId    Aid,Bid=0;

1359:   PetscCitationsRegister(citation,&cited);
1360:   if (ctx->global) {
1361:     EPSSolve_KrylovSchur_Slice(ctx->eps);
1362:     ctx->eps->state = EPS_STATE_SOLVED;
1363:     eps->reason = EPS_CONVERGED_TOL;
1364:     if (ctx->npart>1) {
1365:       /* Gather solution from subsolvers */
1366:       EPSSliceGatherSolution(eps);
1367:     } else {
1368:       eps->nconv = sr->numEigs;
1369:       eps->its   = ctx->eps->its;
1370:       PetscFree(ctx->inertias);
1371:       PetscFree(ctx->shifts);
1372:       EPSSliceGetInertias(ctx->eps,&ctx->nshifts,&ctx->shifts,&ctx->inertias);
1373:     }
1374:   } else {
1375:     if (ctx->npart==1) {
1376:       sr->eigr   = ctx->eps->eigr;
1377:       sr->eigi   = ctx->eps->eigi;
1378:       sr->perm   = ctx->eps->perm;
1379:       sr->errest = ctx->eps->errest;
1380:       sr->V      = ctx->eps->V;
1381:     }
1382:     /* Check that the user did not modify subcomm matrices */
1383:     EPSGetOperators(eps,&A,&B);
1384:     PetscObjectStateGet((PetscObject)A,&Astate);
1385:     PetscObjectGetId((PetscObject)A,&Aid);
1386:     if (B) { 
1387:       PetscObjectStateGet((PetscObject)B,&Bstate);
1388:       PetscObjectGetId((PetscObject)B,&Bid);
1389:     }
1390:     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");
1391:     /* Only with eigenvalues present in the interval ...*/
1392:     if (sr->numEigs==0) {
1393:       eps->reason = EPS_CONVERGED_TOL;
1394:       return(0);
1395:     }
1396:     /* Array of pending shifts */
1397:     sr->maxPend = 100; /* Initial size */
1398:     sr->nPend = 0;
1399:     PetscMalloc1(sr->maxPend,&sr->pending);
1400:     PetscLogObjectMemory((PetscObject)eps,(sr->maxPend)*sizeof(EPS_shift));
1401:     EPSCreateShift(eps,sr->int0,NULL,NULL);
1402:     /* extract first shift */
1403:     sr->sPrev = NULL;
1404:     sr->sPres = sr->pending[--sr->nPend];
1405:     sr->sPres->inertia = sr->inertia0;
1406:     eps->target = sr->sPres->value;
1407:     sr->s0 = sr->sPres;
1408:     sr->indexEig = 0;
1409:     /* Memory reservation for auxiliary variables */
1410:     lds = PetscMin(eps->mpd,eps->ncv);
1411:     PetscCalloc1(lds*lds,&sr->S);
1412:     PetscMalloc1(eps->ncv,&sr->back);
1413:     PetscLogObjectMemory((PetscObject)eps,(sr->numEigs+2*eps->ncv)*sizeof(PetscScalar));
1414:     for (i=0;i<sr->numEigs;i++) {
1415:       sr->eigr[i]   = 0.0;
1416:       sr->eigi[i]   = 0.0;
1417:       sr->errest[i] = 0.0;
1418:       sr->perm[i]   = i;
1419:     }
1420:     /* Vectors for deflation */
1421:     PetscMalloc1(sr->numEigs,&sr->idxDef);
1422:     PetscLogObjectMemory((PetscObject)eps,sr->numEigs*sizeof(PetscInt));
1423:     sr->indexEig = 0;
1424:     /* Main loop */
1425:     while (sr->sPres) {
1426:       /* Search for deflation */
1427:       EPSLookForDeflation(eps);
1428:       /* KrylovSchur */
1429:       EPSKrylovSchur_Slice(eps);

1431:       EPSStoreEigenpairs(eps);
1432:       /* Select new shift */
1433:       if (!sr->sPres->comp[1]) {
1434:         EPSGetNewShiftValue(eps,1,&newS);
1435:         EPSCreateShift(eps,newS,sr->sPres,sr->sPres->neighb[1]);
1436:       }
1437:       if (!sr->sPres->comp[0]) {
1438:         /* Completing earlier interval */
1439:         EPSGetNewShiftValue(eps,0,&newS);
1440:         EPSCreateShift(eps,newS,sr->sPres->neighb[0],sr->sPres);
1441:       }
1442:       /* Preparing for a new search of values */
1443:       EPSExtractShift(eps);
1444:     }

1446:     /* Updating eps values prior to exit */
1447:     PetscFree(sr->S);
1448:     PetscFree(sr->idxDef);
1449:     PetscFree(sr->pending);
1450:     PetscFree(sr->back);
1451:     BVDuplicateResize(eps->V,eps->ncv+1,&sr->Vnext);
1452:     BVSetNumConstraints(sr->Vnext,0);
1453:     BVDestroy(&eps->V);
1454:     eps->V      = sr->Vnext;
1455:     eps->nconv  = sr->indexEig;
1456:     eps->reason = EPS_CONVERGED_TOL;
1457:     eps->its    = sr->itsKs;
1458:     eps->nds    = 0;
1459:   }
1460:   return(0);
1461: }