Actual source code: ks-slice.c
slepc-3.14.2 2021-02-01
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: PetscBool asymm,bsymm,aherm,bherm;
172: PetscMPIInt rank;
173: EPS_SR sr=ctx->sr;
174: PC pc;
175: PCType pctype;
176: KSP ksp;
177: KSPType ksptype;
178: STType sttype;
179: PetscObjectState Astate,Bstate=0;
180: PetscObjectId Aid,Bid=0;
181: MatSolverType stype;
184: EPSGetOperators(eps,&A,&B);
185: if (ctx->npart==1) {
186: if (!ctx->eps) { EPSCreate(((PetscObject)eps)->comm,&ctx->eps); }
187: EPSSetType(ctx->eps,((PetscObject)eps)->type_name);
188: EPSSetOperators(ctx->eps,A,B);
189: a = eps->inta; b = eps->intb;
190: } else {
191: PetscObjectStateGet((PetscObject)A,&Astate);
192: PetscObjectGetId((PetscObject)A,&Aid);
193: MatGetOption(A,MAT_SYMMETRIC,&asymm);
194: MatGetOption(A,MAT_HERMITIAN,&aherm);
195: if (B) {
196: PetscObjectStateGet((PetscObject)B,&Bstate);
197: PetscObjectGetId((PetscObject)B,&Bid);
198: MatGetOption(B,MAT_SYMMETRIC,&bsymm);
199: MatGetOption(B,MAT_HERMITIAN,&bherm);
200: }
201: if (!ctx->subc) {
202: /* Create context for subcommunicators */
203: PetscSubcommCreate(PetscObjectComm((PetscObject)eps),&ctx->subc);
204: PetscSubcommSetNumber(ctx->subc,ctx->npart);
205: PetscSubcommSetType(ctx->subc,PETSC_SUBCOMM_CONTIGUOUS);
206: PetscLogObjectMemory((PetscObject)eps,sizeof(PetscSubcomm));
208: /* Duplicate matrices */
209: MatCreateRedundantMatrix(A,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Ar);
210: ctx->Astate = Astate;
211: ctx->Aid = Aid;
212: MatSetOption(Ar,MAT_SYMMETRIC,asymm);
213: MatSetOption(Ar,MAT_HERMITIAN,aherm);
214: if (B) {
215: MatCreateRedundantMatrix(B,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Br);
216: ctx->Bstate = Bstate;
217: ctx->Bid = Bid;
218: MatSetOption(Br,MAT_SYMMETRIC,bsymm);
219: MatSetOption(Br,MAT_HERMITIAN,bherm);
220: }
221: } else {
222: if (ctx->Astate != Astate || (B && ctx->Bstate != Bstate) || ctx->Aid != Aid || (B && ctx->Bid != Bid)) {
223: EPSGetOperators(ctx->eps,&Ar,&Br);
224: MatCreateRedundantMatrix(A,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Ar);
225: ctx->Astate = Astate;
226: ctx->Aid = Aid;
227: MatSetOption(Ar,MAT_SYMMETRIC,asymm);
228: MatSetOption(Ar,MAT_HERMITIAN,aherm);
229: if (B) {
230: MatCreateRedundantMatrix(B,0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&Br);
231: ctx->Bstate = Bstate;
232: ctx->Bid = Bid;
233: MatSetOption(Br,MAT_SYMMETRIC,bsymm);
234: MatSetOption(Br,MAT_HERMITIAN,bherm);
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: PCFactorGetMatSolverType(pc,&stype);
281: PCFactorGetZeroPivot(pc,&zero);
282: STGetKSP(ctx->eps->st,&ksp);
283: KSPSetType(ksp,ksptype);
284: KSPGetPC(ksp,&pc);
285: PCSetType(pc,pctype);
286: PCFactorSetZeroPivot(pc,zero);
287: if (stype) { PCFactorSetMatSolverType(pc,stype); }
289: EPSSetConvergenceTest(ctx->eps,eps->conv);
290: EPSSetInterval(ctx->eps,a,b);
291: ctx_local = (EPS_KRYLOVSCHUR*)ctx->eps->data;
292: ctx_local->npart = ctx->npart;
293: ctx_local->detect = ctx->detect;
294: ctx_local->global = PETSC_FALSE;
295: ctx_local->eps = eps;
296: ctx_local->subc = ctx->subc;
297: ctx_local->commrank = ctx->commrank;
299: EPSSetDimensions(ctx->eps,ctx->nev,ctx->ncv,ctx->mpd);
300: EPSKrylovSchurSetLocking(ctx->eps,ctx->lock);
302: /* transfer options from eps->V */
303: EPSGetBV(ctx->eps,&V);
304: if (!eps->V) { EPSGetBV(eps,&eps->V); }
305: if (!((PetscObject)(eps->V))->type_name) {
306: BVSetType(V,BVSVEC);
307: } else {
308: BVGetType(eps->V,&type);
309: BVSetType(V,type);
310: }
311: BVGetOrthogonalization(eps->V,&orthog_type,&orthog_ref,&eta,&ob_type);
312: BVSetOrthogonalization(V,orthog_type,orthog_ref,eta,ob_type);
313: ctx->eps->which = eps->which;
314: ctx->eps->max_it = eps->max_it;
315: ctx->eps->tol = eps->tol;
316: ctx->eps->purify = eps->purify;
317: if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
318: EPSSetProblemType(ctx->eps,eps->problem_type);
319: EPSSetUp(ctx->eps);
320: ctx->eps->nconv = 0;
321: ctx->eps->its = 0;
322: for (i=0;i<ctx->eps->ncv;i++) {
323: ctx->eps->eigr[i] = 0.0;
324: ctx->eps->eigi[i] = 0.0;
325: ctx->eps->errest[i] = 0.0;
326: }
327: return(0);
328: }
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=shift;
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: STGetKSP(eps->st,&ksp);
351: KSPGetPC(ksp,&pc);
352: PCFactorGetMatrix(pc,&F);
353: MatGetInertia(F,inertia,zeros,NULL);
354: }
355: if (inertia) { PetscInfo2(eps,"Computed inertia at shift %g: %D\n",(double)nzshift,*inertia); }
356: return(0);
357: }
359: /*
360: Dummy backtransform operation
361: */
362: static PetscErrorCode EPSBackTransform_Skip(EPS eps)
363: {
365: return(0);
366: }
368: PetscErrorCode EPSSetUp_KrylovSchur_Slice(EPS eps)
369: {
370: PetscErrorCode ierr;
371: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data,*ctx_glob;
372: EPS_SR sr,sr_loc,sr_glob;
373: PetscInt nEigs,dssz=1,i,zeros=0,off=0,method,hiteig=0;
374: PetscMPIInt nproc,rank=0,aux;
375: PetscReal r;
376: MPI_Request req;
377: Mat A,B=NULL;
378: DSParallelType ptype;
381: if (ctx->global) {
382: EPSCheckHermitianDefiniteCondition(eps,PETSC_TRUE," with spectrum slicing");
383: EPSCheckSinvertCayleyCondition(eps,PETSC_TRUE," with spectrum slicing");
384: if (eps->inta==eps->intb) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues unless you provide a computational interval with EPSSetInterval()");
385: 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");
386: EPSCheckUnsupportedCondition(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING,PETSC_TRUE," with spectrum slicing");
387: EPSCheckIgnoredCondition(eps,EPS_FEATURE_BALANCE,PETSC_TRUE," with spectrum slicing");
388: if (eps->tol==PETSC_DEFAULT) {
389: #if defined(PETSC_USE_REAL_SINGLE)
390: eps->tol = SLEPC_DEFAULT_TOL;
391: #else
392: /* use tighter tolerance */
393: eps->tol = SLEPC_DEFAULT_TOL*1e-2;
394: #endif
395: }
396: if (eps->max_it==PETSC_DEFAULT) eps->max_it = 100;
397: if (ctx->nev==1) ctx->nev = PetscMin(40,eps->n); /* nev not set, use default value */
398: 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");
399: }
400: eps->ops->backtransform = EPSBackTransform_Skip;
402: /* create spectrum slicing context and initialize it */
403: EPSSliceResetSR(eps);
404: PetscNewLog(eps,&sr);
405: ctx->sr = sr;
406: sr->itsKs = 0;
407: sr->nleap = 0;
408: sr->nMAXCompl = eps->nev/4;
409: sr->iterCompl = eps->max_it/4;
410: sr->sPres = NULL;
411: sr->nS = 0;
413: if (ctx->npart==1 || ctx->global) {
414: /* check presence of ends and finding direction */
415: if ((eps->inta > PETSC_MIN_REAL && !(ctx->subintervals && ctx->subintervals[0]==ctx->subintervals[1])) || eps->intb >= PETSC_MAX_REAL) {
416: sr->int0 = eps->inta;
417: sr->int1 = eps->intb;
418: sr->dir = 1;
419: if (eps->intb >= PETSC_MAX_REAL) { /* Right-open interval */
420: sr->hasEnd = PETSC_FALSE;
421: } else sr->hasEnd = PETSC_TRUE;
422: } else {
423: sr->int0 = eps->intb;
424: sr->int1 = eps->inta;
425: sr->dir = -1;
426: sr->hasEnd = PetscNot(eps->inta <= PETSC_MIN_REAL);
427: }
428: }
429: if (ctx->global) {
430: /* prevent computation of factorization in global eps */
431: STSetTransform(eps->st,PETSC_FALSE);
432: EPSSetDimensions_Default(eps,ctx->nev,&ctx->ncv,&ctx->mpd);
433: /* create subintervals and initialize auxiliary eps for slicing runs */
434: EPSSliceGetEPS(eps);
435: sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
436: if (ctx->npart>1) {
437: if ((sr->dir>0&&ctx->subc->color==0)||(sr->dir<0&&ctx->subc->color==ctx->npart-1)) sr->inertia0 = sr_loc->inertia0;
438: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
439: if (!rank) {
440: MPI_Bcast(&sr->inertia0,1,MPIU_INT,(sr->dir>0)?0:ctx->npart-1,ctx->commrank);
441: }
442: MPI_Bcast(&sr->inertia0,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
443: PetscFree(ctx->nconv_loc);
444: PetscMalloc1(ctx->npart,&ctx->nconv_loc);
445: MPI_Comm_size(((PetscObject)eps)->comm,&nproc);
446: if (sr->dir<0) off = 1;
447: if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
448: PetscMPIIntCast(sr_loc->numEigs,&aux);
449: MPI_Allgather(&aux,1,MPI_INT,ctx->nconv_loc,1,MPI_INT,ctx->commrank);
450: MPI_Allgather(sr_loc->dir==sr->dir?&sr_loc->int0:&sr_loc->int1,1,MPIU_REAL,ctx->subintervals+off,1,MPIU_REAL,ctx->commrank);
451: } else {
452: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
453: if (!rank) {
454: PetscMPIIntCast(sr_loc->numEigs,&aux);
455: MPI_Allgather(&aux,1,MPI_INT,ctx->nconv_loc,1,MPI_INT,ctx->commrank);
456: MPI_Allgather(sr_loc->dir==sr->dir?&sr_loc->int0:&sr_loc->int1,1,MPIU_REAL,ctx->subintervals+off,1,MPIU_REAL,ctx->commrank);
457: }
458: PetscMPIIntCast(ctx->npart,&aux);
459: MPI_Bcast(ctx->nconv_loc,aux,MPI_INT,0,PetscSubcommChild(ctx->subc));
460: MPI_Bcast(ctx->subintervals+off,aux,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
461: }
462: nEigs = 0;
463: for (i=0;i<ctx->npart;i++) nEigs += ctx->nconv_loc[i];
464: } else {
465: nEigs = sr_loc->numEigs;
466: sr->inertia0 = sr_loc->inertia0;
467: sr->dir = sr_loc->dir;
468: }
469: sr->inertia1 = sr->inertia0+sr->dir*nEigs;
470: sr->numEigs = nEigs;
471: eps->nev = nEigs;
472: eps->ncv = nEigs;
473: eps->mpd = nEigs;
474: } else {
475: ctx_glob = (EPS_KRYLOVSCHUR*)ctx->eps->data;
476: sr_glob = ctx_glob->sr;
477: if (ctx->npart>1) {
478: sr->dir = sr_glob->dir;
479: sr->int0 = (sr->dir==1)?eps->inta:eps->intb;
480: sr->int1 = (sr->dir==1)?eps->intb:eps->inta;
481: if ((sr->dir>0&&ctx->subc->color==ctx->npart-1)||(sr->dir<0&&ctx->subc->color==0)) sr->hasEnd = sr_glob->hasEnd;
482: else sr->hasEnd = PETSC_TRUE;
483: }
484: STSetUp(eps->st);
486: /* compute inertia0 */
487: EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,ctx->detect?&zeros:NULL);
488: /* undocumented option to control what to do when an eigenvalue is found:
489: - error out if it's the endpoint of the user-provided interval (or sub-interval)
490: - if it's an endpoint computed internally:
491: + if hiteig=0 error out
492: + else if hiteig=1 the subgroup that hit the eigenvalue does nothing
493: + otherwise the subgroup that hit the eigenvalue perturbs the shift and recomputes inertia
494: */
495: PetscOptionsGetInt(NULL,NULL,"-eps_krylovschur_hiteigenvalue",&hiteig,NULL);
496: if (zeros) { /* error in factorization */
497: 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");
498: 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");
499: else {
500: if (hiteig==1) { /* idle subgroup */
501: sr->inertia0 = -1;
502: } else { /* perturb shift */
503: sr->int0 *= (1.0+SLICE_PTOL);
504: EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,&zeros);
505: if (zeros) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",sr->int1);
506: }
507: }
508: }
509: if (ctx->npart>1) {
510: /* inertia1 is received from neighbour */
511: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
512: if (!rank) {
513: if (sr->inertia0!=-1 && ((sr->dir>0 && ctx->subc->color>0) || (sr->dir<0 && ctx->subc->color<ctx->npart-1))) { /* send inertia0 to neighbour0 */
514: MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
515: MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
516: }
517: if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)|| (sr->dir<0 && ctx->subc->color>0)) { /* receive inertia1 from neighbour1 */
518: MPI_Recv(&(sr->inertia1),1,MPIU_INT,ctx->subc->color+sr->dir,0,ctx->commrank,MPI_STATUS_IGNORE);
519: MPI_Recv(&(sr->int1),1,MPIU_REAL,ctx->subc->color+sr->dir,0,ctx->commrank,MPI_STATUS_IGNORE);
520: }
521: if (sr->inertia0==-1 && !(sr->dir>0 && ctx->subc->color==ctx->npart-1) && !(sr->dir<0 && ctx->subc->color==0)) {
522: sr->inertia0 = sr->inertia1; sr->int0 = sr->int1;
523: MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
524: MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
525: }
526: }
527: if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)||(sr->dir<0 && ctx->subc->color>0)) {
528: MPI_Bcast(&sr->inertia1,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
529: MPI_Bcast(&sr->int1,1,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
530: } else sr_glob->inertia1 = sr->inertia1;
531: }
533: /* last process in eps comm computes inertia1 */
534: if (ctx->npart==1 || ((sr->dir>0 && ctx->subc->color==ctx->npart-1) || (sr->dir<0 && ctx->subc->color==0))) {
535: EPSSliceGetInertia(eps,sr->int1,&sr->inertia1,ctx->detect?&zeros:NULL);
536: if (zeros) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
537: if (!rank && sr->inertia0==-1) {
538: sr->inertia0 = sr->inertia1; sr->int0 = sr->int1;
539: MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
540: MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
541: }
542: if (sr->hasEnd) {
543: sr->dir = -sr->dir; r = sr->int0; sr->int0 = sr->int1; sr->int1 = r;
544: i = sr->inertia0; sr->inertia0 = sr->inertia1; sr->inertia1 = i;
545: }
546: }
548: /* number of eigenvalues in interval */
549: sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
550: if (ctx->npart>1) {
551: /* memory allocate for subinterval eigenpairs */
552: EPSSliceAllocateSolution(eps,1);
553: }
554: dssz = eps->ncv+1;
555: DSGetParallel(ctx->eps->ds,&ptype);
556: DSSetParallel(eps->ds,ptype);
557: DSGetMethod(ctx->eps->ds,&method);
558: DSSetMethod(eps->ds,method);
559: }
560: DSSetType(eps->ds,DSHEP);
561: DSSetCompact(eps->ds,PETSC_TRUE);
562: DSAllocate(eps->ds,dssz);
563: /* keep state of subcomm matrices to check that the user does not modify them */
564: EPSGetOperators(eps,&A,&B);
565: PetscObjectStateGet((PetscObject)A,&ctx->Astate);
566: PetscObjectGetId((PetscObject)A,&ctx->Aid);
567: if (B) {
568: PetscObjectStateGet((PetscObject)B,&ctx->Bstate);
569: PetscObjectGetId((PetscObject)B,&ctx->Bid);
570: } else {
571: ctx->Bstate=0;
572: ctx->Bid=0;
573: }
574: return(0);
575: }
577: static PetscErrorCode EPSSliceGatherEigenVectors(EPS eps)
578: {
579: PetscErrorCode ierr;
580: Vec v,vg,v_loc;
581: IS is1,is2;
582: VecScatter vec_sc;
583: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
584: PetscInt nloc,m0,n0,i,si,idx,*idx1,*idx2,j;
585: PetscScalar *array;
586: EPS_SR sr_loc;
587: BV V_loc;
590: sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
591: V_loc = sr_loc->V;
593: /* Gather parallel eigenvectors */
594: BVGetColumn(eps->V,0,&v);
595: VecGetOwnershipRange(v,&n0,&m0);
596: BVRestoreColumn(eps->V,0,&v);
597: BVGetColumn(ctx->eps->V,0,&v);
598: VecGetLocalSize(v,&nloc);
599: BVRestoreColumn(ctx->eps->V,0,&v);
600: PetscMalloc2(m0-n0,&idx1,m0-n0,&idx2);
601: VecCreateMPI(PetscObjectComm((PetscObject)eps),nloc,PETSC_DECIDE,&vg);
602: idx = -1;
603: for (si=0;si<ctx->npart;si++) {
604: j = 0;
605: for (i=n0;i<m0;i++) {
606: idx1[j] = i;
607: idx2[j++] = i+eps->n*si;
608: }
609: ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx1,PETSC_COPY_VALUES,&is1);
610: ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx2,PETSC_COPY_VALUES,&is2);
611: BVGetColumn(eps->V,0,&v);
612: VecScatterCreate(v,is1,vg,is2,&vec_sc);
613: BVRestoreColumn(eps->V,0,&v);
614: ISDestroy(&is1);
615: ISDestroy(&is2);
616: for (i=0;i<ctx->nconv_loc[si];i++) {
617: BVGetColumn(eps->V,++idx,&v);
618: if (ctx->subc->color==si) {
619: BVGetColumn(V_loc,i,&v_loc);
620: VecGetArray(v_loc,&array);
621: VecPlaceArray(vg,array);
622: }
623: VecScatterBegin(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE);
624: VecScatterEnd(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE);
625: if (ctx->subc->color==si) {
626: VecResetArray(vg);
627: VecRestoreArray(v_loc,&array);
628: BVRestoreColumn(V_loc,i,&v_loc);
629: }
630: BVRestoreColumn(eps->V,idx,&v);
631: }
632: VecScatterDestroy(&vec_sc);
633: }
634: PetscFree2(idx1,idx2);
635: VecDestroy(&vg);
636: return(0);
637: }
639: /*
640: EPSComputeVectors_Slice - Recover Eigenvectors from subcomunicators
641: */
642: PetscErrorCode EPSComputeVectors_Slice(EPS eps)
643: {
644: PetscErrorCode ierr;
645: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
648: if (ctx->global && ctx->npart>1) {
649: EPSComputeVectors(ctx->eps);
650: EPSSliceGatherEigenVectors(eps);
651: }
652: return(0);
653: }
655: #define SWAP(a,b,t) {t=a;a=b;b=t;}
657: static PetscErrorCode EPSSliceGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
658: {
659: PetscErrorCode ierr;
660: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
661: PetscInt i=0,j,tmpi;
662: PetscReal v,tmpr;
663: EPS_shift s;
666: if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
667: if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
668: if (!ctx->sr->s0) { /* EPSSolve not called yet */
669: *n = 2;
670: } else {
671: *n = 1;
672: s = ctx->sr->s0;
673: while (s) {
674: (*n)++;
675: s = s->neighb[1];
676: }
677: }
678: PetscMalloc1(*n,shifts);
679: PetscMalloc1(*n,inertias);
680: if (!ctx->sr->s0) { /* EPSSolve not called yet */
681: (*shifts)[0] = ctx->sr->int0;
682: (*shifts)[1] = ctx->sr->int1;
683: (*inertias)[0] = ctx->sr->inertia0;
684: (*inertias)[1] = ctx->sr->inertia1;
685: } else {
686: s = ctx->sr->s0;
687: while (s) {
688: (*shifts)[i] = s->value;
689: (*inertias)[i++] = s->inertia;
690: s = s->neighb[1];
691: }
692: (*shifts)[i] = ctx->sr->int1;
693: (*inertias)[i] = ctx->sr->inertia1;
694: }
695: /* remove possible duplicate in last position */
696: if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
697: /* sort result */
698: for (i=0;i<*n;i++) {
699: v = (*shifts)[i];
700: for (j=i+1;j<*n;j++) {
701: if (v > (*shifts)[j]) {
702: SWAP((*shifts)[i],(*shifts)[j],tmpr);
703: SWAP((*inertias)[i],(*inertias)[j],tmpi);
704: v = (*shifts)[i];
705: }
706: }
707: }
708: return(0);
709: }
711: static PetscErrorCode EPSSliceGatherSolution(EPS eps)
712: {
713: PetscErrorCode ierr;
714: PetscMPIInt rank,nproc;
715: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
716: PetscInt i,idx,j;
717: PetscInt *perm_loc,off=0,*inertias_loc,ns;
718: PetscScalar *eigr_loc;
719: EPS_SR sr_loc;
720: PetscReal *shifts_loc;
721: PetscMPIInt *disp,*ns_loc,aux;
724: eps->nconv = 0;
725: for (i=0;i<ctx->npart;i++) eps->nconv += ctx->nconv_loc[i];
726: sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
728: /* Gather the shifts used and the inertias computed */
729: EPSSliceGetInertias(ctx->eps,&ns,&shifts_loc,&inertias_loc);
730: if (ctx->sr->dir>0 && shifts_loc[ns-1]==sr_loc->int1 && ctx->subc->color<ctx->npart-1) ns--;
731: if (ctx->sr->dir<0 && shifts_loc[ns-1]==sr_loc->int0 && ctx->subc->color>0) {
732: ns--;
733: for (i=0;i<ns;i++) {
734: inertias_loc[i] = inertias_loc[i+1];
735: shifts_loc[i] = shifts_loc[i+1];
736: }
737: }
738: PetscMalloc1(ctx->npart,&ns_loc);
739: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
740: PetscMPIIntCast(ns,&aux);
741: if (!rank) { MPI_Allgather(&aux,1,MPI_INT,ns_loc,1,MPI_INT,ctx->commrank); }
742: PetscMPIIntCast(ctx->npart,&aux);
743: MPI_Bcast(ns_loc,aux,MPI_INT,0,PetscSubcommChild(ctx->subc));
744: ctx->nshifts = 0;
745: for (i=0;i<ctx->npart;i++) ctx->nshifts += ns_loc[i];
746: PetscFree(ctx->inertias);
747: PetscFree(ctx->shifts);
748: PetscMalloc1(ctx->nshifts,&ctx->inertias);
749: PetscMalloc1(ctx->nshifts,&ctx->shifts);
751: /* Gather eigenvalues (same ranks have fully set of eigenvalues)*/
752: eigr_loc = sr_loc->eigr;
753: perm_loc = sr_loc->perm;
754: MPI_Comm_size(((PetscObject)eps)->comm,&nproc);
755: PetscMalloc1(ctx->npart,&disp);
756: disp[0] = 0;
757: for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ctx->nconv_loc[i-1];
758: if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
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: MPIU_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank);
767: } else { /* subcommunicators with different size */
768: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
769: if (!rank) {
770: PetscMPIIntCast(sr_loc->numEigs,&aux);
771: MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank); /* eigenvalues */
772: MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank); /* perm */
773: for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
774: PetscMPIIntCast(ns,&aux);
775: MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank); /* shifts */
776: MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank); /* inertias */
777: MPIU_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank);
778: }
779: PetscMPIIntCast(eps->nconv,&aux);
780: MPI_Bcast(eps->eigr,aux,MPIU_SCALAR,0,PetscSubcommChild(ctx->subc));
781: MPI_Bcast(eps->perm,aux,MPIU_INT,0,PetscSubcommChild(ctx->subc));
782: MPI_Bcast(ctx->shifts,ctx->nshifts,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
783: PetscMPIIntCast(ctx->nshifts,&aux);
784: MPI_Bcast(ctx->inertias,aux,MPIU_INT,0,PetscSubcommChild(ctx->subc));
785: MPI_Bcast(&eps->its,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
786: }
787: /* Update global array eps->perm */
788: idx = ctx->nconv_loc[0];
789: for (i=1;i<ctx->npart;i++) {
790: off += ctx->nconv_loc[i-1];
791: for (j=0;j<ctx->nconv_loc[i];j++) eps->perm[idx++] += off;
792: }
794: /* Gather parallel eigenvectors */
795: PetscFree(ns_loc);
796: PetscFree(disp);
797: PetscFree(shifts_loc);
798: PetscFree(inertias_loc);
799: return(0);
800: }
802: /*
803: Fills the fields of a shift structure
804: */
805: static PetscErrorCode EPSCreateShift(EPS eps,PetscReal val,EPS_shift neighb0,EPS_shift neighb1)
806: {
807: PetscErrorCode ierr;
808: EPS_shift s,*pending2;
809: PetscInt i;
810: EPS_SR sr;
811: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
814: sr = ctx->sr;
815: if ((neighb0 && val==neighb0->value) || (neighb1 && val==neighb1->value)) {
816: sr->nPend++;
817: return(0);
818: }
819: PetscNewLog(eps,&s);
820: s->value = val;
821: s->neighb[0] = neighb0;
822: if (neighb0) neighb0->neighb[1] = s;
823: s->neighb[1] = neighb1;
824: if (neighb1) neighb1->neighb[0] = s;
825: s->comp[0] = PETSC_FALSE;
826: s->comp[1] = PETSC_FALSE;
827: s->index = -1;
828: s->neigs = 0;
829: s->nconv[0] = s->nconv[1] = 0;
830: s->nsch[0] = s->nsch[1]=0;
831: /* Inserts in the stack of pending shifts */
832: /* If needed, the array is resized */
833: if (sr->nPend >= sr->maxPend) {
834: sr->maxPend *= 2;
835: PetscMalloc1(sr->maxPend,&pending2);
836: PetscLogObjectMemory((PetscObject)eps,sizeof(EPS_shift));
837: for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
838: PetscFree(sr->pending);
839: sr->pending = pending2;
840: }
841: sr->pending[sr->nPend++]=s;
842: return(0);
843: }
845: /* Prepare for Rational Krylov update */
846: static PetscErrorCode EPSPrepareRational(EPS eps)
847: {
848: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
849: PetscErrorCode ierr;
850: PetscInt dir,i,k,ld,nv;
851: PetscScalar *A;
852: EPS_SR sr = ctx->sr;
853: Vec v;
856: DSGetLeadingDimension(eps->ds,&ld);
857: dir = (sr->sPres->neighb[0] == sr->sPrev)?1:-1;
858: dir*=sr->dir;
859: k = 0;
860: for (i=0;i<sr->nS;i++) {
861: if (dir*PetscRealPart(sr->S[i])>0.0) {
862: sr->S[k] = sr->S[i];
863: sr->S[sr->nS+k] = sr->S[sr->nS+i];
864: BVGetColumn(sr->Vnext,k,&v);
865: BVCopyVec(eps->V,eps->nconv+i,v);
866: BVRestoreColumn(sr->Vnext,k,&v);
867: k++;
868: if (k>=sr->nS/2)break;
869: }
870: }
871: /* Copy to DS */
872: DSGetArray(eps->ds,DS_MAT_A,&A);
873: PetscArrayzero(A,ld*ld);
874: for (i=0;i<k;i++) {
875: A[i*(1+ld)] = sr->S[i];
876: A[k+i*ld] = sr->S[sr->nS+i];
877: }
878: sr->nS = k;
879: DSRestoreArray(eps->ds,DS_MAT_A,&A);
880: DSGetDimensions(eps->ds,&nv,NULL,NULL,NULL,NULL);
881: DSSetDimensions(eps->ds,nv,0,0,k);
882: /* Append u to V */
883: BVGetColumn(sr->Vnext,sr->nS,&v);
884: BVCopyVec(eps->V,sr->nv,v);
885: BVRestoreColumn(sr->Vnext,sr->nS,&v);
886: return(0);
887: }
889: /* Provides next shift to be computed */
890: static PetscErrorCode EPSExtractShift(EPS eps)
891: {
892: PetscErrorCode ierr;
893: PetscInt iner,zeros=0;
894: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
895: EPS_SR sr;
896: PetscReal newShift,diam,ptol;
897: EPS_shift sPres;
900: sr = ctx->sr;
901: if (sr->nPend > 0) {
902: if (sr->sPres==sr->pending[sr->nPend-1]) {
903: eps->reason = EPS_CONVERGED_ITERATING;
904: eps->its = 0;
905: sr->nPend--;
906: sr->sPres->rep = PETSC_TRUE;
907: return(0);
908: }
909: sr->sPrev = sr->sPres;
910: sr->sPres = sr->pending[--sr->nPend];
911: sPres = sr->sPres;
912: EPSSliceGetInertia(eps,sPres->value,&iner,ctx->detect?&zeros:NULL);
913: if (zeros) {
914: diam = PetscMin(PetscAbsReal(sPres->neighb[0]->value-sPres->value),PetscAbsReal(sPres->neighb[1]->value-sPres->value));
915: ptol = PetscMin(SLICE_PTOL,diam/2);
916: newShift = sPres->value*(1.0+ptol);
917: if (sr->dir*(sPres->neighb[0] && newShift-sPres->neighb[0]->value) < 0) newShift = (sPres->value+sPres->neighb[0]->value)/2;
918: else if (sPres->neighb[1] && sr->dir*(sPres->neighb[1]->value-newShift) < 0) newShift = (sPres->value+sPres->neighb[1]->value)/2;
919: EPSSliceGetInertia(eps,newShift,&iner,&zeros);
920: if (zeros) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",newShift);
921: sPres->value = newShift;
922: }
923: sr->sPres->inertia = iner;
924: eps->target = sr->sPres->value;
925: eps->reason = EPS_CONVERGED_ITERATING;
926: eps->its = 0;
927: } else sr->sPres = NULL;
928: return(0);
929: }
931: /*
932: Symmetric KrylovSchur adapted to spectrum slicing:
933: Allows searching an specific amount of eigenvalues in the subintervals left and right.
934: Returns whether the search has succeeded
935: */
936: static PetscErrorCode EPSKrylovSchur_Slice(EPS eps)
937: {
938: PetscErrorCode ierr;
939: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
940: PetscInt i,k,l,ld,nv,*iwork,j,count0,count1,iterCompl=0,n0,n1;
941: Mat U,Op;
942: PetscScalar *Q,*A;
943: PetscReal *a,*b,beta,lambda;
944: EPS_shift sPres;
945: PetscBool breakdown,complIterating,sch0,sch1;
946: EPS_SR sr = ctx->sr;
947: char messg[100];
950: /* Spectrum slicing data */
951: sPres = sr->sPres;
952: complIterating =PETSC_FALSE;
953: sch1 = sch0 = PETSC_TRUE;
954: DSGetLeadingDimension(eps->ds,&ld);
955: PetscMalloc1(2*ld,&iwork);
956: count0=0;count1=0; /* Found on both sides */
957: if (!sPres->rep && sr->nS > 0 && (sPres->neighb[0] == sr->sPrev || sPres->neighb[1] == sr->sPrev)) {
958: /* Rational Krylov */
959: DSTranslateRKS(eps->ds,sr->sPrev->value-sPres->value);
960: DSGetDimensions(eps->ds,NULL,NULL,NULL,&l,NULL);
961: DSSetDimensions(eps->ds,l+1,0,0,0);
962: BVSetActiveColumns(eps->V,0,l+1);
963: DSGetMat(eps->ds,DS_MAT_Q,&U);
964: BVMultInPlace(eps->V,U,0,l+1);
965: MatDestroy(&U);
966: } else {
967: /* Get the starting Lanczos vector */
968: EPSGetStartVector(eps,0,NULL);
969: l = 0;
970: }
971: /* Restart loop */
972: while (eps->reason == EPS_CONVERGED_ITERATING) {
973: eps->its++; sr->itsKs++;
974: /* Compute an nv-step Lanczos factorization */
975: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
976: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
977: b = a + ld;
978: STGetOperator(eps->st,&Op);
979: BVMatLanczos(eps->V,Op,a,b,eps->nconv+l,&nv,&breakdown);
980: STRestoreOperator(eps->st,&Op);
981: sr->nv = nv;
982: beta = b[nv-1];
983: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
984: DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
985: if (l==0) {
986: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
987: } else {
988: DSSetState(eps->ds,DS_STATE_RAW);
989: }
990: BVSetActiveColumns(eps->V,eps->nconv,nv);
992: /* Solve projected problem and compute residual norm estimates */
993: if (eps->its == 1 && l > 0) {/* After rational update */
994: DSGetArray(eps->ds,DS_MAT_A,&A);
995: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
996: b = a + ld;
997: k = eps->nconv+l;
998: A[k*ld+k-1] = A[(k-1)*ld+k];
999: A[k*ld+k] = a[k];
1000: for (j=k+1; j< nv; j++) {
1001: A[j*ld+j] = a[j];
1002: A[j*ld+j-1] = b[j-1] ;
1003: A[(j-1)*ld+j] = b[j-1];
1004: }
1005: DSRestoreArray(eps->ds,DS_MAT_A,&A);
1006: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
1007: DSSolve(eps->ds,eps->eigr,NULL);
1008: DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
1009: DSSetCompact(eps->ds,PETSC_TRUE);
1010: } else { /* Restart */
1011: DSSolve(eps->ds,eps->eigr,NULL);
1012: DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
1013: }
1014: DSSynchronize(eps->ds,eps->eigr,NULL);
1016: /* Residual */
1017: EPSKrylovConvergence(eps,PETSC_TRUE,eps->nconv,nv-eps->nconv,beta,0.0,1.0,&k);
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 = nv-k;
1054: if (breakdown) l=0;
1055: if (!ctx->lock && l>0 && eps->reason == EPS_CONVERGED_ITERATING) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
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: PetscSNPrintf(messg,sizeof(messg),"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])?"*":"");
1108: PetscInfo(eps,messg);
1109: if (count0 > sPres->nsch[0] || count1 > sPres->nsch[1]) InertiaMismatch(eps,ctx->detect);
1110: PetscFree(iwork);
1111: return(0);
1112: }
1114: /*
1115: Obtains value of subsequent shift
1116: */
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: *newS = (sPres->value + sPres->neighb[side]->value)/2;
1131: if (PetscAbsReal(sPres->value - *newS)/PetscAbsReal(sPres->value)<=100*PETSC_SQRT_MACHINE_EPSILON) *newS = sPres->value;
1132: } else { /* (Only for side=1). Creating a new interval. */
1133: if (sPres->neigs==0) {/* No value has been accepted*/
1134: if (sPres->neighb[0]) {
1135: /* Multiplying by 10 the previous distance */
1136: *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
1137: sr->nleap++;
1138: /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
1139: if (!sr->hasEnd && sr->nleap > 5) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to compute the wanted eigenvalues with open interval");
1140: } else { /* First shift */
1141: if (eps->nconv != 0) {
1142: /* Unaccepted values give information for next shift */
1143: idxP=0;/* Number of values left from shift */
1144: for (i=0;i<eps->nconv;i++) {
1145: lambda = PetscRealPart(eps->eigr[i]);
1146: if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
1147: else break;
1148: }
1149: /* Avoiding subtraction of eigenvalues (might be the same).*/
1150: if (idxP>0) {
1151: d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[0]))/(idxP+0.3);
1152: } else {
1153: d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[eps->nconv-1]))/(eps->nconv+0.3);
1154: }
1155: *newS = sPres->value + ((sr->dir)*d_prev*eps->nev)/2;
1156: } else { /* No values found, no information for next shift */
1157: SETERRQ(PetscObjectComm((PetscObject)eps),1,"First shift renders no information");
1158: }
1159: }
1160: } else { /* Accepted values found */
1161: sr->nleap = 0;
1162: /* Average distance of values in previous subinterval */
1163: s = sPres->neighb[0];
1164: while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
1165: s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
1166: }
1167: if (s) {
1168: d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
1169: } else { /* First shift. Average distance obtained with values in this shift */
1170: /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
1171: 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)) {
1172: d_prev = PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0])))/(sPres->neigs+0.3);
1173: } else {
1174: d_prev = PetscAbsReal(PetscRealPart(sr->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
1175: }
1176: }
1177: /* Average distance is used for next shift by adding it to value on the right or to shift */
1178: if ((sr->dir)*(PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
1179: *newS = PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(eps->nev))/2;
1180: } else { /* Last accepted value is on the left of shift. Adding to shift */
1181: *newS = sPres->value + ((sr->dir)*d_prev*(eps->nev))/2;
1182: }
1183: }
1184: /* End of interval can not be surpassed */
1185: if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
1186: }/* of neighb[side]==null */
1187: return(0);
1188: }
1190: /*
1191: Function for sorting an array of real values
1192: */
1193: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
1194: {
1195: PetscReal re;
1196: PetscInt i,j,tmp;
1199: if (!prev) for (i=0;i<nr;i++) perm[i] = i;
1200: /* Insertion sort */
1201: for (i=1;i<nr;i++) {
1202: re = PetscRealPart(r[perm[i]]);
1203: j = i-1;
1204: while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
1205: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
1206: }
1207: }
1208: return(0);
1209: }
1211: /* Stores the pairs obtained since the last shift in the global arrays */
1212: static PetscErrorCode EPSStoreEigenpairs(EPS eps)
1213: {
1214: PetscErrorCode ierr;
1215: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1216: PetscReal lambda,err,norm;
1217: PetscInt i,count;
1218: PetscBool iscayley;
1219: EPS_SR sr = ctx->sr;
1220: EPS_shift sPres;
1221: Vec v,w;
1224: sPres = sr->sPres;
1225: sPres->index = sr->indexEig;
1226: count = sr->indexEig;
1227: /* Back-transform */
1228: STBackTransform(eps->st,eps->nconv,eps->eigr,eps->eigi);
1229: PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
1230: /* Sort eigenvalues */
1231: sortRealEigenvalues(eps->eigr,eps->perm,eps->nconv,PETSC_FALSE,sr->dir);
1232: /* Values stored in global array */
1233: for (i=0;i<eps->nconv;i++) {
1234: lambda = PetscRealPart(eps->eigr[eps->perm[i]]);
1235: err = eps->errest[eps->perm[i]];
1237: if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
1238: if (count>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error in Spectrum Slicing");
1239: sr->eigr[count] = lambda;
1240: sr->errest[count] = err;
1241: /* Explicit purification */
1242: BVGetColumn(eps->V,eps->perm[i],&w);
1243: if (eps->purify) {
1244: BVGetColumn(sr->V,count,&v);
1245: STApply(eps->st,w,v);
1246: BVRestoreColumn(sr->V,count,&v);
1247: } else {
1248: BVInsertVec(sr->V,count,w);
1249: }
1250: BVRestoreColumn(eps->V,eps->perm[i],&w);
1251: BVNormColumn(sr->V,count,NORM_2,&norm);
1252: BVScaleColumn(sr->V,count,1.0/norm);
1253: count++;
1254: }
1255: }
1256: sPres->neigs = count - sr->indexEig;
1257: sr->indexEig = count;
1258: /* Global ordering array updating */
1259: sortRealEigenvalues(sr->eigr,sr->perm,count,PETSC_TRUE,sr->dir);
1260: return(0);
1261: }
1263: static PetscErrorCode EPSLookForDeflation(EPS eps)
1264: {
1265: PetscErrorCode ierr;
1266: PetscReal val;
1267: PetscInt i,count0=0,count1=0;
1268: EPS_shift sPres;
1269: PetscInt ini,fin,k,idx0,idx1;
1270: EPS_SR sr;
1271: Vec v;
1272: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1275: sr = ctx->sr;
1276: sPres = sr->sPres;
1278: if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
1279: else ini = 0;
1280: fin = sr->indexEig;
1281: /* Selection of ends for searching new values */
1282: if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
1283: else sPres->ext[0] = sPres->neighb[0]->value;
1284: if (!sPres->neighb[1]) {
1285: if (sr->hasEnd) sPres->ext[1] = sr->int1;
1286: else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
1287: } else sPres->ext[1] = sPres->neighb[1]->value;
1288: /* Selection of values between right and left ends */
1289: for (i=ini;i<fin;i++) {
1290: val=PetscRealPart(sr->eigr[sr->perm[i]]);
1291: /* Values to the right of left shift */
1292: if ((sr->dir)*(val - sPres->ext[1]) < 0) {
1293: if ((sr->dir)*(val - sPres->value) < 0) count0++;
1294: else count1++;
1295: } else break;
1296: }
1297: /* The number of values on each side are found */
1298: if (sPres->neighb[0]) {
1299: sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
1300: if (sPres->nsch[0]<0) InertiaMismatch(eps,ctx->detect);
1301: } else sPres->nsch[0] = 0;
1303: if (sPres->neighb[1]) {
1304: sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
1305: if (sPres->nsch[1]<0) InertiaMismatch(eps,ctx->detect);
1306: } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);
1308: /* Completing vector of indexes for deflation */
1309: idx0 = ini;
1310: idx1 = ini+count0+count1;
1311: k=0;
1312: for (i=idx0;i<idx1;i++) sr->idxDef[k++]=sr->perm[i];
1313: BVDuplicateResize(eps->V,k+eps->ncv+1,&sr->Vnext);
1314: BVSetNumConstraints(sr->Vnext,k);
1315: for (i=0;i<k;i++) {
1316: BVGetColumn(sr->Vnext,-i-1,&v);
1317: BVCopyVec(sr->V,sr->idxDef[i],v);
1318: BVRestoreColumn(sr->Vnext,-i-1,&v);
1319: }
1321: /* For rational Krylov */
1322: if (sr->nS>0 && (sr->sPrev == sr->sPres->neighb[0] || sr->sPrev == sr->sPres->neighb[1])) {
1323: EPSPrepareRational(eps);
1324: }
1325: eps->nconv = 0;
1326: /* Get rid of temporary Vnext */
1327: BVDestroy(&eps->V);
1328: eps->V = sr->Vnext;
1329: sr->Vnext = NULL;
1330: return(0);
1331: }
1333: PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS eps)
1334: {
1335: PetscErrorCode ierr;
1336: PetscInt i,lds,ti;
1337: PetscReal newS;
1338: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1339: EPS_SR sr=ctx->sr;
1340: Mat A,B=NULL;
1341: PetscObjectState Astate,Bstate=0;
1342: PetscObjectId Aid,Bid=0;
1345: PetscCitationsRegister(citation,&cited);
1346: if (ctx->global) {
1347: EPSSolve_KrylovSchur_Slice(ctx->eps);
1348: ctx->eps->state = EPS_STATE_SOLVED;
1349: eps->reason = EPS_CONVERGED_TOL;
1350: if (ctx->npart>1) {
1351: /* Gather solution from subsolvers */
1352: EPSSliceGatherSolution(eps);
1353: } else {
1354: eps->nconv = sr->numEigs;
1355: eps->its = ctx->eps->its;
1356: PetscFree(ctx->inertias);
1357: PetscFree(ctx->shifts);
1358: EPSSliceGetInertias(ctx->eps,&ctx->nshifts,&ctx->shifts,&ctx->inertias);
1359: }
1360: } else {
1361: if (ctx->npart==1) {
1362: sr->eigr = ctx->eps->eigr;
1363: sr->eigi = ctx->eps->eigi;
1364: sr->perm = ctx->eps->perm;
1365: sr->errest = ctx->eps->errest;
1366: sr->V = ctx->eps->V;
1367: }
1368: /* Check that the user did not modify subcomm matrices */
1369: EPSGetOperators(eps,&A,&B);
1370: PetscObjectStateGet((PetscObject)A,&Astate);
1371: PetscObjectGetId((PetscObject)A,&Aid);
1372: if (B) {
1373: PetscObjectStateGet((PetscObject)B,&Bstate);
1374: PetscObjectGetId((PetscObject)B,&Bid);
1375: }
1376: 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");
1377: /* Only with eigenvalues present in the interval ...*/
1378: if (sr->numEigs==0) {
1379: eps->reason = EPS_CONVERGED_TOL;
1380: return(0);
1381: }
1382: /* Array of pending shifts */
1383: sr->maxPend = 100; /* Initial size */
1384: sr->nPend = 0;
1385: PetscMalloc1(sr->maxPend,&sr->pending);
1386: PetscLogObjectMemory((PetscObject)eps,(sr->maxPend)*sizeof(EPS_shift));
1387: EPSCreateShift(eps,sr->int0,NULL,NULL);
1388: /* extract first shift */
1389: sr->sPrev = NULL;
1390: sr->sPres = sr->pending[--sr->nPend];
1391: sr->sPres->inertia = sr->inertia0;
1392: eps->target = sr->sPres->value;
1393: sr->s0 = sr->sPres;
1394: sr->indexEig = 0;
1395: /* Memory reservation for auxiliary variables */
1396: lds = PetscMin(eps->mpd,eps->ncv);
1397: PetscCalloc1(lds*lds,&sr->S);
1398: PetscMalloc1(eps->ncv,&sr->back);
1399: PetscLogObjectMemory((PetscObject)eps,(sr->numEigs+2*eps->ncv)*sizeof(PetscScalar));
1400: for (i=0;i<sr->numEigs;i++) {
1401: sr->eigr[i] = 0.0;
1402: sr->eigi[i] = 0.0;
1403: sr->errest[i] = 0.0;
1404: sr->perm[i] = i;
1405: }
1406: /* Vectors for deflation */
1407: PetscMalloc1(sr->numEigs,&sr->idxDef);
1408: PetscLogObjectMemory((PetscObject)eps,sr->numEigs*sizeof(PetscInt));
1409: sr->indexEig = 0;
1410: /* Main loop */
1411: while (sr->sPres) {
1412: /* Search for deflation */
1413: EPSLookForDeflation(eps);
1414: /* KrylovSchur */
1415: EPSKrylovSchur_Slice(eps);
1417: EPSStoreEigenpairs(eps);
1418: /* Select new shift */
1419: if (!sr->sPres->comp[1]) {
1420: EPSGetNewShiftValue(eps,1,&newS);
1421: EPSCreateShift(eps,newS,sr->sPres,sr->sPres->neighb[1]);
1422: }
1423: if (!sr->sPres->comp[0]) {
1424: /* Completing earlier interval */
1425: EPSGetNewShiftValue(eps,0,&newS);
1426: EPSCreateShift(eps,newS,sr->sPres->neighb[0],sr->sPres);
1427: }
1428: /* Preparing for a new search of values */
1429: EPSExtractShift(eps);
1430: }
1432: /* Updating eps values prior to exit */
1433: PetscFree(sr->S);
1434: PetscFree(sr->idxDef);
1435: PetscFree(sr->pending);
1436: PetscFree(sr->back);
1437: BVDuplicateResize(eps->V,eps->ncv+1,&sr->Vnext);
1438: BVSetNumConstraints(sr->Vnext,0);
1439: BVDestroy(&eps->V);
1440: eps->V = sr->Vnext;
1441: eps->nconv = sr->indexEig;
1442: eps->reason = EPS_CONVERGED_TOL;
1443: eps->its = sr->itsKs;
1444: eps->nds = 0;
1445: if (sr->dir<0) {
1446: for (i=0;i<eps->nconv/2;i++) {
1447: ti = sr->perm[i]; sr->perm[i] = sr->perm[eps->nconv-1-i]; sr->perm[eps->nconv-1-i] = ti;
1448: }
1449: }
1450: }
1451: return(0);
1452: }