Actual source code: ks-slice.c
slepc-3.16.0 2021-09-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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),PETSC_ERR_PLIB,"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: /* Reset auxiliary EPS */
85: EPSSliceResetSR(ctx->eps);
86: EPSReset(ctx->eps);
87: EPSSliceResetSR(eps);
88: PetscFree(ctx->inertias);
89: PetscFree(ctx->shifts);
90: return(0);
91: }
93: PetscErrorCode EPSDestroy_KrylovSchur_Slice(EPS eps)
94: {
95: PetscErrorCode ierr;
96: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
99: if (!ctx->global) return(0);
100: /* Destroy auxiliary EPS */
101: EPSReset_KrylovSchur_Slice(eps);
102: EPSDestroy(&ctx->eps);
103: if (ctx->npart>1) {
104: PetscSubcommDestroy(&ctx->subc);
105: if (ctx->commset) {
106: MPI_Comm_free(&ctx->commrank);
107: ctx->commset = PETSC_FALSE;
108: }
109: ISDestroy(&ctx->isrow);
110: ISDestroy(&ctx->iscol);
111: MatDestroyMatrices(1,&ctx->submata);
112: MatDestroyMatrices(1,&ctx->submatb);
113: }
114: PetscFree(ctx->subintervals);
115: PetscFree(ctx->nconv_loc);
116: return(0);
117: }
119: /*
120: EPSSliceAllocateSolution - Allocate memory storage for common variables such
121: as eigenvalues and eigenvectors. The argument extra is used for methods
122: that require a working basis slightly larger than ncv.
123: */
124: static PetscErrorCode EPSSliceAllocateSolution(EPS eps,PetscInt extra)
125: {
126: PetscErrorCode ierr;
127: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
128: PetscReal eta;
129: PetscInt k;
130: PetscLogDouble cnt;
131: BVType type;
132: BVOrthogType orthog_type;
133: BVOrthogRefineType orthog_ref;
134: BVOrthogBlockType ob_type;
135: Mat matrix;
136: Vec t;
137: EPS_SR sr = ctx->sr;
140: /* allocate space for eigenvalues and friends */
141: k = PetscMax(1,sr->numEigs);
142: PetscFree4(sr->eigr,sr->eigi,sr->errest,sr->perm);
143: PetscMalloc4(k,&sr->eigr,k,&sr->eigi,k,&sr->errest,k,&sr->perm);
144: cnt = 2*k*sizeof(PetscScalar) + 2*k*sizeof(PetscReal) + k*sizeof(PetscInt);
145: PetscLogObjectMemory((PetscObject)eps,cnt);
147: /* allocate sr->V and transfer options from eps->V */
148: BVDestroy(&sr->V);
149: BVCreate(PetscObjectComm((PetscObject)eps),&sr->V);
150: PetscLogObjectParent((PetscObject)eps,(PetscObject)sr->V);
151: if (!eps->V) { EPSGetBV(eps,&eps->V); }
152: if (!((PetscObject)(eps->V))->type_name) {
153: BVSetType(sr->V,BVSVEC);
154: } else {
155: BVGetType(eps->V,&type);
156: BVSetType(sr->V,type);
157: }
158: STMatCreateVecsEmpty(eps->st,&t,NULL);
159: BVSetSizesFromVec(sr->V,t,k);
160: VecDestroy(&t);
161: EPS_SetInnerProduct(eps);
162: BVGetMatrix(eps->V,&matrix,NULL);
163: BVSetMatrix(sr->V,matrix,PETSC_FALSE);
164: BVGetOrthogonalization(eps->V,&orthog_type,&orthog_ref,&eta,&ob_type);
165: BVSetOrthogonalization(sr->V,orthog_type,orthog_ref,eta,ob_type);
166: return(0);
167: }
169: static PetscErrorCode EPSSliceGetEPS(EPS eps)
170: {
171: PetscErrorCode ierr;
172: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data,*ctx_local;
173: BV V;
174: BVType type;
175: PetscReal eta;
176: BVOrthogType orthog_type;
177: BVOrthogRefineType orthog_ref;
178: BVOrthogBlockType ob_type;
179: PetscInt i;
180: PetscReal h,a,b;
181: PetscRandom rand;
182: EPS_SR sr=ctx->sr;
185: if (!ctx->eps) { EPSKrylovSchurGetChildEPS(eps,&ctx->eps); }
187: /* Determine subintervals */
188: if (ctx->npart==1) {
189: a = eps->inta; b = eps->intb;
190: } else {
191: if (!ctx->subintset) { /* uniform distribution if no set by user */
192: if (!sr->hasEnd) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Global interval must be bounded for splitting it in uniform subintervals");
193: h = (eps->intb-eps->inta)/ctx->npart;
194: a = eps->inta+ctx->subc->color*h;
195: b = (ctx->subc->color==ctx->npart-1)?eps->intb:eps->inta+(ctx->subc->color+1)*h;
196: PetscFree(ctx->subintervals);
197: PetscMalloc1(ctx->npart+1,&ctx->subintervals);
198: for (i=0;i<ctx->npart;i++) ctx->subintervals[i] = eps->inta+h*i;
199: ctx->subintervals[ctx->npart] = eps->intb;
200: } else {
201: a = ctx->subintervals[ctx->subc->color];
202: b = ctx->subintervals[ctx->subc->color+1];
203: }
204: }
205: EPSSetInterval(ctx->eps,a,b);
206: EPSSetConvergenceTest(ctx->eps,eps->conv);
207: EPSSetDimensions(ctx->eps,ctx->nev,ctx->ncv,ctx->mpd);
208: EPSKrylovSchurSetLocking(ctx->eps,ctx->lock);
210: ctx_local = (EPS_KRYLOVSCHUR*)ctx->eps->data;
211: ctx_local->detect = ctx->detect;
213: /* transfer options from eps->V */
214: EPSGetBV(ctx->eps,&V);
215: BVGetRandomContext(V,&rand); /* make sure the random context is available when duplicating */
216: if (!eps->V) { EPSGetBV(eps,&eps->V); }
217: if (!((PetscObject)(eps->V))->type_name) {
218: BVSetType(V,BVSVEC);
219: } else {
220: BVGetType(eps->V,&type);
221: BVSetType(V,type);
222: }
223: BVGetOrthogonalization(eps->V,&orthog_type,&orthog_ref,&eta,&ob_type);
224: BVSetOrthogonalization(V,orthog_type,orthog_ref,eta,ob_type);
226: ctx->eps->which = eps->which;
227: ctx->eps->max_it = eps->max_it;
228: ctx->eps->tol = eps->tol;
229: ctx->eps->purify = eps->purify;
230: if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
231: EPSSetProblemType(ctx->eps,eps->problem_type);
232: EPSSetUp(ctx->eps);
233: ctx->eps->nconv = 0;
234: ctx->eps->its = 0;
235: for (i=0;i<ctx->eps->ncv;i++) {
236: ctx->eps->eigr[i] = 0.0;
237: ctx->eps->eigi[i] = 0.0;
238: ctx->eps->errest[i] = 0.0;
239: }
240: return(0);
241: }
243: static PetscErrorCode EPSSliceGetInertia(EPS eps,PetscReal shift,PetscInt *inertia,PetscInt *zeros)
244: {
246: KSP ksp,kspr;
247: PC pc;
248: Mat F;
249: PetscReal nzshift=shift;
250: PetscBool flg;
253: if (shift >= PETSC_MAX_REAL) { /* Right-open interval */
254: if (inertia) *inertia = eps->n;
255: } else if (shift <= PETSC_MIN_REAL) {
256: if (inertia) *inertia = 0;
257: if (zeros) *zeros = 0;
258: } else {
259: /* If the shift is zero, perturb it to a very small positive value.
260: The goal is that the nonzero pattern is the same in all cases and reuse
261: the symbolic factorizations */
262: nzshift = (shift==0.0)? 10.0/PETSC_MAX_REAL: shift;
263: STSetShift(eps->st,nzshift);
264: STGetKSP(eps->st,&ksp);
265: KSPGetPC(ksp,&pc);
266: PetscObjectTypeCompare((PetscObject)pc,PCREDUNDANT,&flg);
267: if (flg) {
268: PCRedundantGetKSP(pc,&kspr);
269: KSPGetPC(kspr,&pc);
270: }
271: PCFactorGetMatrix(pc,&F);
272: MatGetInertia(F,inertia,zeros,NULL);
273: }
274: if (inertia) { PetscInfo2(eps,"Computed inertia at shift %g: %D\n",(double)nzshift,*inertia); }
275: return(0);
276: }
278: /*
279: Dummy backtransform operation
280: */
281: static PetscErrorCode EPSBackTransform_Skip(EPS eps)
282: {
284: return(0);
285: }
287: PetscErrorCode EPSSetUp_KrylovSchur_Slice(EPS eps)
288: {
289: PetscErrorCode ierr;
290: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data,*ctx_glob;
291: EPS_SR sr,sr_loc,sr_glob;
292: PetscInt nEigs,dssz=1,i,zeros=0,off=0,method,hiteig=0;
293: PetscMPIInt nproc,rank=0,aux;
294: PetscReal r;
295: MPI_Request req;
296: Mat A,B=NULL;
297: DSParallelType ptype;
300: if (ctx->global) {
301: EPSCheckHermitianDefiniteCondition(eps,PETSC_TRUE," with spectrum slicing");
302: EPSCheckSinvertCayleyCondition(eps,PETSC_TRUE," with spectrum slicing");
303: 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()");
304: 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");
305: if (eps->nds) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing not supported in combination with deflation space");
306: EPSCheckUnsupportedCondition(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING,PETSC_TRUE," with spectrum slicing");
307: EPSCheckIgnoredCondition(eps,EPS_FEATURE_BALANCE,PETSC_TRUE," with spectrum slicing");
308: if (eps->tol==PETSC_DEFAULT) {
309: #if defined(PETSC_USE_REAL_SINGLE)
310: eps->tol = SLEPC_DEFAULT_TOL;
311: #else
312: /* use tighter tolerance */
313: eps->tol = SLEPC_DEFAULT_TOL*1e-2;
314: #endif
315: }
316: if (eps->max_it==PETSC_DEFAULT) eps->max_it = 100;
317: if (ctx->nev==1) ctx->nev = PetscMin(40,eps->n); /* nev not set, use default value */
318: 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");
319: }
320: eps->ops->backtransform = EPSBackTransform_Skip;
322: /* create spectrum slicing context and initialize it */
323: EPSSliceResetSR(eps);
324: PetscNewLog(eps,&sr);
325: ctx->sr = sr;
326: sr->itsKs = 0;
327: sr->nleap = 0;
328: sr->nMAXCompl = eps->nev/4;
329: sr->iterCompl = eps->max_it/4;
330: sr->sPres = NULL;
331: sr->nS = 0;
333: if (ctx->npart==1 || ctx->global) {
334: /* check presence of ends and finding direction */
335: if ((eps->inta > PETSC_MIN_REAL && !(ctx->subintervals && ctx->subintervals[0]==ctx->subintervals[1])) || eps->intb >= PETSC_MAX_REAL) {
336: sr->int0 = eps->inta;
337: sr->int1 = eps->intb;
338: sr->dir = 1;
339: if (eps->intb >= PETSC_MAX_REAL) { /* Right-open interval */
340: sr->hasEnd = PETSC_FALSE;
341: } else sr->hasEnd = PETSC_TRUE;
342: } else {
343: sr->int0 = eps->intb;
344: sr->int1 = eps->inta;
345: sr->dir = -1;
346: sr->hasEnd = PetscNot(eps->inta <= PETSC_MIN_REAL);
347: }
348: }
349: if (ctx->global) {
350: EPSSetDimensions_Default(eps,ctx->nev,&ctx->ncv,&ctx->mpd);
351: /* create subintervals and initialize auxiliary eps for slicing runs */
352: EPSKrylovSchurGetChildEPS(eps,&ctx->eps);
353: /* prevent computation of factorization in global eps */
354: STSetTransform(eps->st,PETSC_FALSE);
355: EPSSliceGetEPS(eps);
356: sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
357: if (ctx->npart>1) {
358: if ((sr->dir>0&&ctx->subc->color==0)||(sr->dir<0&&ctx->subc->color==ctx->npart-1)) sr->inertia0 = sr_loc->inertia0;
359: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
360: if (!rank) {
361: MPI_Bcast(&sr->inertia0,1,MPIU_INT,(sr->dir>0)?0:ctx->npart-1,ctx->commrank);
362: }
363: MPI_Bcast(&sr->inertia0,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
364: PetscFree(ctx->nconv_loc);
365: PetscMalloc1(ctx->npart,&ctx->nconv_loc);
366: MPI_Comm_size(((PetscObject)eps)->comm,&nproc);
367: if (sr->dir<0) off = 1;
368: if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
369: PetscMPIIntCast(sr_loc->numEigs,&aux);
370: MPI_Allgather(&aux,1,MPI_INT,ctx->nconv_loc,1,MPI_INT,ctx->commrank);
371: MPI_Allgather(sr_loc->dir==sr->dir?&sr_loc->int0:&sr_loc->int1,1,MPIU_REAL,ctx->subintervals+off,1,MPIU_REAL,ctx->commrank);
372: } else {
373: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
374: if (!rank) {
375: PetscMPIIntCast(sr_loc->numEigs,&aux);
376: MPI_Allgather(&aux,1,MPI_INT,ctx->nconv_loc,1,MPI_INT,ctx->commrank);
377: MPI_Allgather(sr_loc->dir==sr->dir?&sr_loc->int0:&sr_loc->int1,1,MPIU_REAL,ctx->subintervals+off,1,MPIU_REAL,ctx->commrank);
378: }
379: PetscMPIIntCast(ctx->npart,&aux);
380: MPI_Bcast(ctx->nconv_loc,aux,MPI_INT,0,PetscSubcommChild(ctx->subc));
381: MPI_Bcast(ctx->subintervals+off,aux,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
382: }
383: nEigs = 0;
384: for (i=0;i<ctx->npart;i++) nEigs += ctx->nconv_loc[i];
385: } else {
386: nEigs = sr_loc->numEigs;
387: sr->inertia0 = sr_loc->inertia0;
388: sr->dir = sr_loc->dir;
389: }
390: sr->inertia1 = sr->inertia0+sr->dir*nEigs;
391: sr->numEigs = nEigs;
392: eps->nev = nEigs;
393: eps->ncv = nEigs;
394: eps->mpd = nEigs;
395: } else {
396: ctx_glob = (EPS_KRYLOVSCHUR*)ctx->eps->data;
397: sr_glob = ctx_glob->sr;
398: if (ctx->npart>1) {
399: sr->dir = sr_glob->dir;
400: sr->int0 = (sr->dir==1)?eps->inta:eps->intb;
401: sr->int1 = (sr->dir==1)?eps->intb:eps->inta;
402: if ((sr->dir>0&&ctx->subc->color==ctx->npart-1)||(sr->dir<0&&ctx->subc->color==0)) sr->hasEnd = sr_glob->hasEnd;
403: else sr->hasEnd = PETSC_TRUE;
404: }
405: /* sets first shift */
406: STSetShift(eps->st,(sr->int0==0.0)?10.0/PETSC_MAX_REAL:sr->int0);
407: STSetUp(eps->st);
409: /* compute inertia0 */
410: EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,ctx->detect?&zeros:NULL);
411: /* undocumented option to control what to do when an eigenvalue is found:
412: - error out if it's the endpoint of the user-provided interval (or sub-interval)
413: - if it's an endpoint computed internally:
414: + if hiteig=0 error out
415: + else if hiteig=1 the subgroup that hit the eigenvalue does nothing
416: + otherwise the subgroup that hit the eigenvalue perturbs the shift and recomputes inertia
417: */
418: PetscOptionsGetInt(NULL,NULL,"-eps_krylovschur_hiteigenvalue",&hiteig,NULL);
419: if (zeros) { /* error in factorization */
420: 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");
421: 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");
422: else {
423: if (hiteig==1) { /* idle subgroup */
424: sr->inertia0 = -1;
425: } else { /* perturb shift */
426: sr->int0 *= (1.0+SLICE_PTOL);
427: EPSSliceGetInertia(eps,sr->int0,&sr->inertia0,&zeros);
428: if (zeros) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",sr->int1);
429: }
430: }
431: }
432: if (ctx->npart>1) {
433: /* inertia1 is received from neighbour */
434: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
435: if (!rank) {
436: if (sr->inertia0!=-1 && ((sr->dir>0 && ctx->subc->color>0) || (sr->dir<0 && ctx->subc->color<ctx->npart-1))) { /* send inertia0 to neighbour0 */
437: MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
438: MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
439: }
440: if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)|| (sr->dir<0 && ctx->subc->color>0)) { /* receive inertia1 from neighbour1 */
441: MPI_Recv(&(sr->inertia1),1,MPIU_INT,ctx->subc->color+sr->dir,0,ctx->commrank,MPI_STATUS_IGNORE);
442: MPI_Recv(&(sr->int1),1,MPIU_REAL,ctx->subc->color+sr->dir,0,ctx->commrank,MPI_STATUS_IGNORE);
443: }
444: if (sr->inertia0==-1 && !(sr->dir>0 && ctx->subc->color==ctx->npart-1) && !(sr->dir<0 && ctx->subc->color==0)) {
445: sr->inertia0 = sr->inertia1; sr->int0 = sr->int1;
446: MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
447: MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
448: }
449: }
450: if ((sr->dir>0 && ctx->subc->color<ctx->npart-1)||(sr->dir<0 && ctx->subc->color>0)) {
451: MPI_Bcast(&sr->inertia1,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
452: MPI_Bcast(&sr->int1,1,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
453: } else sr_glob->inertia1 = sr->inertia1;
454: }
456: /* last process in eps comm computes inertia1 */
457: if (ctx->npart==1 || ((sr->dir>0 && ctx->subc->color==ctx->npart-1) || (sr->dir<0 && ctx->subc->color==0))) {
458: EPSSliceGetInertia(eps,sr->int1,&sr->inertia1,ctx->detect?&zeros:NULL);
459: if (zeros) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
460: if (!rank && sr->inertia0==-1) {
461: sr->inertia0 = sr->inertia1; sr->int0 = sr->int1;
462: MPI_Isend(&(sr->inertia0),1,MPIU_INT,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
463: MPI_Isend(&(sr->int0),1,MPIU_REAL,ctx->subc->color-sr->dir,0,ctx->commrank,&req);
464: }
465: if (sr->hasEnd) {
466: sr->dir = -sr->dir; r = sr->int0; sr->int0 = sr->int1; sr->int1 = r;
467: i = sr->inertia0; sr->inertia0 = sr->inertia1; sr->inertia1 = i;
468: }
469: }
471: /* number of eigenvalues in interval */
472: sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
473: if (ctx->npart>1) {
474: /* memory allocate for subinterval eigenpairs */
475: EPSSliceAllocateSolution(eps,1);
476: }
477: dssz = eps->ncv+1;
478: DSGetParallel(ctx->eps->ds,&ptype);
479: DSSetParallel(eps->ds,ptype);
480: DSGetMethod(ctx->eps->ds,&method);
481: DSSetMethod(eps->ds,method);
482: }
483: DSSetType(eps->ds,DSHEP);
484: DSSetCompact(eps->ds,PETSC_TRUE);
485: DSAllocate(eps->ds,dssz);
486: /* keep state of subcomm matrices to check that the user does not modify them */
487: EPSGetOperators(eps,&A,&B);
488: PetscObjectStateGet((PetscObject)A,&ctx->Astate);
489: PetscObjectGetId((PetscObject)A,&ctx->Aid);
490: if (B) {
491: PetscObjectStateGet((PetscObject)B,&ctx->Bstate);
492: PetscObjectGetId((PetscObject)B,&ctx->Bid);
493: } else {
494: ctx->Bstate=0;
495: ctx->Bid=0;
496: }
497: return(0);
498: }
500: static PetscErrorCode EPSSliceGatherEigenVectors(EPS eps)
501: {
502: PetscErrorCode ierr;
503: Vec v,vg,v_loc;
504: IS is1,is2;
505: VecScatter vec_sc;
506: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
507: PetscInt nloc,m0,n0,i,si,idx,*idx1,*idx2,j;
508: PetscScalar *array;
509: EPS_SR sr_loc;
510: BV V_loc;
513: sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
514: V_loc = sr_loc->V;
516: /* Gather parallel eigenvectors */
517: BVGetColumn(eps->V,0,&v);
518: VecGetOwnershipRange(v,&n0,&m0);
519: BVRestoreColumn(eps->V,0,&v);
520: BVGetColumn(ctx->eps->V,0,&v);
521: VecGetLocalSize(v,&nloc);
522: BVRestoreColumn(ctx->eps->V,0,&v);
523: PetscMalloc2(m0-n0,&idx1,m0-n0,&idx2);
524: VecCreateMPI(PetscObjectComm((PetscObject)eps),nloc,PETSC_DECIDE,&vg);
525: idx = -1;
526: for (si=0;si<ctx->npart;si++) {
527: j = 0;
528: for (i=n0;i<m0;i++) {
529: idx1[j] = i;
530: idx2[j++] = i+eps->n*si;
531: }
532: ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx1,PETSC_COPY_VALUES,&is1);
533: ISCreateGeneral(PetscObjectComm((PetscObject)eps),(m0-n0),idx2,PETSC_COPY_VALUES,&is2);
534: BVGetColumn(eps->V,0,&v);
535: VecScatterCreate(v,is1,vg,is2,&vec_sc);
536: BVRestoreColumn(eps->V,0,&v);
537: ISDestroy(&is1);
538: ISDestroy(&is2);
539: for (i=0;i<ctx->nconv_loc[si];i++) {
540: BVGetColumn(eps->V,++idx,&v);
541: if (ctx->subc->color==si) {
542: BVGetColumn(V_loc,i,&v_loc);
543: VecGetArray(v_loc,&array);
544: VecPlaceArray(vg,array);
545: }
546: VecScatterBegin(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE);
547: VecScatterEnd(vec_sc,vg,v,INSERT_VALUES,SCATTER_REVERSE);
548: if (ctx->subc->color==si) {
549: VecResetArray(vg);
550: VecRestoreArray(v_loc,&array);
551: BVRestoreColumn(V_loc,i,&v_loc);
552: }
553: BVRestoreColumn(eps->V,idx,&v);
554: }
555: VecScatterDestroy(&vec_sc);
556: }
557: PetscFree2(idx1,idx2);
558: VecDestroy(&vg);
559: return(0);
560: }
562: /*
563: EPSComputeVectors_Slice - Recover Eigenvectors from subcomunicators
564: */
565: PetscErrorCode EPSComputeVectors_Slice(EPS eps)
566: {
567: PetscErrorCode ierr;
568: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
571: if (ctx->global && ctx->npart>1) {
572: EPSComputeVectors(ctx->eps);
573: EPSSliceGatherEigenVectors(eps);
574: }
575: return(0);
576: }
578: #define SWAP(a,b,t) {t=a;a=b;b=t;}
580: static PetscErrorCode EPSSliceGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
581: {
582: PetscErrorCode ierr;
583: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
584: PetscInt i=0,j,tmpi;
585: PetscReal v,tmpr;
586: EPS_shift s;
589: if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
590: if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
591: if (!ctx->sr->s0) { /* EPSSolve not called yet */
592: *n = 2;
593: } else {
594: *n = 1;
595: s = ctx->sr->s0;
596: while (s) {
597: (*n)++;
598: s = s->neighb[1];
599: }
600: }
601: PetscMalloc1(*n,shifts);
602: PetscMalloc1(*n,inertias);
603: if (!ctx->sr->s0) { /* EPSSolve not called yet */
604: (*shifts)[0] = ctx->sr->int0;
605: (*shifts)[1] = ctx->sr->int1;
606: (*inertias)[0] = ctx->sr->inertia0;
607: (*inertias)[1] = ctx->sr->inertia1;
608: } else {
609: s = ctx->sr->s0;
610: while (s) {
611: (*shifts)[i] = s->value;
612: (*inertias)[i++] = s->inertia;
613: s = s->neighb[1];
614: }
615: (*shifts)[i] = ctx->sr->int1;
616: (*inertias)[i] = ctx->sr->inertia1;
617: }
618: /* remove possible duplicate in last position */
619: if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
620: /* sort result */
621: for (i=0;i<*n;i++) {
622: v = (*shifts)[i];
623: for (j=i+1;j<*n;j++) {
624: if (v > (*shifts)[j]) {
625: SWAP((*shifts)[i],(*shifts)[j],tmpr);
626: SWAP((*inertias)[i],(*inertias)[j],tmpi);
627: v = (*shifts)[i];
628: }
629: }
630: }
631: return(0);
632: }
634: static PetscErrorCode EPSSliceGatherSolution(EPS eps)
635: {
636: PetscErrorCode ierr;
637: PetscMPIInt rank,nproc;
638: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
639: PetscInt i,idx,j;
640: PetscInt *perm_loc,off=0,*inertias_loc,ns;
641: PetscScalar *eigr_loc;
642: EPS_SR sr_loc;
643: PetscReal *shifts_loc;
644: PetscMPIInt *disp,*ns_loc,aux;
647: eps->nconv = 0;
648: for (i=0;i<ctx->npart;i++) eps->nconv += ctx->nconv_loc[i];
649: sr_loc = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;
651: /* Gather the shifts used and the inertias computed */
652: EPSSliceGetInertias(ctx->eps,&ns,&shifts_loc,&inertias_loc);
653: if (ctx->sr->dir>0 && shifts_loc[ns-1]==sr_loc->int1 && ctx->subc->color<ctx->npart-1) ns--;
654: if (ctx->sr->dir<0 && shifts_loc[ns-1]==sr_loc->int0 && ctx->subc->color>0) {
655: ns--;
656: for (i=0;i<ns;i++) {
657: inertias_loc[i] = inertias_loc[i+1];
658: shifts_loc[i] = shifts_loc[i+1];
659: }
660: }
661: PetscMalloc1(ctx->npart,&ns_loc);
662: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
663: PetscMPIIntCast(ns,&aux);
664: if (!rank) { MPI_Allgather(&aux,1,MPI_INT,ns_loc,1,MPI_INT,ctx->commrank); }
665: PetscMPIIntCast(ctx->npart,&aux);
666: MPI_Bcast(ns_loc,aux,MPI_INT,0,PetscSubcommChild(ctx->subc));
667: ctx->nshifts = 0;
668: for (i=0;i<ctx->npart;i++) ctx->nshifts += ns_loc[i];
669: PetscFree(ctx->inertias);
670: PetscFree(ctx->shifts);
671: PetscMalloc1(ctx->nshifts,&ctx->inertias);
672: PetscMalloc1(ctx->nshifts,&ctx->shifts);
674: /* Gather eigenvalues (same ranks have fully set of eigenvalues)*/
675: eigr_loc = sr_loc->eigr;
676: perm_loc = sr_loc->perm;
677: MPI_Comm_size(((PetscObject)eps)->comm,&nproc);
678: PetscMalloc1(ctx->npart,&disp);
679: disp[0] = 0;
680: for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ctx->nconv_loc[i-1];
681: if (nproc%ctx->npart==0) { /* subcommunicators with the same size */
682: PetscMPIIntCast(sr_loc->numEigs,&aux);
683: MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank); /* eigenvalues */
684: MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank); /* perm */
685: for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
686: PetscMPIIntCast(ns,&aux);
687: MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank); /* shifts */
688: MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank); /* inertias */
689: MPIU_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank);
690: } else { /* subcommunicators with different size */
691: MPI_Comm_rank(PetscSubcommChild(ctx->subc),&rank);
692: if (!rank) {
693: PetscMPIIntCast(sr_loc->numEigs,&aux);
694: MPI_Allgatherv(eigr_loc,aux,MPIU_SCALAR,eps->eigr,ctx->nconv_loc,disp,MPIU_SCALAR,ctx->commrank); /* eigenvalues */
695: MPI_Allgatherv(perm_loc,aux,MPIU_INT,eps->perm,ctx->nconv_loc,disp,MPIU_INT,ctx->commrank); /* perm */
696: for (i=1;i<ctx->npart;i++) disp[i] = disp[i-1]+ns_loc[i-1];
697: PetscMPIIntCast(ns,&aux);
698: MPI_Allgatherv(shifts_loc,aux,MPIU_REAL,ctx->shifts,ns_loc,disp,MPIU_REAL,ctx->commrank); /* shifts */
699: MPI_Allgatherv(inertias_loc,aux,MPIU_INT,ctx->inertias,ns_loc,disp,MPIU_INT,ctx->commrank); /* inertias */
700: MPIU_Allreduce(&sr_loc->itsKs,&eps->its,1,MPIU_INT,MPI_SUM,ctx->commrank);
701: }
702: PetscMPIIntCast(eps->nconv,&aux);
703: MPI_Bcast(eps->eigr,aux,MPIU_SCALAR,0,PetscSubcommChild(ctx->subc));
704: MPI_Bcast(eps->perm,aux,MPIU_INT,0,PetscSubcommChild(ctx->subc));
705: MPI_Bcast(ctx->shifts,ctx->nshifts,MPIU_REAL,0,PetscSubcommChild(ctx->subc));
706: PetscMPIIntCast(ctx->nshifts,&aux);
707: MPI_Bcast(ctx->inertias,aux,MPIU_INT,0,PetscSubcommChild(ctx->subc));
708: MPI_Bcast(&eps->its,1,MPIU_INT,0,PetscSubcommChild(ctx->subc));
709: }
710: /* Update global array eps->perm */
711: idx = ctx->nconv_loc[0];
712: for (i=1;i<ctx->npart;i++) {
713: off += ctx->nconv_loc[i-1];
714: for (j=0;j<ctx->nconv_loc[i];j++) eps->perm[idx++] += off;
715: }
717: /* Gather parallel eigenvectors */
718: PetscFree(ns_loc);
719: PetscFree(disp);
720: PetscFree(shifts_loc);
721: PetscFree(inertias_loc);
722: return(0);
723: }
725: /*
726: Fills the fields of a shift structure
727: */
728: static PetscErrorCode EPSCreateShift(EPS eps,PetscReal val,EPS_shift neighb0,EPS_shift neighb1)
729: {
730: PetscErrorCode ierr;
731: EPS_shift s,*pending2;
732: PetscInt i;
733: EPS_SR sr;
734: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
737: sr = ctx->sr;
738: if ((neighb0 && val==neighb0->value) || (neighb1 && val==neighb1->value)) {
739: sr->nPend++;
740: return(0);
741: }
742: PetscNewLog(eps,&s);
743: s->value = val;
744: s->neighb[0] = neighb0;
745: if (neighb0) neighb0->neighb[1] = s;
746: s->neighb[1] = neighb1;
747: if (neighb1) neighb1->neighb[0] = s;
748: s->comp[0] = PETSC_FALSE;
749: s->comp[1] = PETSC_FALSE;
750: s->index = -1;
751: s->neigs = 0;
752: s->nconv[0] = s->nconv[1] = 0;
753: s->nsch[0] = s->nsch[1]=0;
754: /* Inserts in the stack of pending shifts */
755: /* If needed, the array is resized */
756: if (sr->nPend >= sr->maxPend) {
757: sr->maxPend *= 2;
758: PetscMalloc1(sr->maxPend,&pending2);
759: PetscLogObjectMemory((PetscObject)eps,sr->maxPend*sizeof(EPS_shift*));
760: for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
761: PetscFree(sr->pending);
762: sr->pending = pending2;
763: }
764: sr->pending[sr->nPend++]=s;
765: return(0);
766: }
768: /* Prepare for Rational Krylov update */
769: static PetscErrorCode EPSPrepareRational(EPS eps)
770: {
771: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
772: PetscErrorCode ierr;
773: PetscInt dir,i,k,ld,nv;
774: PetscScalar *A;
775: EPS_SR sr = ctx->sr;
776: Vec v;
779: DSGetLeadingDimension(eps->ds,&ld);
780: dir = (sr->sPres->neighb[0] == sr->sPrev)?1:-1;
781: dir*=sr->dir;
782: k = 0;
783: for (i=0;i<sr->nS;i++) {
784: if (dir*PetscRealPart(sr->S[i])>0.0) {
785: sr->S[k] = sr->S[i];
786: sr->S[sr->nS+k] = sr->S[sr->nS+i];
787: BVGetColumn(sr->Vnext,k,&v);
788: BVCopyVec(eps->V,eps->nconv+i,v);
789: BVRestoreColumn(sr->Vnext,k,&v);
790: k++;
791: if (k>=sr->nS/2)break;
792: }
793: }
794: /* Copy to DS */
795: DSGetArray(eps->ds,DS_MAT_A,&A);
796: PetscArrayzero(A,ld*ld);
797: for (i=0;i<k;i++) {
798: A[i*(1+ld)] = sr->S[i];
799: A[k+i*ld] = sr->S[sr->nS+i];
800: }
801: sr->nS = k;
802: DSRestoreArray(eps->ds,DS_MAT_A,&A);
803: DSGetDimensions(eps->ds,&nv,NULL,NULL,NULL);
804: DSSetDimensions(eps->ds,nv,0,k);
805: /* Append u to V */
806: BVGetColumn(sr->Vnext,sr->nS,&v);
807: BVCopyVec(eps->V,sr->nv,v);
808: BVRestoreColumn(sr->Vnext,sr->nS,&v);
809: return(0);
810: }
812: /* Provides next shift to be computed */
813: static PetscErrorCode EPSExtractShift(EPS eps)
814: {
815: PetscErrorCode ierr;
816: PetscInt iner,zeros=0;
817: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
818: EPS_SR sr;
819: PetscReal newShift,diam,ptol;
820: EPS_shift sPres;
823: sr = ctx->sr;
824: if (sr->nPend > 0) {
825: if (sr->sPres==sr->pending[sr->nPend-1]) {
826: eps->reason = EPS_CONVERGED_ITERATING;
827: eps->its = 0;
828: sr->nPend--;
829: sr->sPres->rep = PETSC_TRUE;
830: return(0);
831: }
832: sr->sPrev = sr->sPres;
833: sr->sPres = sr->pending[--sr->nPend];
834: sPres = sr->sPres;
835: EPSSliceGetInertia(eps,sPres->value,&iner,ctx->detect?&zeros:NULL);
836: if (zeros) {
837: diam = PetscMin(PetscAbsReal(sPres->neighb[0]->value-sPres->value),PetscAbsReal(sPres->neighb[1]->value-sPres->value));
838: ptol = PetscMin(SLICE_PTOL,diam/2);
839: newShift = sPres->value*(1.0+ptol);
840: if (sr->dir*(sPres->neighb[0] && newShift-sPres->neighb[0]->value) < 0) newShift = (sPres->value+sPres->neighb[0]->value)/2;
841: else if (sPres->neighb[1] && sr->dir*(sPres->neighb[1]->value-newShift) < 0) newShift = (sPres->value+sPres->neighb[1]->value)/2;
842: EPSSliceGetInertia(eps,newShift,&iner,&zeros);
843: if (zeros) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",newShift);
844: sPres->value = newShift;
845: }
846: sr->sPres->inertia = iner;
847: eps->target = sr->sPres->value;
848: eps->reason = EPS_CONVERGED_ITERATING;
849: eps->its = 0;
850: } else sr->sPres = NULL;
851: return(0);
852: }
854: /*
855: Symmetric KrylovSchur adapted to spectrum slicing:
856: Allows searching an specific amount of eigenvalues in the subintervals left and right.
857: Returns whether the search has succeeded
858: */
859: static PetscErrorCode EPSKrylovSchur_Slice(EPS eps)
860: {
861: PetscErrorCode ierr;
862: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
863: PetscInt i,k,l,ld,nv,*iwork,j,count0,count1,iterCompl=0,n0,n1;
864: Mat U,Op;
865: PetscScalar *Q,*A;
866: PetscReal *a,*b,beta,lambda;
867: EPS_shift sPres;
868: PetscBool breakdown,complIterating,sch0,sch1;
869: EPS_SR sr = ctx->sr;
870: char messg[100];
873: /* Spectrum slicing data */
874: sPres = sr->sPres;
875: complIterating =PETSC_FALSE;
876: sch1 = sch0 = PETSC_TRUE;
877: DSGetLeadingDimension(eps->ds,&ld);
878: PetscMalloc1(2*ld,&iwork);
879: count0=0;count1=0; /* Found on both sides */
880: if (!sPres->rep && sr->nS > 0 && (sPres->neighb[0] == sr->sPrev || sPres->neighb[1] == sr->sPrev)) {
881: /* Rational Krylov */
882: DSTranslateRKS(eps->ds,sr->sPrev->value-sPres->value);
883: DSGetDimensions(eps->ds,NULL,NULL,&l,NULL);
884: DSSetDimensions(eps->ds,l+1,0,0);
885: BVSetActiveColumns(eps->V,0,l+1);
886: DSGetMat(eps->ds,DS_MAT_Q,&U);
887: BVMultInPlace(eps->V,U,0,l+1);
888: MatDestroy(&U);
889: } else {
890: /* Get the starting Lanczos vector */
891: EPSGetStartVector(eps,0,NULL);
892: l = 0;
893: }
894: /* Restart loop */
895: while (eps->reason == EPS_CONVERGED_ITERATING) {
896: eps->its++; sr->itsKs++;
897: /* Compute an nv-step Lanczos factorization */
898: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
899: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
900: b = a + ld;
901: STGetOperator(eps->st,&Op);
902: BVMatLanczos(eps->V,Op,a,b,eps->nconv+l,&nv,&breakdown);
903: STRestoreOperator(eps->st,&Op);
904: sr->nv = nv;
905: beta = b[nv-1];
906: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
907: DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l);
908: if (l==0) {
909: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
910: } else {
911: DSSetState(eps->ds,DS_STATE_RAW);
912: }
913: BVSetActiveColumns(eps->V,eps->nconv,nv);
915: /* Solve projected problem and compute residual norm estimates */
916: if (eps->its == 1 && l > 0) {/* After rational update */
917: DSGetArray(eps->ds,DS_MAT_A,&A);
918: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
919: b = a + ld;
920: k = eps->nconv+l;
921: A[k*ld+k-1] = A[(k-1)*ld+k];
922: A[k*ld+k] = a[k];
923: for (j=k+1; j< nv; j++) {
924: A[j*ld+j] = a[j];
925: A[j*ld+j-1] = b[j-1] ;
926: A[(j-1)*ld+j] = b[j-1];
927: }
928: DSRestoreArray(eps->ds,DS_MAT_A,&A);
929: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
930: DSSolve(eps->ds,eps->eigr,NULL);
931: DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
932: DSSetCompact(eps->ds,PETSC_TRUE);
933: } else { /* Restart */
934: DSSolve(eps->ds,eps->eigr,NULL);
935: DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
936: }
937: DSSynchronize(eps->ds,eps->eigr,NULL);
939: /* Residual */
940: EPSKrylovConvergence(eps,PETSC_TRUE,eps->nconv,nv-eps->nconv,beta,0.0,1.0,&k);
941: /* Checking values obtained for completing */
942: for (i=0;i<k;i++) {
943: sr->back[i]=eps->eigr[i];
944: }
945: STBackTransform(eps->st,k,sr->back,eps->eigi);
946: count0=count1=0;
947: for (i=0;i<k;i++) {
948: lambda = PetscRealPart(sr->back[i]);
949: if (((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0)) count0++;
950: if (((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0)) count1++;
951: }
952: if (k>eps->nev && eps->ncv-k<5) eps->reason = EPS_CONVERGED_TOL;
953: else {
954: /* Checks completion */
955: if ((!sch0||count0 >= sPres->nsch[0]) && (!sch1 ||count1 >= sPres->nsch[1])) {
956: eps->reason = EPS_CONVERGED_TOL;
957: } else {
958: if (!complIterating && eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
959: if (complIterating) {
960: if (--iterCompl <= 0) eps->reason = EPS_DIVERGED_ITS;
961: } else if (k >= eps->nev) {
962: n0 = sPres->nsch[0]-count0;
963: n1 = sPres->nsch[1]-count1;
964: if (sr->iterCompl>0 && ((n0>0 && n0<= sr->nMAXCompl)||(n1>0&&n1<=sr->nMAXCompl))) {
965: /* Iterating for completion*/
966: complIterating = PETSC_TRUE;
967: if (n0 >sr->nMAXCompl)sch0 = PETSC_FALSE;
968: if (n1 >sr->nMAXCompl)sch1 = PETSC_FALSE;
969: iterCompl = sr->iterCompl;
970: } else eps->reason = EPS_CONVERGED_TOL;
971: }
972: }
973: }
974: /* Update l */
975: if (eps->reason == EPS_CONVERGED_ITERATING) l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
976: else l = nv-k;
977: if (breakdown) l=0;
978: if (!ctx->lock && l>0 && eps->reason == EPS_CONVERGED_ITERATING) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
980: if (eps->reason == EPS_CONVERGED_ITERATING) {
981: if (breakdown) {
982: /* Start a new Lanczos factorization */
983: PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
984: EPSGetStartVector(eps,k,&breakdown);
985: if (breakdown) {
986: eps->reason = EPS_DIVERGED_BREAKDOWN;
987: PetscInfo(eps,"Unable to generate more start vectors\n");
988: }
989: } else {
990: /* Prepare the Rayleigh quotient for restart */
991: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
992: DSGetArray(eps->ds,DS_MAT_Q,&Q);
993: b = a + ld;
994: for (i=k;i<k+l;i++) {
995: a[i] = PetscRealPart(eps->eigr[i]);
996: b[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
997: }
998: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
999: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
1000: }
1001: }
1002: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1003: DSGetMat(eps->ds,DS_MAT_Q,&U);
1004: BVMultInPlace(eps->V,U,eps->nconv,k+l);
1005: MatDestroy(&U);
1007: /* Normalize u and append it to V */
1008: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
1009: BVCopyColumn(eps->V,nv,k+l);
1010: }
1011: eps->nconv = k;
1012: if (eps->reason != EPS_CONVERGED_ITERATING) {
1013: /* Store approximated values for next shift */
1014: DSGetArray(eps->ds,DS_MAT_Q,&Q);
1015: sr->nS = l;
1016: for (i=0;i<l;i++) {
1017: sr->S[i] = eps->eigr[i+k];/* Diagonal elements */
1018: sr->S[i+l] = Q[nv-1+(i+k)*ld]*beta; /* Out of diagonal elements */
1019: }
1020: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
1021: }
1022: }
1023: /* Check for completion */
1024: for (i=0;i< eps->nconv; i++) {
1025: if ((sr->dir)*PetscRealPart(eps->eigr[i])>0) sPres->nconv[1]++;
1026: else sPres->nconv[0]++;
1027: }
1028: sPres->comp[0] = PetscNot(count0 < sPres->nsch[0]);
1029: sPres->comp[1] = PetscNot(count1 < sPres->nsch[1]);
1030: 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])?"*":"");
1031: PetscInfo(eps,messg);
1032: if (count0 > sPres->nsch[0] || count1 > sPres->nsch[1]) InertiaMismatch(eps,ctx->detect);
1033: PetscFree(iwork);
1034: return(0);
1035: }
1037: /*
1038: Obtains value of subsequent shift
1039: */
1040: static PetscErrorCode EPSGetNewShiftValue(EPS eps,PetscInt side,PetscReal *newS)
1041: {
1042: PetscReal lambda,d_prev;
1043: PetscInt i,idxP;
1044: EPS_SR sr;
1045: EPS_shift sPres,s;
1046: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1049: sr = ctx->sr;
1050: sPres = sr->sPres;
1051: if (sPres->neighb[side]) {
1052: /* Completing a previous interval */
1053: *newS = (sPres->value + sPres->neighb[side]->value)/2;
1054: if (PetscAbsReal(sPres->value - *newS)/PetscAbsReal(sPres->value)<=100*PETSC_SQRT_MACHINE_EPSILON) *newS = sPres->value;
1055: } else { /* (Only for side=1). Creating a new interval. */
1056: if (sPres->neigs==0) {/* No value has been accepted*/
1057: if (sPres->neighb[0]) {
1058: /* Multiplying by 10 the previous distance */
1059: *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
1060: sr->nleap++;
1061: /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
1062: if (!sr->hasEnd && sr->nleap > 5) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Unable to compute the wanted eigenvalues with open interval");
1063: } else { /* First shift */
1064: if (eps->nconv != 0) {
1065: /* Unaccepted values give information for next shift */
1066: idxP=0;/* Number of values left from shift */
1067: for (i=0;i<eps->nconv;i++) {
1068: lambda = PetscRealPart(eps->eigr[i]);
1069: if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
1070: else break;
1071: }
1072: /* Avoiding subtraction of eigenvalues (might be the same).*/
1073: if (idxP>0) {
1074: d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[0]))/(idxP+0.3);
1075: } else {
1076: d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[eps->nconv-1]))/(eps->nconv+0.3);
1077: }
1078: *newS = sPres->value + ((sr->dir)*d_prev*eps->nev)/2;
1079: } else { /* No values found, no information for next shift */
1080: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"First shift renders no information");
1081: }
1082: }
1083: } else { /* Accepted values found */
1084: sr->nleap = 0;
1085: /* Average distance of values in previous subinterval */
1086: s = sPres->neighb[0];
1087: while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
1088: s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
1089: }
1090: if (s) {
1091: d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
1092: } else { /* First shift. Average distance obtained with values in this shift */
1093: /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
1094: 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)) {
1095: d_prev = PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0])))/(sPres->neigs+0.3);
1096: } else {
1097: d_prev = PetscAbsReal(PetscRealPart(sr->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
1098: }
1099: }
1100: /* Average distance is used for next shift by adding it to value on the right or to shift */
1101: if ((sr->dir)*(PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
1102: *newS = PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(eps->nev))/2;
1103: } else { /* Last accepted value is on the left of shift. Adding to shift */
1104: *newS = sPres->value + ((sr->dir)*d_prev*(eps->nev))/2;
1105: }
1106: }
1107: /* End of interval can not be surpassed */
1108: if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
1109: }/* of neighb[side]==null */
1110: return(0);
1111: }
1113: /*
1114: Function for sorting an array of real values
1115: */
1116: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
1117: {
1118: PetscReal re;
1119: PetscInt i,j,tmp;
1122: if (!prev) for (i=0;i<nr;i++) perm[i] = i;
1123: /* Insertion sort */
1124: for (i=1;i<nr;i++) {
1125: re = PetscRealPart(r[perm[i]]);
1126: j = i-1;
1127: while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
1128: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
1129: }
1130: }
1131: return(0);
1132: }
1134: /* Stores the pairs obtained since the last shift in the global arrays */
1135: static PetscErrorCode EPSStoreEigenpairs(EPS eps)
1136: {
1137: PetscErrorCode ierr;
1138: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1139: PetscReal lambda,err,norm;
1140: PetscInt i,count;
1141: PetscBool iscayley;
1142: EPS_SR sr = ctx->sr;
1143: EPS_shift sPres;
1144: Vec v,w;
1147: sPres = sr->sPres;
1148: sPres->index = sr->indexEig;
1149: count = sr->indexEig;
1150: /* Back-transform */
1151: STBackTransform(eps->st,eps->nconv,eps->eigr,eps->eigi);
1152: PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
1153: /* Sort eigenvalues */
1154: sortRealEigenvalues(eps->eigr,eps->perm,eps->nconv,PETSC_FALSE,sr->dir);
1155: /* Values stored in global array */
1156: for (i=0;i<eps->nconv;i++) {
1157: lambda = PetscRealPart(eps->eigr[eps->perm[i]]);
1158: err = eps->errest[eps->perm[i]];
1160: if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
1161: if (count>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Unexpected error in Spectrum Slicing");
1162: sr->eigr[count] = lambda;
1163: sr->errest[count] = err;
1164: /* Explicit purification */
1165: BVGetColumn(eps->V,eps->perm[i],&w);
1166: if (eps->purify) {
1167: BVGetColumn(sr->V,count,&v);
1168: STApply(eps->st,w,v);
1169: BVRestoreColumn(sr->V,count,&v);
1170: } else {
1171: BVInsertVec(sr->V,count,w);
1172: }
1173: BVRestoreColumn(eps->V,eps->perm[i],&w);
1174: BVNormColumn(sr->V,count,NORM_2,&norm);
1175: BVScaleColumn(sr->V,count,1.0/norm);
1176: count++;
1177: }
1178: }
1179: sPres->neigs = count - sr->indexEig;
1180: sr->indexEig = count;
1181: /* Global ordering array updating */
1182: sortRealEigenvalues(sr->eigr,sr->perm,count,PETSC_TRUE,sr->dir);
1183: return(0);
1184: }
1186: static PetscErrorCode EPSLookForDeflation(EPS eps)
1187: {
1188: PetscErrorCode ierr;
1189: PetscReal val;
1190: PetscInt i,count0=0,count1=0;
1191: EPS_shift sPres;
1192: PetscInt ini,fin,k,idx0,idx1;
1193: EPS_SR sr;
1194: Vec v;
1195: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1198: sr = ctx->sr;
1199: sPres = sr->sPres;
1201: if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
1202: else ini = 0;
1203: fin = sr->indexEig;
1204: /* Selection of ends for searching new values */
1205: if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
1206: else sPres->ext[0] = sPres->neighb[0]->value;
1207: if (!sPres->neighb[1]) {
1208: if (sr->hasEnd) sPres->ext[1] = sr->int1;
1209: else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
1210: } else sPres->ext[1] = sPres->neighb[1]->value;
1211: /* Selection of values between right and left ends */
1212: for (i=ini;i<fin;i++) {
1213: val=PetscRealPart(sr->eigr[sr->perm[i]]);
1214: /* Values to the right of left shift */
1215: if ((sr->dir)*(val - sPres->ext[1]) < 0) {
1216: if ((sr->dir)*(val - sPres->value) < 0) count0++;
1217: else count1++;
1218: } else break;
1219: }
1220: /* The number of values on each side are found */
1221: if (sPres->neighb[0]) {
1222: sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
1223: if (sPres->nsch[0]<0) InertiaMismatch(eps,ctx->detect);
1224: } else sPres->nsch[0] = 0;
1226: if (sPres->neighb[1]) {
1227: sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
1228: if (sPres->nsch[1]<0) InertiaMismatch(eps,ctx->detect);
1229: } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);
1231: /* Completing vector of indexes for deflation */
1232: idx0 = ini;
1233: idx1 = ini+count0+count1;
1234: k=0;
1235: for (i=idx0;i<idx1;i++) sr->idxDef[k++]=sr->perm[i];
1236: BVDuplicateResize(eps->V,k+eps->ncv+1,&sr->Vnext);
1237: BVSetNumConstraints(sr->Vnext,k);
1238: for (i=0;i<k;i++) {
1239: BVGetColumn(sr->Vnext,-i-1,&v);
1240: BVCopyVec(sr->V,sr->idxDef[i],v);
1241: BVRestoreColumn(sr->Vnext,-i-1,&v);
1242: }
1244: /* For rational Krylov */
1245: if (sr->nS>0 && (sr->sPrev == sr->sPres->neighb[0] || sr->sPrev == sr->sPres->neighb[1])) {
1246: EPSPrepareRational(eps);
1247: }
1248: eps->nconv = 0;
1249: /* Get rid of temporary Vnext */
1250: BVDestroy(&eps->V);
1251: eps->V = sr->Vnext;
1252: sr->Vnext = NULL;
1253: return(0);
1254: }
1256: PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS eps)
1257: {
1258: PetscErrorCode ierr;
1259: PetscInt i,lds,ti;
1260: PetscReal newS;
1261: EPS_KRYLOVSCHUR *ctx=(EPS_KRYLOVSCHUR*)eps->data;
1262: EPS_SR sr=ctx->sr;
1263: Mat A,B=NULL;
1264: PetscObjectState Astate,Bstate=0;
1265: PetscObjectId Aid,Bid=0;
1268: PetscCitationsRegister(citation,&cited);
1269: if (ctx->global) {
1270: EPSSolve_KrylovSchur_Slice(ctx->eps);
1271: ctx->eps->state = EPS_STATE_SOLVED;
1272: eps->reason = EPS_CONVERGED_TOL;
1273: if (ctx->npart>1) {
1274: /* Gather solution from subsolvers */
1275: EPSSliceGatherSolution(eps);
1276: } else {
1277: eps->nconv = sr->numEigs;
1278: eps->its = ctx->eps->its;
1279: PetscFree(ctx->inertias);
1280: PetscFree(ctx->shifts);
1281: EPSSliceGetInertias(ctx->eps,&ctx->nshifts,&ctx->shifts,&ctx->inertias);
1282: }
1283: } else {
1284: if (ctx->npart==1) {
1285: sr->eigr = ctx->eps->eigr;
1286: sr->eigi = ctx->eps->eigi;
1287: sr->perm = ctx->eps->perm;
1288: sr->errest = ctx->eps->errest;
1289: sr->V = ctx->eps->V;
1290: }
1291: /* Check that the user did not modify subcomm matrices */
1292: EPSGetOperators(eps,&A,&B);
1293: PetscObjectStateGet((PetscObject)A,&Astate);
1294: PetscObjectGetId((PetscObject)A,&Aid);
1295: if (B) {
1296: PetscObjectStateGet((PetscObject)B,&Bstate);
1297: PetscObjectGetId((PetscObject)B,&Bid);
1298: }
1299: 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");
1300: /* Only with eigenvalues present in the interval ...*/
1301: if (sr->numEigs==0) {
1302: eps->reason = EPS_CONVERGED_TOL;
1303: return(0);
1304: }
1305: /* Array of pending shifts */
1306: sr->maxPend = 100; /* Initial size */
1307: sr->nPend = 0;
1308: PetscMalloc1(sr->maxPend,&sr->pending);
1309: PetscLogObjectMemory((PetscObject)eps,sr->maxPend*sizeof(EPS_shift*));
1310: EPSCreateShift(eps,sr->int0,NULL,NULL);
1311: /* extract first shift */
1312: sr->sPrev = NULL;
1313: sr->sPres = sr->pending[--sr->nPend];
1314: sr->sPres->inertia = sr->inertia0;
1315: eps->target = sr->sPres->value;
1316: sr->s0 = sr->sPres;
1317: sr->indexEig = 0;
1318: /* Memory reservation for auxiliary variables */
1319: lds = PetscMin(eps->mpd,eps->ncv);
1320: PetscCalloc1(lds*lds,&sr->S);
1321: PetscMalloc1(eps->ncv,&sr->back);
1322: PetscLogObjectMemory((PetscObject)eps,(lds*lds+eps->ncv)*sizeof(PetscScalar));
1323: for (i=0;i<sr->numEigs;i++) {
1324: sr->eigr[i] = 0.0;
1325: sr->eigi[i] = 0.0;
1326: sr->errest[i] = 0.0;
1327: sr->perm[i] = i;
1328: }
1329: /* Vectors for deflation */
1330: PetscMalloc1(sr->numEigs,&sr->idxDef);
1331: PetscLogObjectMemory((PetscObject)eps,sr->numEigs*sizeof(PetscInt));
1332: sr->indexEig = 0;
1333: /* Main loop */
1334: while (sr->sPres) {
1335: /* Search for deflation */
1336: EPSLookForDeflation(eps);
1337: /* KrylovSchur */
1338: EPSKrylovSchur_Slice(eps);
1340: EPSStoreEigenpairs(eps);
1341: /* Select new shift */
1342: if (!sr->sPres->comp[1]) {
1343: EPSGetNewShiftValue(eps,1,&newS);
1344: EPSCreateShift(eps,newS,sr->sPres,sr->sPres->neighb[1]);
1345: }
1346: if (!sr->sPres->comp[0]) {
1347: /* Completing earlier interval */
1348: EPSGetNewShiftValue(eps,0,&newS);
1349: EPSCreateShift(eps,newS,sr->sPres->neighb[0],sr->sPres);
1350: }
1351: /* Preparing for a new search of values */
1352: EPSExtractShift(eps);
1353: }
1355: /* Updating eps values prior to exit */
1356: PetscFree(sr->S);
1357: PetscFree(sr->idxDef);
1358: PetscFree(sr->pending);
1359: PetscFree(sr->back);
1360: BVDuplicateResize(eps->V,eps->ncv+1,&sr->Vnext);
1361: BVSetNumConstraints(sr->Vnext,0);
1362: BVDestroy(&eps->V);
1363: eps->V = sr->Vnext;
1364: eps->nconv = sr->indexEig;
1365: eps->reason = EPS_CONVERGED_TOL;
1366: eps->its = sr->itsKs;
1367: eps->nds = 0;
1368: if (sr->dir<0) {
1369: for (i=0;i<eps->nconv/2;i++) {
1370: ti = sr->perm[i]; sr->perm[i] = sr->perm[eps->nconv-1-i]; sr->perm[eps->nconv-1-i] = ti;
1371: }
1372: }
1373: }
1374: return(0);
1375: }