Actual source code: ks-slice.c
slepc-3.13.4 2020-09-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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=shift;
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: STGetKSP(eps->st,&ksp);
338: KSPGetPC(ksp,&pc);
339: PCFactorGetMatrix(pc,&F);
340: MatGetInertia(F,inertia,zeros,NULL);
341: }
342: if (inertia) { PetscInfo2(eps,"Computed inertia at shift %g: %D\n",(double)nzshift,*inertia); }
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,hiteig=0;
362: PetscMPIInt nproc,rank=0,aux;
363: PetscReal r;
364: MPI_Request req;
365: Mat A,B=NULL;
366: SlepcSC sc;
367: DSParallelType ptype;
370: if (ctx->global) {
371: 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");
372: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
373: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs cannot be used with spectrum slicing");
374: PetscObjectTypeCompareAny((PetscObject)eps->st,&issinv,STSINVERT,STCAYLEY,"");
375: if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
376: if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL*1e-2; /* use tighter tolerance */
377: if (eps->max_it==PETSC_DEFAULT) eps->max_it = 100;
378: if (ctx->nev==1) ctx->nev = PetscMin(40,eps->n); /* nev not set, use default value */
379: 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");
380: }
381: eps->ops->backtransform = EPSBackTransform_Skip;
383: /* create spectrum slicing context and initialize it */
384: EPSSliceResetSR(eps);
385: PetscNewLog(eps,&sr);
386: ctx->sr = sr;
387: sr->itsKs = 0;
388: sr->nleap = 0;
389: sr->nMAXCompl = eps->nev/4;
390: sr->iterCompl = eps->max_it/4;
391: sr->sPres = NULL;
392: sr->nS = 0;
394: if (ctx->npart==1 || ctx->global) {
395: /* check presence of ends and finding direction */
396: if ((eps->inta > PETSC_MIN_REAL && !(ctx->subintervals && ctx->subintervals[0]==ctx->subintervals[1])) || eps->intb >= PETSC_MAX_REAL) {
397: sr->int0 = eps->inta;
398: sr->int1 = eps->intb;
399: sr->dir = 1;
400: if (eps->intb >= PETSC_MAX_REAL) { /* Right-open interval */
401: sr->hasEnd = PETSC_FALSE;
402: } else sr->hasEnd = PETSC_TRUE;
403: } else {
404: sr->int0 = eps->intb;
405: sr->int1 = eps->inta;
406: sr->dir = -1;
407: sr->hasEnd = PetscNot(eps->inta <= PETSC_MIN_REAL);
408: }
409: }
410: if (ctx->global) {
411: /* prevent computation of factorization in global eps */
412: STSetTransform(eps->st,PETSC_FALSE);
413: EPSSetDimensions_Default(eps,ctx->nev,&ctx->ncv,&ctx->mpd);
414: /* create subintervals and initialize auxiliary eps for slicing runs */
415: EPSSliceGetEPS(eps);
416: sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
417: if (ctx->npart>1) {
418: if ((sr->dir>0&&ctx->subc->color==0)||(sr->dir<0&&ctx->subc->color==ctx->npart-1)) sr->inertia0 = sr_loc->inertia0;
419: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
420: if (!rank) {
421: MPI_Bcast(&sr->inertia0,1,MPIU_INT,(sr->dir>0)?0:ctx->npart-1,ctx->commrank);
422: }
423: MPI_Bcast(&sr->inertia0,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
424: PetscFree(ctx->nconv_loc);
425: PetscMalloc1(ctx->npart,&ctx->nconv_loc);
426: MPI_Comm_size(((PetscObject)eps)->comm,&nproc);
427: if (sr->dir<0) off = 1;
428: if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
429: PetscMPIIntCast(sr_loc->numEigs,&aux);
430: MPI_Allgather(&aux,1,MPI_INT,ctx->nconv_loc,1,MPI_INT,ctx->commrank);
431: MPI_Allgather(sr_loc->dir==sr->dir?&sr_loc->int0:&sr_loc->int1,1,MPIU_REAL,ctx->subintervals+off,1,MPIU_REAL,ctx->commrank);
432: } else {
433: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
434: if (!rank) {
435: PetscMPIIntCast(sr_loc->numEigs,&aux);
436: MPI_Allgather(&aux,1,MPI_INT,ctx->nconv_loc,1,MPI_INT,ctx->commrank);
437: MPI_Allgather(sr_loc->dir==sr->dir?&sr_loc->int0:&sr_loc->int1,1,MPIU_REAL,ctx->subintervals+off,1,MPIU_REAL,ctx->commrank);
438: }
439: PetscMPIIntCast(ctx->npart,&aux);
440: MPI_Bcast(ctx->nconv_loc,aux,MPI_INT,0,PetscSubcommChild(ctx->subc));
441: MPI_Bcast(ctx->subintervals+off,aux,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
442: }
443: nEigs = 0;
444: for (i=0;i<ctx->npart;i++) nEigs += ctx->nconv_loc[i];
445: } else {
446: nEigs = sr_loc->numEigs;
447: sr->inertia0 = sr_loc->inertia0;
448: sr->dir = sr_loc->dir;
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: }
465: STSetUp(eps->st);
467: /* compute inertia0 */
468: EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,ctx->detect?&zeros:NULL);
469: /* undocumented option to control what to do when an eigenvalue is found:
470: - error out if it's the endpoint of the user-provided interval (or sub-interval)
471: - if it's an endpoint computed internally:
472: + if hiteig=0 error out
473: + else if hiteig=1 the subgroup that hit the eigenvalue does nothing
474: + otherwise the subgroup that hit the eigenvalue perturbs the shift and recomputes inertia
475: */
476: PetscOptionsGetInt(NULL,NULL,"-eps_krylovschur_hiteigenvalue",&hiteig,NULL);
477: if (zeros) { /* error in factorization */
478: 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");
479: else if (ctx_glob->subintset && !hiteig) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
480: else {
481: if (hiteig==1) { /* idle subgroup */
482: sr->inertia0 = -1;
483: } else { /* perturb shift */
484: sr->int0 *= (1.0+SLICE_PTOL);
485: EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,&zeros);
486: if (zeros) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",sr->int1);
487: }
488: }
489: }
490: if (ctx->npart>1) {
491: /* inertia1 is received from neighbour */
492: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
493: if (!rank) {
494: if (sr->inertia0!=-1 && ((sr->dir>0 && ctx->subc->color>0) || (sr->dir<0 && ctx->subc->color<ctx->npart-1))) { /* send inertia0 to neighbour0 */
495: MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
496: MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
497: }
498: if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)|| (sr->dir<0 && ctx->subc->color>0)) { /* receive inertia1 from neighbour1 */
499: MPI_Recv(&(sr->inertia1),1,MPIU_INT,ctx->subc->color+sr->dir,0,ctx->commrank,MPI_STATUS_IGNORE);
500: MPI_Recv(&(sr->int1),1,MPIU_REAL,ctx->subc->color+sr->dir,0,ctx->commrank,MPI_STATUS_IGNORE);
501: }
502: if (sr->inertia0==-1 && !(sr->dir>0 && ctx->subc->color==ctx->npart-1) && !(sr->dir<0 && ctx->subc->color==0)) {
503: sr->inertia0 = sr->inertia1; sr->int0 = sr->int1;
504: MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
505: MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
506: }
507: }
508: if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)||(sr->dir<0 && ctx->subc->color>0)) {
509: MPI_Bcast(&sr->inertia1,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
510: MPI_Bcast(&sr->int1,1,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
511: } else sr_glob->inertia1 = sr->inertia1;
512: }
514: /* last process in eps comm computes inertia1 */
515: if (ctx->npart==1 || ((sr->dir>0 && ctx->subc->color==ctx->npart-1) || (sr->dir<0 && ctx->subc->color==0))) {
516: EPSSliceGetInertia(eps,sr->int1,&sr->inertia1,ctx->detect?&zeros:NULL);
517: if (zeros) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
518: if (!rank && sr->inertia0==-1) {
519: sr->inertia0 = sr->inertia1; sr->int0 = sr->int1;
520: MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
521: MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
522: }
523: if (sr->hasEnd) {
524: sr->dir = -sr->dir; r = sr->int0; sr->int0 = sr->int1; sr->int1 = r;
525: i = sr->inertia0; sr->inertia0 = sr->inertia1; sr->inertia1 = i;
526: }
527: }
529: /* number of eigenvalues in interval */
530: sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
531: if (ctx->npart>1) {
532: /* memory allocate for subinterval eigenpairs */
533: EPSSliceAllocateSolution(eps,1);
534: }
535: dssz = eps->ncv+1;
536: if (sr->numEigs>0) {
537: DSGetSlepcSC(eps->ds,&sc);
538: sc->rg = NULL;
539: sc->comparison = SlepcCompareLargestMagnitude;
540: sc->comparisonctx = NULL;
541: sc->map = NULL;
542: sc->mapobj = NULL;
543: }
544: DSGetParallel(ctx->eps->ds,&ptype);
545: DSSetParallel(eps->ds,ptype);
546: DSGetMethod(ctx->eps->ds,&method);
547: DSSetMethod(eps->ds,method);
548: }
549: DSSetType(eps->ds,DSHEP);
550: DSSetCompact(eps->ds,PETSC_TRUE);
551: DSAllocate(eps->ds,dssz);
552: /* keep state of subcomm matrices to check that the user does not modify them */
553: EPSGetOperators(eps,&A,&B);
554: PetscObjectStateGet((PetscObject)A,&ctx->Astate);
555: PetscObjectGetId((PetscObject)A,&ctx->Aid);
556: if (B) {
557: PetscObjectStateGet((PetscObject)B,&ctx->Bstate);
558: PetscObjectGetId((PetscObject)B,&ctx->Bid);
559: } else {
560: ctx->Bstate=0;
561: ctx->Bid=0;
562: }
563: return(0);
564: }
566: static PetscErrorCode EPSSliceGatherEigenVectors(EPS eps)
567: {
568: PetscErrorCode ierr;
569: Vec v,vg,v_loc;
570: IS is1,is2;
571: VecScatter vec_sc;
572: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
573: PetscInt nloc,m0,n0,i,si,idx,*idx1,*idx2,j;
574: PetscScalar *array;
575: EPS_SR sr_loc;
576: BV V_loc;
579: sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
580: V_loc = sr_loc->V;
582: /* Gather parallel eigenvectors */
583: BVGetColumn(eps->V,0,&v);
584: VecGetOwnershipRange(v,&n0,&m0);
585: BVRestoreColumn(eps->V,0,&v);
586: BVGetColumn(ctx->eps->V,0,&v);
587: VecGetLocalSize(v,&nloc);
588: BVRestoreColumn(ctx->eps->V,0,&v);
589: PetscMalloc2(m0-n0,&idx1,m0-n0,&idx2);
590: VecCreateMPI(PetscObjectComm((PetscObject)eps),nloc,PETSC_DECIDE,&vg);
591: idx = -1;
592: for (si=0;si<ctx->npart;si++) {
593: j = 0;
594: for (i=n0;i<m0;i++) {
595: idx1[j] = i;
596: idx2[j++] = i+eps->n*si;
597: }
598: ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx1,PETSC_COPY_VALUES,&is1);
599: ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx2,PETSC_COPY_VALUES,&is2);
600: BVGetColumn(eps->V,0,&v);
601: VecScatterCreate(v,is1,vg,is2,&vec_sc);
602: BVRestoreColumn(eps->V,0,&v);
603: ISDestroy(&is1);
604: ISDestroy(&is2);
605: for (i=0;i<ctx->nconv_loc[si];i++) {
606: BVGetColumn(eps->V,++idx,&v);
607: if (ctx->subc->color==si) {
608: BVGetColumn(V_loc,i,&v_loc);
609: VecGetArray(v_loc,&array);
610: VecPlaceArray(vg,array);
611: }
612: VecScatterBegin(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE);
613: VecScatterEnd(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE);
614: if (ctx->subc->color==si) {
615: VecResetArray(vg);
616: VecRestoreArray(v_loc,&array);
617: BVRestoreColumn(V_loc,i,&v_loc);
618: }
619: BVRestoreColumn(eps->V,idx,&v);
620: }
621: VecScatterDestroy(&vec_sc);
622: }
623: PetscFree2(idx1,idx2);
624: VecDestroy(&vg);
625: return(0);
626: }
628: /*
629: EPSComputeVectors_Slice - Recover Eigenvectors from subcomunicators
630: */
631: PetscErrorCode EPSComputeVectors_Slice(EPS eps)
632: {
633: PetscErrorCode ierr;
634: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
637: if (ctx->global && ctx->npart>1) {
638: EPSComputeVectors(ctx->eps);
639: EPSSliceGatherEigenVectors(eps);
640: }
641: return(0);
642: }
644: #define SWAP(a,b,t) {t=a;a=b;b=t;}
646: static PetscErrorCode EPSSliceGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
647: {
648: PetscErrorCode ierr;
649: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
650: PetscInt i=0,j,tmpi;
651: PetscReal v,tmpr;
652: EPS_shift s;
655: if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
656: if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
657: if (!ctx->sr->s0) { /* EPSSolve not called yet */
658: *n = 2;
659: } else {
660: *n = 1;
661: s = ctx->sr->s0;
662: while (s) {
663: (*n)++;
664: s = s->neighb[1];
665: }
666: }
667: PetscMalloc1(*n,shifts);
668: PetscMalloc1(*n,inertias);
669: if (!ctx->sr->s0) { /* EPSSolve not called yet */
670: (*shifts)[0] = ctx->sr->int0;
671: (*shifts)[1] = ctx->sr->int1;
672: (*inertias)[0] = ctx->sr->inertia0;
673: (*inertias)[1] = ctx->sr->inertia1;
674: } else {
675: s = ctx->sr->s0;
676: while (s) {
677: (*shifts)[i] = s->value;
678: (*inertias)[i++] = s->inertia;
679: s = s->neighb[1];
680: }
681: (*shifts)[i] = ctx->sr->int1;
682: (*inertias)[i] = ctx->sr->inertia1;
683: }
684: /* remove possible duplicate in last position */
685: if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
686: /* sort result */
687: for (i=0;i<*n;i++) {
688: v = (*shifts)[i];
689: for (j=i+1;j<*n;j++) {
690: if (v > (*shifts)[j]) {
691: SWAP((*shifts)[i],(*shifts)[j],tmpr);
692: SWAP((*inertias)[i],(*inertias)[j],tmpi);
693: v = (*shifts)[i];
694: }
695: }
696: }
697: return(0);
698: }
700: static PetscErrorCode EPSSliceGatherSolution(EPS eps)
701: {
702: PetscErrorCode ierr;
703: PetscMPIInt rank,nproc;
704: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
705: PetscInt i,idx,j;
706: PetscInt *perm_loc,off=0,*inertias_loc,ns;
707: PetscScalar *eigr_loc;
708: EPS_SR sr_loc;
709: PetscReal *shifts_loc;
710: PetscMPIInt *disp,*ns_loc,aux;
713: eps->nconv = 0;
714: for (i=0;i<ctx->npart;i++) eps->nconv += ctx->nconv_loc[i];
715: sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
717: /* Gather the shifts used and the inertias computed */
718: EPSSliceGetInertias(ctx->eps,&ns,&shifts_loc,&inertias_loc);
719: if (ctx->sr->dir>0 && shifts_loc[ns-1]==sr_loc->int1 && ctx->subc->color<ctx->npart-1) ns--;
720: if (ctx->sr->dir<0 && shifts_loc[ns-1]==sr_loc->int0 && ctx->subc->color>0) {
721: ns--;
722: for (i=0;i<ns;i++) {
723: inertias_loc[i] = inertias_loc[i+1];
724: shifts_loc[i] = shifts_loc[i+1];
725: }
726: }
727: PetscMalloc1(ctx->npart,&ns_loc);
728: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
729: PetscMPIIntCast(ns,&aux);
730: if (!rank) { MPI_Allgather(&aux,1,MPI_INT,ns_loc,1,MPI_INT,ctx->commrank); }
731: PetscMPIIntCast(ctx->npart,&aux);
732: MPI_Bcast(ns_loc,aux,MPI_INT,0,PetscSubcommChild(ctx->subc));
733: ctx->nshifts = 0;
734: for (i=0;i<ctx->npart;i++) ctx->nshifts += ns_loc[i];
735: PetscFree(ctx->inertias);
736: PetscFree(ctx->shifts);
737: PetscMalloc1(ctx->nshifts,&ctx->inertias);
738: PetscMalloc1(ctx->nshifts,&ctx->shifts);
740: /* Gather eigenvalues (same ranks have fully set of eigenvalues)*/
741: eigr_loc = sr_loc->eigr;
742: perm_loc = sr_loc->perm;
743: MPI_Comm_size(((PetscObject)eps)->comm,&nproc);
744: PetscMalloc1(ctx->npart,&disp);
745: disp[0] = 0;
746: for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ctx->nconv_loc[i-1];
747: if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
748: PetscMPIIntCast(sr_loc->numEigs,&aux);
749: MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank); /* eigenvalues */
750: MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank); /* perm */
751: for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
752: PetscMPIIntCast(ns,&aux);
753: MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank); /* shifts */
754: MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank); /* inertias */
755: MPI_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank);
756: } else { /* subcommunicators with different size */
757: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
758: if (!rank) {
759: PetscMPIIntCast(sr_loc->numEigs,&aux);
760: MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank); /* eigenvalues */
761: MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank); /* perm */
762: for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
763: PetscMPIIntCast(ns,&aux);
764: MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank); /* shifts */
765: MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank); /* inertias */
766: MPI_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank);
767: }
768: PetscMPIIntCast(eps->nconv,&aux);
769: MPI_Bcast(eps->eigr,aux,MPIU_SCALAR,0,PetscSubcommChild(ctx->subc));
770: MPI_Bcast(eps->perm,aux,MPIU_INT,0,PetscSubcommChild(ctx->subc));
771: MPI_Bcast(ctx->shifts,ctx->nshifts,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
772: PetscMPIIntCast(ctx->nshifts,&aux);
773: MPI_Bcast(ctx->inertias,aux,MPIU_INT,0,PetscSubcommChild(ctx->subc));
774: MPI_Bcast(&eps->its,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
775: }
776: /* Update global array eps->perm */
777: idx = ctx->nconv_loc[0];
778: for (i=1;i<ctx->npart;i++) {
779: off += ctx->nconv_loc[i-1];
780: for (j=0;j<ctx->nconv_loc[i];j++) eps->perm[idx++] += off;
781: }
783: /* Gather parallel eigenvectors */
784: PetscFree(ns_loc);
785: PetscFree(disp);
786: PetscFree(shifts_loc);
787: PetscFree(inertias_loc);
788: return(0);
789: }
791: /*
792: Fills the fields of a shift structure
793: */
794: static PetscErrorCode EPSCreateShift(EPS eps,PetscReal val,EPS_shift neighb0,EPS_shift neighb1)
795: {
796: PetscErrorCode ierr;
797: EPS_shift s,*pending2;
798: PetscInt i;
799: EPS_SR sr;
800: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
803: sr = ctx->sr;
804: if ((neighb0 && val==neighb0->value) || (neighb1 && val==neighb1->value)) {
805: sr->nPend++;
806: return(0);
807: }
808: PetscNewLog(eps,&s);
809: s->value = val;
810: s->neighb[0] = neighb0;
811: if (neighb0) neighb0->neighb[1] = s;
812: s->neighb[1] = neighb1;
813: if (neighb1) neighb1->neighb[0] = s;
814: s->comp[0] = PETSC_FALSE;
815: s->comp[1] = PETSC_FALSE;
816: s->index = -1;
817: s->neigs = 0;
818: s->nconv[0] = s->nconv[1] = 0;
819: s->nsch[0] = s->nsch[1]=0;
820: /* Inserts in the stack of pending shifts */
821: /* If needed, the array is resized */
822: if (sr->nPend >= sr->maxPend) {
823: sr->maxPend *= 2;
824: PetscMalloc1(sr->maxPend,&pending2);
825: PetscLogObjectMemory((PetscObject)eps,sizeof(EPS_shift));
826: for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
827: PetscFree(sr->pending);
828: sr->pending = pending2;
829: }
830: sr->pending[sr->nPend++]=s;
831: return(0);
832: }
834: /* Prepare for Rational Krylov update */
835: static PetscErrorCode EPSPrepareRational(EPS eps)
836: {
837: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
838: PetscErrorCode ierr;
839: PetscInt dir,i,k,ld,nv;
840: PetscScalar *A;
841: EPS_SR sr = ctx->sr;
842: Vec v;
845: DSGetLeadingDimension(eps->ds,&ld);
846: dir = (sr->sPres->neighb[0] == sr->sPrev)?1:-1;
847: dir*=sr->dir;
848: k = 0;
849: for (i=0;i<sr->nS;i++) {
850: if (dir*PetscRealPart(sr->S[i])>0.0) {
851: sr->S[k] = sr->S[i];
852: sr->S[sr->nS+k] = sr->S[sr->nS+i];
853: BVGetColumn(sr->Vnext,k,&v);
854: BVCopyVec(eps->V,eps->nconv+i,v);
855: BVRestoreColumn(sr->Vnext,k,&v);
856: k++;
857: if (k>=sr->nS/2)break;
858: }
859: }
860: /* Copy to DS */
861: DSGetArray(eps->ds,DS_MAT_A,&A);
862: PetscArrayzero(A,ld*ld);
863: for (i=0;i<k;i++) {
864: A[i*(1+ld)] = sr->S[i];
865: A[k+i*ld] = sr->S[sr->nS+i];
866: }
867: sr->nS = k;
868: DSRestoreArray(eps->ds,DS_MAT_A,&A);
869: DSGetDimensions(eps->ds,&nv,NULL,NULL,NULL,NULL);
870: DSSetDimensions(eps->ds,nv,0,0,k);
871: /* Append u to V */
872: BVGetColumn(sr->Vnext,sr->nS,&v);
873: BVCopyVec(eps->V,sr->nv,v);
874: BVRestoreColumn(sr->Vnext,sr->nS,&v);
875: return(0);
876: }
878: /* Provides next shift to be computed */
879: static PetscErrorCode EPSExtractShift(EPS eps)
880: {
881: PetscErrorCode ierr;
882: PetscInt iner,zeros=0;
883: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
884: EPS_SR sr;
885: PetscReal newShift,diam,ptol;
886: EPS_shift sPres;
889: sr = ctx->sr;
890: if (sr->nPend > 0) {
891: if (sr->sPres==sr->pending[sr->nPend-1]) {
892: eps->reason = EPS_CONVERGED_ITERATING;
893: eps->its = 0;
894: sr->nPend--;
895: sr->sPres->rep = PETSC_TRUE;
896: return(0);
897: }
898: sr->sPrev = sr->sPres;
899: sr->sPres = sr->pending[--sr->nPend];
900: sPres = sr->sPres;
901: EPSSliceGetInertia(eps,sPres->value,&iner,ctx->detect?&zeros:NULL);
902: if (zeros) {
903: diam = PetscMin(PetscAbsReal(sPres->neighb[0]->value-sPres->value),PetscAbsReal(sPres->neighb[1]->value-sPres->value));
904: ptol = PetscMin(SLICE_PTOL,diam/2);
905: newShift = sPres->value*(1.0+ptol);
906: if (sr->dir*(sPres->neighb[0] && newShift-sPres->neighb[0]->value) < 0) newShift = (sPres->value+sPres->neighb[0]->value)/2;
907: else if (sPres->neighb[1] && sr->dir*(sPres->neighb[1]->value-newShift) < 0) newShift = (sPres->value+sPres->neighb[1]->value)/2;
908: EPSSliceGetInertia(eps,newShift,&iner,&zeros);
909: if (zeros) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",newShift);
910: sPres->value = newShift;
911: }
912: sr->sPres->inertia = iner;
913: eps->target = sr->sPres->value;
914: eps->reason = EPS_CONVERGED_ITERATING;
915: eps->its = 0;
916: } else sr->sPres = NULL;
917: return(0);
918: }
920: /*
921: Symmetric KrylovSchur adapted to spectrum slicing:
922: Allows searching an specific amount of eigenvalues in the subintervals left and right.
923: Returns whether the search has succeeded
924: */
925: static PetscErrorCode EPSKrylovSchur_Slice(EPS eps)
926: {
927: PetscErrorCode ierr;
928: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
929: PetscInt i,k,l,ld,nv,*iwork,j,count0,count1,iterCompl=0,n0,n1;
930: Mat U,Op;
931: PetscScalar *Q,*A;
932: PetscReal *a,*b,beta,lambda;
933: EPS_shift sPres;
934: PetscBool breakdown,complIterating,sch0,sch1;
935: EPS_SR sr = ctx->sr;
936: #define MSGLEN 100
937: char messg[MSGLEN];
940: /* Spectrum slicing data */
941: sPres = sr->sPres;
942: complIterating =PETSC_FALSE;
943: sch1 = sch0 = PETSC_TRUE;
944: DSGetLeadingDimension(eps->ds,&ld);
945: PetscMalloc1(2*ld,&iwork);
946: count0=0;count1=0; /* Found on both sides */
947: if (!sPres->rep && sr->nS > 0 && (sPres->neighb[0] == sr->sPrev || sPres->neighb[1] == sr->sPrev)) {
948: /* Rational Krylov */
949: DSTranslateRKS(eps->ds,sr->sPrev->value-sPres->value);
950: DSGetDimensions(eps->ds,NULL,NULL,NULL,&l,NULL);
951: DSSetDimensions(eps->ds,l+1,0,0,0);
952: BVSetActiveColumns(eps->V,0,l+1);
953: DSGetMat(eps->ds,DS_MAT_Q,&U);
954: BVMultInPlace(eps->V,U,0,l+1);
955: MatDestroy(&U);
956: } else {
957: /* Get the starting Lanczos vector */
958: EPSGetStartVector(eps,0,NULL);
959: l = 0;
960: }
961: /* Restart loop */
962: while (eps->reason == EPS_CONVERGED_ITERATING) {
963: eps->its++; sr->itsKs++;
964: /* Compute an nv-step Lanczos factorization */
965: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
966: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
967: b = a + ld;
968: STGetOperator(eps->st,&Op);
969: BVMatLanczos(eps->V,Op,a,b,eps->nconv+l,&nv,&breakdown);
970: STRestoreOperator(eps->st,&Op);
971: sr->nv = nv;
972: beta = b[nv-1];
973: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
974: DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
975: if (l==0) {
976: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
977: } else {
978: DSSetState(eps->ds,DS_STATE_RAW);
979: }
980: BVSetActiveColumns(eps->V,eps->nconv,nv);
982: /* Solve projected problem and compute residual norm estimates */
983: if (eps->its == 1 && l > 0) {/* After rational update */
984: DSGetArray(eps->ds,DS_MAT_A,&A);
985: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
986: b = a + ld;
987: k = eps->nconv+l;
988: A[k*ld+k-1] = A[(k-1)*ld+k];
989: A[k*ld+k] = a[k];
990: for (j=k+1; j< nv; j++) {
991: A[j*ld+j] = a[j];
992: A[j*ld+j-1] = b[j-1] ;
993: A[(j-1)*ld+j] = b[j-1];
994: }
995: DSRestoreArray(eps->ds,DS_MAT_A,&A);
996: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
997: DSSolve(eps->ds,eps->eigr,NULL);
998: DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
999: DSSetCompact(eps->ds,PETSC_TRUE);
1000: } else { /* Restart */
1001: DSSolve(eps->ds,eps->eigr,NULL);
1002: DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
1003: }
1004: DSSynchronize(eps->ds,eps->eigr,NULL);
1006: /* Residual */
1007: EPSKrylovConvergence(eps,PETSC_TRUE,eps->nconv,nv-eps->nconv,beta,0.0,1.0,&k);
1008: /* Checking values obtained for completing */
1009: for (i=0;i<k;i++) {
1010: sr->back[i]=eps->eigr[i];
1011: }
1012: STBackTransform(eps->st,k,sr->back,eps->eigi);
1013: count0=count1=0;
1014: for (i=0;i<k;i++) {
1015: lambda = PetscRealPart(sr->back[i]);
1016: if (((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0)) count0++;
1017: if (((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0)) count1++;
1018: }
1019: if (k>eps->nev && eps->ncv-k<5) eps->reason = EPS_CONVERGED_TOL;
1020: else {
1021: /* Checks completion */
1022: if ((!sch0||count0 >= sPres->nsch[0]) && (!sch1 ||count1 >= sPres->nsch[1])) {
1023: eps->reason = EPS_CONVERGED_TOL;
1024: } else {
1025: if (!complIterating && eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
1026: if (complIterating) {
1027: if (--iterCompl <= 0) eps->reason = EPS_DIVERGED_ITS;
1028: } else if (k >= eps->nev) {
1029: n0 = sPres->nsch[0]-count0;
1030: n1 = sPres->nsch[1]-count1;
1031: if (sr->iterCompl>0 && ((n0>0 && n0<= sr->nMAXCompl)||(n1>0&&n1<=sr->nMAXCompl))) {
1032: /* Iterating for completion*/
1033: complIterating = PETSC_TRUE;
1034: if (n0 >sr->nMAXCompl)sch0 = PETSC_FALSE;
1035: if (n1 >sr->nMAXCompl)sch1 = PETSC_FALSE;
1036: iterCompl = sr->iterCompl;
1037: } else eps->reason = EPS_CONVERGED_TOL;
1038: }
1039: }
1040: }
1041: /* Update l */
1042: if (eps->reason == EPS_CONVERGED_ITERATING) l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
1043: else l = nv-k;
1044: if (breakdown) l=0;
1045: if (!ctx->lock && l>0 && eps->reason == EPS_CONVERGED_ITERATING) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
1047: if (eps->reason == EPS_CONVERGED_ITERATING) {
1048: if (breakdown) {
1049: /* Start a new Lanczos factorization */
1050: PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
1051: EPSGetStartVector(eps,k,&breakdown);
1052: if (breakdown) {
1053: eps->reason = EPS_DIVERGED_BREAKDOWN;
1054: PetscInfo(eps,"Unable to generate more start vectors\n");
1055: }
1056: } else {
1057: /* Prepare the Rayleigh quotient for restart */
1058: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
1059: DSGetArray(eps->ds,DS_MAT_Q,&Q);
1060: b = a + ld;
1061: for (i=k;i<k+l;i++) {
1062: a[i] = PetscRealPart(eps->eigr[i]);
1063: b[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
1064: }
1065: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
1066: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
1067: }
1068: }
1069: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1070: DSGetMat(eps->ds,DS_MAT_Q,&U);
1071: BVMultInPlace(eps->V,U,eps->nconv,k+l);
1072: MatDestroy(&U);
1074: /* Normalize u and append it to V */
1075: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
1076: BVCopyColumn(eps->V,nv,k+l);
1077: }
1078: eps->nconv = k;
1079: if (eps->reason != EPS_CONVERGED_ITERATING) {
1080: /* Store approximated values for next shift */
1081: DSGetArray(eps->ds,DS_MAT_Q,&Q);
1082: sr->nS = l;
1083: for (i=0;i<l;i++) {
1084: sr->S[i] = eps->eigr[i+k];/* Diagonal elements */
1085: sr->S[i+l] = Q[nv-1+(i+k)*ld]*beta; /* Out of diagonal elements */
1086: }
1087: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
1088: }
1089: }
1090: /* Check for completion */
1091: for (i=0;i< eps->nconv; i++) {
1092: if ((sr->dir)*PetscRealPart(eps->eigr[i])>0) sPres->nconv[1]++;
1093: else sPres->nconv[0]++;
1094: }
1095: sPres->comp[0] = PetscNot(count0 < sPres->nsch[0]);
1096: sPres->comp[1] = PetscNot(count1 < sPres->nsch[1]);
1097: PetscSNPrintf(messg,MSGLEN,"Lanczos: %D evals in [%g,%g]%s and %D evals in [%g,%g]%s\n",count0,(double)(sr->dir==1)?sPres->ext[0]:sPres->value,(double)(sr->dir==1)?sPres->value:sPres->ext[0],(sPres->comp[0])?"*":"",count1,(double)(sr->dir==1)?sPres->value:sPres->ext[1],(double)(sr->dir==1)?sPres->ext[1]:sPres->value,(sPres->comp[1])?"*":"");
1098: PetscInfo(eps,messg);
1099: if (count0 > sPres->nsch[0] || count1 > sPres->nsch[1]) InertiaMismatch(eps,ctx->detect);
1100: PetscFree(iwork);
1101: return(0);
1102: }
1104: /*
1105: Obtains value of subsequent shift
1106: */
1107: static PetscErrorCode EPSGetNewShiftValue(EPS eps,PetscInt side,PetscReal *newS)
1108: {
1109: PetscReal lambda,d_prev;
1110: PetscInt i,idxP;
1111: EPS_SR sr;
1112: EPS_shift sPres,s;
1113: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1116: sr = ctx->sr;
1117: sPres = sr->sPres;
1118: if (sPres->neighb[side]) {
1119: /* Completing a previous interval */
1120: *newS = (sPres->value + sPres->neighb[side]->value)/2;
1121: if (PetscAbsReal(sPres->value - *newS)/PetscAbsReal(sPres->value)<=100*PETSC_SQRT_MACHINE_EPSILON) *newS = sPres->value;
1122: } else { /* (Only for side=1). Creating a new interval. */
1123: if (sPres->neigs==0) {/* No value has been accepted*/
1124: if (sPres->neighb[0]) {
1125: /* Multiplying by 10 the previous distance */
1126: *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
1127: sr->nleap++;
1128: /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
1129: if (!sr->hasEnd && sr->nleap > 5) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to compute the wanted eigenvalues with open interval");
1130: } else { /* First shift */
1131: if (eps->nconv != 0) {
1132: /* Unaccepted values give information for next shift */
1133: idxP=0;/* Number of values left from shift */
1134: for (i=0;i<eps->nconv;i++) {
1135: lambda = PetscRealPart(eps->eigr[i]);
1136: if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
1137: else break;
1138: }
1139: /* Avoiding subtraction of eigenvalues (might be the same).*/
1140: if (idxP>0) {
1141: d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[0]))/(idxP+0.3);
1142: } else {
1143: d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[eps->nconv-1]))/(eps->nconv+0.3);
1144: }
1145: *newS = sPres->value + ((sr->dir)*d_prev*eps->nev)/2;
1146: } else { /* No values found, no information for next shift */
1147: SETERRQ(PetscObjectComm((PetscObject)eps),1,"First shift renders no information");
1148: }
1149: }
1150: } else { /* Accepted values found */
1151: sr->nleap = 0;
1152: /* Average distance of values in previous subinterval */
1153: s = sPres->neighb[0];
1154: while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
1155: s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
1156: }
1157: if (s) {
1158: d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
1159: } else { /* First shift. Average distance obtained with values in this shift */
1160: /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
1161: 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)) {
1162: d_prev = PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0])))/(sPres->neigs+0.3);
1163: } else {
1164: d_prev = PetscAbsReal(PetscRealPart(sr->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
1165: }
1166: }
1167: /* Average distance is used for next shift by adding it to value on the right or to shift */
1168: if ((sr->dir)*(PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
1169: *newS = PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(eps->nev))/2;
1170: } else { /* Last accepted value is on the left of shift. Adding to shift */
1171: *newS = sPres->value + ((sr->dir)*d_prev*(eps->nev))/2;
1172: }
1173: }
1174: /* End of interval can not be surpassed */
1175: if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
1176: }/* of neighb[side]==null */
1177: return(0);
1178: }
1180: /*
1181: Function for sorting an array of real values
1182: */
1183: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
1184: {
1185: PetscReal re;
1186: PetscInt i,j,tmp;
1189: if (!prev) for (i=0;i<nr;i++) perm[i] = i;
1190: /* Insertion sort */
1191: for (i=1;i<nr;i++) {
1192: re = PetscRealPart(r[perm[i]]);
1193: j = i-1;
1194: while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
1195: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
1196: }
1197: }
1198: return(0);
1199: }
1201: /* Stores the pairs obtained since the last shift in the global arrays */
1202: static PetscErrorCode EPSStoreEigenpairs(EPS eps)
1203: {
1204: PetscErrorCode ierr;
1205: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1206: PetscReal lambda,err,norm;
1207: PetscInt i,count;
1208: PetscBool iscayley;
1209: EPS_SR sr = ctx->sr;
1210: EPS_shift sPres;
1211: Vec v,w;
1214: sPres = sr->sPres;
1215: sPres->index = sr->indexEig;
1216: count = sr->indexEig;
1217: /* Back-transform */
1218: STBackTransform(eps->st,eps->nconv,eps->eigr,eps->eigi);
1219: PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
1220: /* Sort eigenvalues */
1221: sortRealEigenvalues(eps->eigr,eps->perm,eps->nconv,PETSC_FALSE,sr->dir);
1222: /* Values stored in global array */
1223: for (i=0;i<eps->nconv;i++) {
1224: lambda = PetscRealPart(eps->eigr[eps->perm[i]]);
1225: err = eps->errest[eps->perm[i]];
1227: if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
1228: if (count>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error in Spectrum Slicing");
1229: sr->eigr[count] = lambda;
1230: sr->errest[count] = err;
1231: /* Explicit purification */
1232: BVGetColumn(eps->V,eps->perm[i],&w);
1233: if (eps->purify) {
1234: BVGetColumn(sr->V,count,&v);
1235: STApply(eps->st,w,v);
1236: BVRestoreColumn(sr->V,count,&v);
1237: } else {
1238: BVInsertVec(sr->V,count,w);
1239: }
1240: BVRestoreColumn(eps->V,eps->perm[i],&w);
1241: BVNormColumn(sr->V,count,NORM_2,&norm);
1242: BVScaleColumn(sr->V,count,1.0/norm);
1243: count++;
1244: }
1245: }
1246: sPres->neigs = count - sr->indexEig;
1247: sr->indexEig = count;
1248: /* Global ordering array updating */
1249: sortRealEigenvalues(sr->eigr,sr->perm,count,PETSC_TRUE,sr->dir);
1250: return(0);
1251: }
1253: static PetscErrorCode EPSLookForDeflation(EPS eps)
1254: {
1255: PetscErrorCode ierr;
1256: PetscReal val;
1257: PetscInt i,count0=0,count1=0;
1258: EPS_shift sPres;
1259: PetscInt ini,fin,k,idx0,idx1;
1260: EPS_SR sr;
1261: Vec v;
1262: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1265: sr = ctx->sr;
1266: sPres = sr->sPres;
1268: if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
1269: else ini = 0;
1270: fin = sr->indexEig;
1271: /* Selection of ends for searching new values */
1272: if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
1273: else sPres->ext[0] = sPres->neighb[0]->value;
1274: if (!sPres->neighb[1]) {
1275: if (sr->hasEnd) sPres->ext[1] = sr->int1;
1276: else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
1277: } else sPres->ext[1] = sPres->neighb[1]->value;
1278: /* Selection of values between right and left ends */
1279: for (i=ini;i<fin;i++) {
1280: val=PetscRealPart(sr->eigr[sr->perm[i]]);
1281: /* Values to the right of left shift */
1282: if ((sr->dir)*(val - sPres->ext[1]) < 0) {
1283: if ((sr->dir)*(val - sPres->value) < 0) count0++;
1284: else count1++;
1285: } else break;
1286: }
1287: /* The number of values on each side are found */
1288: if (sPres->neighb[0]) {
1289: sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
1290: if (sPres->nsch[0]<0) InertiaMismatch(eps,ctx->detect);
1291: } else sPres->nsch[0] = 0;
1293: if (sPres->neighb[1]) {
1294: sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
1295: if (sPres->nsch[1]<0) InertiaMismatch(eps,ctx->detect);
1296: } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);
1298: /* Completing vector of indexes for deflation */
1299: idx0 = ini;
1300: idx1 = ini+count0+count1;
1301: k=0;
1302: for (i=idx0;i<idx1;i++) sr->idxDef[k++]=sr->perm[i];
1303: BVDuplicateResize(eps->V,k+eps->ncv+1,&sr->Vnext);
1304: BVSetNumConstraints(sr->Vnext,k);
1305: for (i=0;i<k;i++) {
1306: BVGetColumn(sr->Vnext,-i-1,&v);
1307: BVCopyVec(sr->V,sr->idxDef[i],v);
1308: BVRestoreColumn(sr->Vnext,-i-1,&v);
1309: }
1311: /* For rational Krylov */
1312: if (sr->nS>0 && (sr->sPrev == sr->sPres->neighb[0] || sr->sPrev == sr->sPres->neighb[1])) {
1313: EPSPrepareRational(eps);
1314: }
1315: eps->nconv = 0;
1316: /* Get rid of temporary Vnext */
1317: BVDestroy(&eps->V);
1318: eps->V = sr->Vnext;
1319: sr->Vnext = NULL;
1320: return(0);
1321: }
1323: PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS eps)
1324: {
1325: PetscErrorCode ierr;
1326: PetscInt i,lds,ti;
1327: PetscReal newS;
1328: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1329: EPS_SR sr=ctx->sr;
1330: Mat A,B=NULL;
1331: PetscObjectState Astate,Bstate=0;
1332: PetscObjectId Aid,Bid=0;
1335: PetscCitationsRegister(citation,&cited);
1336: if (ctx->global) {
1337: EPSSolve_KrylovSchur_Slice(ctx->eps);
1338: ctx->eps->state = EPS_STATE_SOLVED;
1339: eps->reason = EPS_CONVERGED_TOL;
1340: if (ctx->npart>1) {
1341: /* Gather solution from subsolvers */
1342: EPSSliceGatherSolution(eps);
1343: } else {
1344: eps->nconv = sr->numEigs;
1345: eps->its = ctx->eps->its;
1346: PetscFree(ctx->inertias);
1347: PetscFree(ctx->shifts);
1348: EPSSliceGetInertias(ctx->eps,&ctx->nshifts,&ctx->shifts,&ctx->inertias);
1349: }
1350: } else {
1351: if (ctx->npart==1) {
1352: sr->eigr = ctx->eps->eigr;
1353: sr->eigi = ctx->eps->eigi;
1354: sr->perm = ctx->eps->perm;
1355: sr->errest = ctx->eps->errest;
1356: sr->V = ctx->eps->V;
1357: }
1358: /* Check that the user did not modify subcomm matrices */
1359: EPSGetOperators(eps,&A,&B);
1360: PetscObjectStateGet((PetscObject)A,&Astate);
1361: PetscObjectGetId((PetscObject)A,&Aid);
1362: if (B) {
1363: PetscObjectStateGet((PetscObject)B,&Bstate);
1364: PetscObjectGetId((PetscObject)B,&Bid);
1365: }
1366: 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");
1367: /* Only with eigenvalues present in the interval ...*/
1368: if (sr->numEigs==0) {
1369: eps->reason = EPS_CONVERGED_TOL;
1370: return(0);
1371: }
1372: /* Array of pending shifts */
1373: sr->maxPend = 100; /* Initial size */
1374: sr->nPend = 0;
1375: PetscMalloc1(sr->maxPend,&sr->pending);
1376: PetscLogObjectMemory((PetscObject)eps,(sr->maxPend)*sizeof(EPS_shift));
1377: EPSCreateShift(eps,sr->int0,NULL,NULL);
1378: /* extract first shift */
1379: sr->sPrev = NULL;
1380: sr->sPres = sr->pending[--sr->nPend];
1381: sr->sPres->inertia = sr->inertia0;
1382: eps->target = sr->sPres->value;
1383: sr->s0 = sr->sPres;
1384: sr->indexEig = 0;
1385: /* Memory reservation for auxiliary variables */
1386: lds = PetscMin(eps->mpd,eps->ncv);
1387: PetscCalloc1(lds*lds,&sr->S);
1388: PetscMalloc1(eps->ncv,&sr->back);
1389: PetscLogObjectMemory((PetscObject)eps,(sr->numEigs+2*eps->ncv)*sizeof(PetscScalar));
1390: for (i=0;i<sr->numEigs;i++) {
1391: sr->eigr[i] = 0.0;
1392: sr->eigi[i] = 0.0;
1393: sr->errest[i] = 0.0;
1394: sr->perm[i] = i;
1395: }
1396: /* Vectors for deflation */
1397: PetscMalloc1(sr->numEigs,&sr->idxDef);
1398: PetscLogObjectMemory((PetscObject)eps,sr->numEigs*sizeof(PetscInt));
1399: sr->indexEig = 0;
1400: /* Main loop */
1401: while (sr->sPres) {
1402: /* Search for deflation */
1403: EPSLookForDeflation(eps);
1404: /* KrylovSchur */
1405: EPSKrylovSchur_Slice(eps);
1407: EPSStoreEigenpairs(eps);
1408: /* Select new shift */
1409: if (!sr->sPres->comp[1]) {
1410: EPSGetNewShiftValue(eps,1,&newS);
1411: EPSCreateShift(eps,newS,sr->sPres,sr->sPres->neighb[1]);
1412: }
1413: if (!sr->sPres->comp[0]) {
1414: /* Completing earlier interval */
1415: EPSGetNewShiftValue(eps,0,&newS);
1416: EPSCreateShift(eps,newS,sr->sPres->neighb[0],sr->sPres);
1417: }
1418: /* Preparing for a new search of values */
1419: EPSExtractShift(eps);
1420: }
1422: /* Updating eps values prior to exit */
1423: PetscFree(sr->S);
1424: PetscFree(sr->idxDef);
1425: PetscFree(sr->pending);
1426: PetscFree(sr->back);
1427: BVDuplicateResize(eps->V,eps->ncv+1,&sr->Vnext);
1428: BVSetNumConstraints(sr->Vnext,0);
1429: BVDestroy(&eps->V);
1430: eps->V = sr->Vnext;
1431: eps->nconv = sr->indexEig;
1432: eps->reason = EPS_CONVERGED_TOL;
1433: eps->its = sr->itsKs;
1434: eps->nds = 0;
1435: if (sr->dir<0) {
1436: for (i=0;i<eps->nconv/2;i++) {
1437: ti = sr->perm[i]; sr->perm[i] = sr->perm[eps->nconv-1-i]; sr->perm[eps->nconv-1-i] = ti;
1438: }
1439: }
1440: }
1441: return(0);
1442: }