Actual source code: ciss.c
slepc-3.18.0 2022-10-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: "ciss"
13: Method: Contour Integral Spectral Slicing
15: Algorithm:
17: Contour integral based on Sakurai-Sugiura method to construct a
18: subspace, with various eigenpair extractions (Rayleigh-Ritz,
19: explicit moment).
21: Based on code contributed by Y. Maeda, T. Sakurai.
23: References:
25: [1] T. Sakurai and H. Sugiura, "A projection method for generalized
26: eigenvalue problems", J. Comput. Appl. Math. 159:119-128, 2003.
28: [2] T. Sakurai and H. Tadano, "CIRR: a Rayleigh-Ritz type method with
29: contour integral for generalized eigenvalue problems", Hokkaido
30: Math. J. 36:745-757, 2007.
31: */
33: #include <slepc/private/epsimpl.h>
34: #include <slepc/private/slepccontour.h>
35: #include <slepcblaslapack.h>
37: typedef struct {
38: /* user parameters */
39: PetscInt N; /* number of integration points (32) */
40: PetscInt L; /* block size (16) */
41: PetscInt M; /* moment degree (N/4 = 4) */
42: PetscReal delta; /* threshold of singular value (1e-12) */
43: PetscInt L_max; /* maximum number of columns of the source matrix V */
44: PetscReal spurious_threshold; /* discard spurious eigenpairs */
45: PetscBool isreal; /* A and B are real */
46: PetscInt npart; /* number of partitions */
47: PetscInt refine_inner;
48: PetscInt refine_blocksize;
49: EPSCISSQuadRule quad;
50: EPSCISSExtraction extraction;
51: PetscBool usest;
52: /* private data */
53: SlepcContourData contour;
54: PetscReal *sigma; /* threshold for numerical rank */
55: PetscScalar *weight;
56: PetscScalar *omega;
57: PetscScalar *pp;
58: BV V;
59: BV S;
60: BV pV;
61: BV Y;
62: PetscBool useconj;
63: PetscBool usest_set; /* whether the user set the usest flag or not */
64: PetscObjectId rgid;
65: PetscObjectState rgstate;
66: } EPS_CISS;
68: /*
69: Set up KSP solvers for every integration point, only called if !ctx->usest
70: */
71: static PetscErrorCode EPSCISSSetUp(EPS eps,Mat A,Mat B,Mat Pa,Mat Pb)
72: {
73: EPS_CISS *ctx = (EPS_CISS*)eps->data;
74: SlepcContourData contour;
75: PetscInt i,p_id,nsplit;
76: Mat Amat,Pmat;
77: MatStructure str,strp;
79: if (!ctx->contour || !ctx->contour->ksp) EPSCISSGetKSPs(eps,NULL,NULL);
80: contour = ctx->contour;
81: STGetMatStructure(eps->st,&str);
82: STGetSplitPreconditionerInfo(eps->st,&nsplit,&strp);
83: for (i=0;i<contour->npoints;i++) {
84: p_id = i*contour->subcomm->n + contour->subcomm->color;
85: MatDuplicate(A,MAT_COPY_VALUES,&Amat);
86: if (B) MatAXPY(Amat,-ctx->omega[p_id],B,str);
87: else MatShift(Amat,-ctx->omega[p_id]);
88: if (nsplit) {
89: MatDuplicate(Pa,MAT_COPY_VALUES,&Pmat);
90: if (Pb) MatAXPY(Pmat,-ctx->omega[p_id],Pb,strp);
91: else MatShift(Pmat,-ctx->omega[p_id]);
92: } else Pmat = Amat;
93: EPS_KSPSetOperators(contour->ksp[i],Amat,Amat);
94: MatDestroy(&Amat);
95: if (nsplit) MatDestroy(&Pmat);
96: }
97: return 0;
98: }
100: /*
101: Y_i = (A-z_i B)^{-1}BV for every integration point, Y=[Y_i] is in the context
102: */
103: static PetscErrorCode EPSCISSSolve(EPS eps,Mat B,BV V,PetscInt L_start,PetscInt L_end)
104: {
105: EPS_CISS *ctx = (EPS_CISS*)eps->data;
106: SlepcContourData contour;
107: PetscInt i,p_id;
108: Mat MV,BMV=NULL,MC;
109: KSP ksp;
111: if (!ctx->contour || !ctx->contour->ksp) EPSCISSGetKSPs(eps,NULL,NULL);
112: contour = ctx->contour;
113: BVSetActiveColumns(V,L_start,L_end);
114: BVGetMat(V,&MV);
115: for (i=0;i<contour->npoints;i++) {
116: p_id = i*contour->subcomm->n + contour->subcomm->color;
117: if (ctx->usest) {
118: STSetShift(eps->st,ctx->omega[p_id]);
119: STGetKSP(eps->st,&ksp);
120: } else ksp = contour->ksp[i];
121: BVSetActiveColumns(ctx->Y,i*ctx->L+L_start,i*ctx->L+L_end);
122: BVGetMat(ctx->Y,&MC);
123: if (B) {
124: if (!i) {
125: MatProductCreate(B,MV,NULL,&BMV);
126: MatProductSetType(BMV,MATPRODUCT_AB);
127: MatProductSetFromOptions(BMV);
128: MatProductSymbolic(BMV);
129: }
130: MatProductNumeric(BMV);
131: KSPMatSolve(ksp,BMV,MC);
132: } else KSPMatSolve(ksp,MV,MC);
133: BVRestoreMat(ctx->Y,&MC);
134: if (ctx->usest && i<contour->npoints-1) KSPReset(ksp);
135: }
136: MatDestroy(&BMV);
137: BVRestoreMat(V,&MV);
138: return 0;
139: }
141: static PetscErrorCode rescale_eig(EPS eps,PetscInt nv)
142: {
143: EPS_CISS *ctx = (EPS_CISS*)eps->data;
144: PetscInt i;
145: PetscScalar center;
146: PetscReal radius,a,b,c,d,rgscale;
147: #if defined(PETSC_USE_COMPLEX)
148: PetscReal start_ang,end_ang,vscale,theta;
149: #endif
150: PetscBool isring,isellipse,isinterval;
152: PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
153: PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring);
154: PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval);
155: RGGetScale(eps->rg,&rgscale);
156: if (isinterval) {
157: RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d);
158: if (c==d) {
159: for (i=0;i<nv;i++) {
160: #if defined(PETSC_USE_COMPLEX)
161: eps->eigr[i] = PetscRealPart(eps->eigr[i]);
162: #else
163: eps->eigi[i] = 0;
164: #endif
165: }
166: }
167: }
168: if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
169: if (isellipse) {
170: RGEllipseGetParameters(eps->rg,¢er,&radius,NULL);
171: for (i=0;i<nv;i++) eps->eigr[i] = rgscale*(center + radius*eps->eigr[i]);
172: } else if (isinterval) {
173: RGIntervalGetEndpoints(eps->rg,&a,&b,&c,&d);
174: if (ctx->quad == EPS_CISS_QUADRULE_CHEBYSHEV) {
175: for (i=0;i<nv;i++) {
176: if (c==d) eps->eigr[i] = ((b-a)*(eps->eigr[i]+1.0)/2.0+a)*rgscale;
177: if (a==b) {
178: #if defined(PETSC_USE_COMPLEX)
179: eps->eigr[i] = ((d-c)*(eps->eigr[i]+1.0)/2.0+c)*rgscale*PETSC_i;
180: #else
181: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Integration points on a vertical line require complex arithmetic");
182: #endif
183: }
184: }
185: } else {
186: center = (b+a)/2.0+(d+c)/2.0*PETSC_PI;
187: radius = PetscSqrtReal(PetscPowRealInt((b-a)/2.0,2)+PetscPowRealInt((d-c)/2.0,2));
188: for (i=0;i<nv;i++) eps->eigr[i] = center + radius*eps->eigr[i];
189: }
190: } else if (isring) { /* only supported in complex scalars */
191: #if defined(PETSC_USE_COMPLEX)
192: RGRingGetParameters(eps->rg,¢er,&radius,&vscale,&start_ang,&end_ang,NULL);
193: if (ctx->quad == EPS_CISS_QUADRULE_CHEBYSHEV) {
194: for (i=0;i<nv;i++) {
195: theta = (start_ang*2.0+(end_ang-start_ang)*(PetscRealPart(eps->eigr[i])+1.0))*PETSC_PI;
196: eps->eigr[i] = rgscale*center + (rgscale*radius+PetscImaginaryPart(eps->eigr[i]))*PetscCMPLX(PetscCosReal(theta),vscale*PetscSinReal(theta));
197: }
198: } else {
199: for (i=0;i<nv;i++) eps->eigr[i] = rgscale*(center + radius*eps->eigr[i]);
200: }
201: #endif
202: }
203: }
204: return 0;
205: }
207: PetscErrorCode EPSSetUp_CISS(EPS eps)
208: {
209: EPS_CISS *ctx = (EPS_CISS*)eps->data;
210: SlepcContourData contour;
211: PetscBool istrivial,isring,isellipse,isinterval,flg;
212: PetscReal c,d;
213: PetscInt nsplit;
214: PetscRandom rand;
215: PetscObjectId id;
216: PetscObjectState state;
217: Mat A[2],Psplit[2];
218: Vec v0;
220: if (eps->ncv==PETSC_DEFAULT) {
221: eps->ncv = ctx->L_max*ctx->M;
222: if (eps->ncv>eps->n) {
223: eps->ncv = eps->n;
224: ctx->L_max = eps->ncv/ctx->M;
226: }
227: } else {
228: EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
229: ctx->L_max = eps->ncv/ctx->M;
230: if (!ctx->L_max) {
231: ctx->L_max = 1;
232: eps->ncv = ctx->L_max*ctx->M;
233: }
234: }
235: ctx->L = PetscMin(ctx->L,ctx->L_max);
236: if (eps->max_it==PETSC_DEFAULT) eps->max_it = 5;
237: if (eps->mpd==PETSC_DEFAULT) eps->mpd = eps->ncv;
238: if (!eps->which) eps->which = EPS_ALL;
240: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_EXTRACTION | EPS_FEATURE_STOPPING | EPS_FEATURE_TWOSIDED);
242: /* check region */
243: RGIsTrivial(eps->rg,&istrivial);
245: RGGetComplement(eps->rg,&flg);
247: PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
248: PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring);
249: PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval);
252: /* if the region has changed, then reset contour data */
253: PetscObjectGetId((PetscObject)eps->rg,&id);
254: PetscObjectStateGet((PetscObject)eps->rg,&state);
255: if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
256: SlepcContourDataDestroy(&ctx->contour);
257: PetscInfo(eps,"Resetting the contour data structure due to a change of region\n");
258: ctx->rgid = id; ctx->rgstate = state;
259: }
261: #if !defined(PETSC_USE_COMPLEX)
263: #endif
264: if (isinterval) {
265: RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d);
266: #if !defined(PETSC_USE_COMPLEX)
268: #endif
269: if (!ctx->quad && c==d) ctx->quad = EPS_CISS_QUADRULE_CHEBYSHEV;
270: }
271: if (!ctx->quad) ctx->quad = EPS_CISS_QUADRULE_TRAPEZOIDAL;
273: /* create contour data structure */
274: if (!ctx->contour) {
275: RGCanUseConjugates(eps->rg,ctx->isreal,&ctx->useconj);
276: SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)eps,&ctx->contour);
277: }
279: EPSAllocateSolution(eps,0);
280: BVGetRandomContext(eps->V,&rand); /* make sure the random context is available when duplicating */
281: if (ctx->weight) PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
282: PetscMalloc4(ctx->N,&ctx->weight,ctx->N+1,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma);
284: /* allocate basis vectors */
285: BVDestroy(&ctx->S);
286: BVDuplicateResize(eps->V,ctx->L*ctx->M,&ctx->S);
287: BVDestroy(&ctx->V);
288: BVDuplicateResize(eps->V,ctx->L,&ctx->V);
290: STGetMatrix(eps->st,0,&A[0]);
291: MatIsShell(A[0],&flg);
293: if (eps->isgeneralized) STGetMatrix(eps->st,1,&A[1]);
295: if (!ctx->usest_set) ctx->usest = (ctx->npart>1)? PETSC_FALSE: PETSC_TRUE;
298: /* check if a user-defined split preconditioner has been set */
299: STGetSplitPreconditionerInfo(eps->st,&nsplit,NULL);
300: if (nsplit) {
301: STGetSplitPreconditionerTerm(eps->st,0,&Psplit[0]);
302: if (eps->isgeneralized) STGetSplitPreconditionerTerm(eps->st,1,&Psplit[1]);
303: }
305: contour = ctx->contour;
306: SlepcContourRedundantMat(contour,eps->isgeneralized?2:1,A,nsplit?Psplit:NULL);
307: if (contour->pA) {
308: BVGetColumn(ctx->V,0,&v0);
309: SlepcContourScatterCreate(contour,v0);
310: BVRestoreColumn(ctx->V,0,&v0);
311: BVDestroy(&ctx->pV);
312: BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV);
313: BVSetSizesFromVec(ctx->pV,contour->xsub,eps->n);
314: BVSetFromOptions(ctx->pV);
315: BVResize(ctx->pV,ctx->L,PETSC_FALSE);
316: }
318: EPSCheckDefinite(eps);
319: EPSCheckSinvertCondition(eps,ctx->usest," (with the usest flag set)");
321: BVDestroy(&ctx->Y);
322: if (contour->pA) {
323: BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y);
324: BVSetSizesFromVec(ctx->Y,contour->xsub,eps->n);
325: BVSetFromOptions(ctx->Y);
326: BVResize(ctx->Y,contour->npoints*ctx->L,PETSC_FALSE);
327: } else BVDuplicateResize(eps->V,contour->npoints*ctx->L,&ctx->Y);
329: if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) DSSetType(eps->ds,DSGNHEP);
330: else if (eps->isgeneralized) {
331: if (eps->ishermitian && eps->ispositive) DSSetType(eps->ds,DSGHEP);
332: else DSSetType(eps->ds,DSGNHEP);
333: } else {
334: if (eps->ishermitian) DSSetType(eps->ds,DSHEP);
335: else DSSetType(eps->ds,DSNHEP);
336: }
337: DSAllocate(eps->ds,eps->ncv);
339: #if !defined(PETSC_USE_COMPLEX)
340: EPSSetWorkVecs(eps,3);
341: if (!eps->ishermitian) PetscInfo(eps,"Warning: complex eigenvalues are not calculated exactly without --with-scalar-type=complex in PETSc\n");
342: #else
343: EPSSetWorkVecs(eps,2);
344: #endif
345: return 0;
346: }
348: PetscErrorCode EPSSetUpSort_CISS(EPS eps)
349: {
350: SlepcSC sc;
352: /* fill sorting criterion context */
353: eps->sc->comparison = SlepcCompareSmallestReal;
354: eps->sc->comparisonctx = NULL;
355: eps->sc->map = NULL;
356: eps->sc->mapobj = NULL;
358: /* fill sorting criterion for DS */
359: DSGetSlepcSC(eps->ds,&sc);
360: sc->comparison = SlepcCompareLargestMagnitude;
361: sc->comparisonctx = NULL;
362: sc->map = NULL;
363: sc->mapobj = NULL;
364: return 0;
365: }
367: PetscErrorCode EPSSolve_CISS(EPS eps)
368: {
369: EPS_CISS *ctx = (EPS_CISS*)eps->data;
370: SlepcContourData contour = ctx->contour;
371: Mat A,B,X,M,pA,pB,T,J,Pa=NULL,Pb=NULL;
372: BV V;
373: PetscInt i,j,ld,nmat,L_add=0,nv=0,L_base=ctx->L,inner,nlocal,*inside,nsplit;
374: PetscScalar *Mu,*H0,*H1=NULL,*rr,*temp;
375: PetscReal error,max_error,norm;
376: PetscBool *fl1;
377: Vec si,si1=NULL,w[3];
378: PetscRandom rand;
379: #if defined(PETSC_USE_COMPLEX)
380: PetscBool isellipse;
381: PetscReal est_eig,eta;
382: #else
383: PetscReal normi;
384: #endif
386: w[0] = eps->work[0];
387: #if defined(PETSC_USE_COMPLEX)
388: w[1] = NULL;
389: #else
390: w[1] = eps->work[2];
391: #endif
392: w[2] = eps->work[1];
393: VecGetLocalSize(w[0],&nlocal);
394: DSGetLeadingDimension(eps->ds,&ld);
395: RGComputeQuadrature(eps->rg,ctx->quad==EPS_CISS_QUADRULE_CHEBYSHEV?RG_QUADRULE_CHEBYSHEV:RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight);
396: STGetNumMatrices(eps->st,&nmat);
397: STGetMatrix(eps->st,0,&A);
398: if (nmat>1) STGetMatrix(eps->st,1,&B);
399: else B = NULL;
400: J = (contour->pA && nmat>1)? contour->pA[1]: B;
401: V = contour->pA? ctx->pV: ctx->V;
402: if (!ctx->usest) {
403: T = contour->pA? contour->pA[0]: A;
404: STGetSplitPreconditionerInfo(eps->st,&nsplit,NULL);
405: if (nsplit) {
406: if (contour->pA) {
407: Pa = contour->pP[0];
408: if (nsplit>1) Pb = contour->pP[1];
409: } else {
410: STGetSplitPreconditionerTerm(eps->st,0,&Pa);
411: if (nsplit>1) STGetSplitPreconditionerTerm(eps->st,1,&Pb);
412: }
413: }
414: EPSCISSSetUp(eps,T,J,Pa,Pb);
415: }
416: BVSetActiveColumns(ctx->V,0,ctx->L);
417: BVSetRandomSign(ctx->V);
418: BVGetRandomContext(ctx->V,&rand);
420: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
421: EPSCISSSolve(eps,J,V,0,ctx->L);
422: #if defined(PETSC_USE_COMPLEX)
423: PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
424: if (isellipse) {
425: BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig);
426: PetscInfo(eps,"Estimated eigenvalue count: %f\n",(double)est_eig);
427: eta = PetscPowReal(10.0,-PetscLog10Real(eps->tol)/ctx->N);
428: L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
429: if (L_add>ctx->L_max-ctx->L) {
430: PetscInfo(eps,"Number of eigenvalues inside the contour path may be too large\n");
431: L_add = ctx->L_max-ctx->L;
432: }
433: }
434: #endif
435: if (L_add>0) {
436: PetscInfo(eps,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by Estimate #Eig\n",ctx->L,ctx->L+L_add);
437: BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints);
438: BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
439: BVSetRandomSign(ctx->V);
440: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
441: ctx->L += L_add;
442: EPSCISSSolve(eps,J,V,ctx->L-L_add,ctx->L);
443: }
444: PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
445: for (i=0;i<ctx->refine_blocksize;i++) {
446: BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj);
447: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
448: PetscLogEventBegin(EPS_CISS_SVD,eps,0,0,0);
449: SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
450: PetscLogEventEnd(EPS_CISS_SVD,eps,0,0,0);
451: if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
452: L_add = L_base;
453: if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
454: PetscInfo(eps,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by SVD(H0)\n",ctx->L,ctx->L+L_add);
455: BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints);
456: BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
457: BVSetRandomSign(ctx->V);
458: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
459: ctx->L += L_add;
460: EPSCISSSolve(eps,J,V,ctx->L-L_add,ctx->L);
461: if (L_add) {
462: PetscFree2(Mu,H0);
463: PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
464: }
465: }
466: if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1);
468: while (eps->reason == EPS_CONVERGED_ITERATING) {
469: eps->its++;
470: for (inner=0;inner<=ctx->refine_inner;inner++) {
471: if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
472: BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj);
473: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
474: PetscLogEventBegin(EPS_CISS_SVD,eps,0,0,0);
475: SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
476: PetscLogEventEnd(EPS_CISS_SVD,eps,0,0,0);
477: break;
478: } else {
479: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
480: BVSetActiveColumns(ctx->S,0,ctx->L);
481: BVSetActiveColumns(ctx->V,0,ctx->L);
482: BVCopy(ctx->S,ctx->V);
483: BVSVDAndRank(ctx->S,ctx->M,ctx->L,ctx->delta,BV_SVD_METHOD_REFINE,H0,ctx->sigma,&nv);
484: if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
485: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
486: EPSCISSSolve(eps,J,V,0,ctx->L);
487: } else break;
488: }
489: }
490: eps->nconv = 0;
491: if (nv == 0) eps->reason = EPS_CONVERGED_TOL;
492: else {
493: DSSetDimensions(eps->ds,nv,0,0);
494: DSSetState(eps->ds,DS_STATE_RAW);
496: if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
497: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
498: CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1);
499: DSGetArray(eps->ds,DS_MAT_A,&temp);
500: for (j=0;j<nv;j++) {
501: for (i=0;i<nv;i++) {
502: temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
503: }
504: }
505: DSRestoreArray(eps->ds,DS_MAT_A,&temp);
506: DSGetArray(eps->ds,DS_MAT_B,&temp);
507: for (j=0;j<nv;j++) {
508: for (i=0;i<nv;i++) {
509: temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
510: }
511: }
512: DSRestoreArray(eps->ds,DS_MAT_B,&temp);
513: } else {
514: BVSetActiveColumns(ctx->S,0,nv);
515: DSGetMat(eps->ds,DS_MAT_A,&pA);
516: MatZeroEntries(pA);
517: BVMatProject(ctx->S,A,ctx->S,pA);
518: DSRestoreMat(eps->ds,DS_MAT_A,&pA);
519: if (B) {
520: DSGetMat(eps->ds,DS_MAT_B,&pB);
521: MatZeroEntries(pB);
522: BVMatProject(ctx->S,B,ctx->S,pB);
523: DSRestoreMat(eps->ds,DS_MAT_B,&pB);
524: }
525: }
527: DSSolve(eps->ds,eps->eigr,eps->eigi);
528: DSSynchronize(eps->ds,eps->eigr,eps->eigi);
530: PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr);
531: rescale_eig(eps,nv);
532: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
533: DSGetMat(eps->ds,DS_MAT_X,&X);
534: SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1);
535: DSRestoreMat(eps->ds,DS_MAT_X,&X);
536: RGCheckInside(eps->rg,nv,eps->eigr,eps->eigi,inside);
537: for (i=0;i<nv;i++) {
538: if (fl1[i] && inside[i]>=0) {
539: rr[i] = 1.0;
540: eps->nconv++;
541: } else rr[i] = 0.0;
542: }
543: DSSort(eps->ds,eps->eigr,eps->eigi,rr,NULL,&eps->nconv);
544: DSSynchronize(eps->ds,eps->eigr,eps->eigi);
545: rescale_eig(eps,nv);
546: PetscFree3(fl1,inside,rr);
547: BVSetActiveColumns(eps->V,0,nv);
548: if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
549: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
550: BVSetActiveColumns(ctx->S,0,ctx->L);
551: BVCopy(ctx->S,ctx->V);
552: BVSetActiveColumns(ctx->S,0,nv);
553: }
554: BVCopy(ctx->S,eps->V);
556: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
557: DSGetMat(eps->ds,DS_MAT_X,&X);
558: BVMultInPlace(ctx->S,X,0,eps->nconv);
559: if (eps->ishermitian) BVMultInPlace(eps->V,X,0,eps->nconv);
560: DSRestoreMat(eps->ds,DS_MAT_X,&X);
561: max_error = 0.0;
562: for (i=0;i<eps->nconv;i++) {
563: BVGetColumn(ctx->S,i,&si);
564: #if !defined(PETSC_USE_COMPLEX)
565: if (eps->eigi[i]!=0.0) BVGetColumn(ctx->S,i+1,&si1);
566: #endif
567: EPSComputeResidualNorm_Private(eps,PETSC_FALSE,eps->eigr[i],eps->eigi[i],si,si1,w,&error);
568: if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) { /* vector is not normalized */
569: VecNorm(si,NORM_2,&norm);
570: #if !defined(PETSC_USE_COMPLEX)
571: if (eps->eigi[i]!=0.0) {
572: VecNorm(si1,NORM_2,&normi);
573: norm = SlepcAbsEigenvalue(norm,normi);
574: }
575: #endif
576: error /= norm;
577: }
578: (*eps->converged)(eps,eps->eigr[i],eps->eigi[i],error,&error,eps->convergedctx);
579: BVRestoreColumn(ctx->S,i,&si);
580: #if !defined(PETSC_USE_COMPLEX)
581: if (eps->eigi[i]!=0.0) {
582: BVRestoreColumn(ctx->S,i+1,&si1);
583: i++;
584: }
585: #endif
586: max_error = PetscMax(max_error,error);
587: }
589: if (max_error <= eps->tol) eps->reason = EPS_CONVERGED_TOL;
590: else if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
591: else {
592: if (eps->nconv > ctx->L) nv = eps->nconv;
593: else if (ctx->L > nv) nv = ctx->L;
594: nv = PetscMin(nv,ctx->L*ctx->M);
595: MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M);
596: MatSetRandom(M,rand);
597: BVSetActiveColumns(ctx->S,0,nv);
598: BVMultInPlace(ctx->S,M,0,ctx->L);
599: MatDestroy(&M);
600: BVSetActiveColumns(ctx->S,0,ctx->L);
601: BVSetActiveColumns(ctx->V,0,ctx->L);
602: BVCopy(ctx->S,ctx->V);
603: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
604: EPSCISSSolve(eps,J,V,0,ctx->L);
605: }
606: }
607: }
608: if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) PetscFree(H1);
609: PetscFree2(Mu,H0);
610: return 0;
611: }
613: PetscErrorCode EPSComputeVectors_CISS(EPS eps)
614: {
615: EPS_CISS *ctx = (EPS_CISS*)eps->data;
616: PetscInt n;
617: Mat Z,B=NULL;
619: if (eps->ishermitian) {
620: if (eps->isgeneralized && !eps->ispositive) EPSComputeVectors_Indefinite(eps);
621: else EPSComputeVectors_Hermitian(eps);
622: if (eps->isgeneralized && eps->ispositive && ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
623: /* normalize to have unit B-norm */
624: STGetMatrix(eps->st,1,&B);
625: BVSetMatrix(eps->V,B,PETSC_FALSE);
626: BVNormalize(eps->V,NULL);
627: BVSetMatrix(eps->V,NULL,PETSC_FALSE);
628: }
629: return 0;
630: }
631: DSGetDimensions(eps->ds,&n,NULL,NULL,NULL);
632: BVSetActiveColumns(eps->V,0,n);
634: /* right eigenvectors */
635: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
637: /* V = V * Z */
638: DSGetMat(eps->ds,DS_MAT_X,&Z);
639: BVMultInPlace(eps->V,Z,0,n);
640: DSRestoreMat(eps->ds,DS_MAT_X,&Z);
641: BVSetActiveColumns(eps->V,0,eps->nconv);
643: /* normalize */
644: if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) BVNormalize(eps->V,NULL);
645: return 0;
646: }
648: static PetscErrorCode EPSCISSSetSizes_CISS(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
649: {
650: EPS_CISS *ctx = (EPS_CISS*)eps->data;
651: PetscInt oN,oL,oM,oLmax,onpart;
652: PetscMPIInt size;
654: oN = ctx->N;
655: if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
656: if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
657: } else {
660: if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
661: }
662: oL = ctx->L;
663: if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
664: ctx->L = 16;
665: } else {
667: ctx->L = bs;
668: }
669: oM = ctx->M;
670: if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
671: ctx->M = ctx->N/4;
672: } else {
675: ctx->M = ms;
676: }
677: onpart = ctx->npart;
678: if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
679: ctx->npart = 1;
680: } else {
681: MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size);
683: ctx->npart = npart;
684: }
685: oLmax = ctx->L_max;
686: if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
687: ctx->L_max = 64;
688: } else {
690: ctx->L_max = PetscMax(bsmax,ctx->L);
691: }
692: if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
693: SlepcContourDataDestroy(&ctx->contour);
694: PetscInfo(eps,"Resetting the contour data structure due to a change of parameters\n");
695: eps->state = EPS_STATE_INITIAL;
696: }
697: ctx->isreal = realmats;
698: if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) eps->state = EPS_STATE_INITIAL;
699: return 0;
700: }
702: /*@
703: EPSCISSSetSizes - Sets the values of various size parameters in the CISS solver.
705: Logically Collective on eps
707: Input Parameters:
708: + eps - the eigenproblem solver context
709: . ip - number of integration points
710: . bs - block size
711: . ms - moment size
712: . npart - number of partitions when splitting the communicator
713: . bsmax - max block size
714: - realmats - A and B are real
716: Options Database Keys:
717: + -eps_ciss_integration_points - Sets the number of integration points
718: . -eps_ciss_blocksize - Sets the block size
719: . -eps_ciss_moments - Sets the moment size
720: . -eps_ciss_partitions - Sets the number of partitions
721: . -eps_ciss_maxblocksize - Sets the maximum block size
722: - -eps_ciss_realmats - A and B are real
724: Note:
725: The default number of partitions is 1. This means the internal KSP object is shared
726: among all processes of the EPS communicator. Otherwise, the communicator is split
727: into npart communicators, so that npart KSP solves proceed simultaneously.
729: Level: advanced
731: .seealso: EPSCISSGetSizes()
732: @*/
733: PetscErrorCode EPSCISSSetSizes(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
734: {
742: PetscTryMethod(eps,"EPSCISSSetSizes_C",(EPS,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(eps,ip,bs,ms,npart,bsmax,realmats));
743: return 0;
744: }
746: static PetscErrorCode EPSCISSGetSizes_CISS(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
747: {
748: EPS_CISS *ctx = (EPS_CISS*)eps->data;
750: if (ip) *ip = ctx->N;
751: if (bs) *bs = ctx->L;
752: if (ms) *ms = ctx->M;
753: if (npart) *npart = ctx->npart;
754: if (bsmax) *bsmax = ctx->L_max;
755: if (realmats) *realmats = ctx->isreal;
756: return 0;
757: }
759: /*@
760: EPSCISSGetSizes - Gets the values of various size parameters in the CISS solver.
762: Not Collective
764: Input Parameter:
765: . eps - the eigenproblem solver context
767: Output Parameters:
768: + ip - number of integration points
769: . bs - block size
770: . ms - moment size
771: . npart - number of partitions when splitting the communicator
772: . bsmax - max block size
773: - realmats - A and B are real
775: Level: advanced
777: .seealso: EPSCISSSetSizes()
778: @*/
779: PetscErrorCode EPSCISSGetSizes(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
780: {
782: PetscUseMethod(eps,"EPSCISSGetSizes_C",(EPS,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(eps,ip,bs,ms,npart,bsmax,realmats));
783: return 0;
784: }
786: static PetscErrorCode EPSCISSSetThreshold_CISS(EPS eps,PetscReal delta,PetscReal spur)
787: {
788: EPS_CISS *ctx = (EPS_CISS*)eps->data;
790: if (delta == PETSC_DEFAULT) {
791: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
792: } else {
794: ctx->delta = delta;
795: }
796: if (spur == PETSC_DEFAULT) {
797: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
798: } else {
800: ctx->spurious_threshold = spur;
801: }
802: return 0;
803: }
805: /*@
806: EPSCISSSetThreshold - Sets the values of various threshold parameters in
807: the CISS solver.
809: Logically Collective on eps
811: Input Parameters:
812: + eps - the eigenproblem solver context
813: . delta - threshold for numerical rank
814: - spur - spurious threshold (to discard spurious eigenpairs)
816: Options Database Keys:
817: + -eps_ciss_delta - Sets the delta
818: - -eps_ciss_spurious_threshold - Sets the spurious threshold
820: Level: advanced
822: .seealso: EPSCISSGetThreshold()
823: @*/
824: PetscErrorCode EPSCISSSetThreshold(EPS eps,PetscReal delta,PetscReal spur)
825: {
829: PetscTryMethod(eps,"EPSCISSSetThreshold_C",(EPS,PetscReal,PetscReal),(eps,delta,spur));
830: return 0;
831: }
833: static PetscErrorCode EPSCISSGetThreshold_CISS(EPS eps,PetscReal *delta,PetscReal *spur)
834: {
835: EPS_CISS *ctx = (EPS_CISS*)eps->data;
837: if (delta) *delta = ctx->delta;
838: if (spur) *spur = ctx->spurious_threshold;
839: return 0;
840: }
842: /*@
843: EPSCISSGetThreshold - Gets the values of various threshold parameters
844: in the CISS solver.
846: Not Collective
848: Input Parameter:
849: . eps - the eigenproblem solver context
851: Output Parameters:
852: + delta - threshold for numerical rank
853: - spur - spurious threshold (to discard spurious eigenpairs)
855: Level: advanced
857: .seealso: EPSCISSSetThreshold()
858: @*/
859: PetscErrorCode EPSCISSGetThreshold(EPS eps,PetscReal *delta,PetscReal *spur)
860: {
862: PetscUseMethod(eps,"EPSCISSGetThreshold_C",(EPS,PetscReal*,PetscReal*),(eps,delta,spur));
863: return 0;
864: }
866: static PetscErrorCode EPSCISSSetRefinement_CISS(EPS eps,PetscInt inner,PetscInt blsize)
867: {
868: EPS_CISS *ctx = (EPS_CISS*)eps->data;
870: if (inner == PETSC_DEFAULT) {
871: ctx->refine_inner = 0;
872: } else {
874: ctx->refine_inner = inner;
875: }
876: if (blsize == PETSC_DEFAULT) {
877: ctx->refine_blocksize = 0;
878: } else {
880: ctx->refine_blocksize = blsize;
881: }
882: return 0;
883: }
885: /*@
886: EPSCISSSetRefinement - Sets the values of various refinement parameters
887: in the CISS solver.
889: Logically Collective on eps
891: Input Parameters:
892: + eps - the eigenproblem solver context
893: . inner - number of iterative refinement iterations (inner loop)
894: - blsize - number of iterative refinement iterations (blocksize loop)
896: Options Database Keys:
897: + -eps_ciss_refine_inner - Sets number of inner iterations
898: - -eps_ciss_refine_blocksize - Sets number of blocksize iterations
900: Level: advanced
902: .seealso: EPSCISSGetRefinement()
903: @*/
904: PetscErrorCode EPSCISSSetRefinement(EPS eps,PetscInt inner,PetscInt blsize)
905: {
909: PetscTryMethod(eps,"EPSCISSSetRefinement_C",(EPS,PetscInt,PetscInt),(eps,inner,blsize));
910: return 0;
911: }
913: static PetscErrorCode EPSCISSGetRefinement_CISS(EPS eps,PetscInt *inner,PetscInt *blsize)
914: {
915: EPS_CISS *ctx = (EPS_CISS*)eps->data;
917: if (inner) *inner = ctx->refine_inner;
918: if (blsize) *blsize = ctx->refine_blocksize;
919: return 0;
920: }
922: /*@
923: EPSCISSGetRefinement - Gets the values of various refinement parameters
924: in the CISS solver.
926: Not Collective
928: Input Parameter:
929: . eps - the eigenproblem solver context
931: Output Parameters:
932: + inner - number of iterative refinement iterations (inner loop)
933: - blsize - number of iterative refinement iterations (blocksize loop)
935: Level: advanced
937: .seealso: EPSCISSSetRefinement()
938: @*/
939: PetscErrorCode EPSCISSGetRefinement(EPS eps, PetscInt *inner, PetscInt *blsize)
940: {
942: PetscUseMethod(eps,"EPSCISSGetRefinement_C",(EPS,PetscInt*,PetscInt*),(eps,inner,blsize));
943: return 0;
944: }
946: static PetscErrorCode EPSCISSSetUseST_CISS(EPS eps,PetscBool usest)
947: {
948: EPS_CISS *ctx = (EPS_CISS*)eps->data;
950: ctx->usest = usest;
951: ctx->usest_set = PETSC_TRUE;
952: eps->state = EPS_STATE_INITIAL;
953: return 0;
954: }
956: /*@
957: EPSCISSSetUseST - Sets a flag indicating that the CISS solver will
958: use the ST object for the linear solves.
960: Logically Collective on eps
962: Input Parameters:
963: + eps - the eigenproblem solver context
964: - usest - boolean flag to use the ST object or not
966: Options Database Keys:
967: . -eps_ciss_usest <bool> - whether the ST object will be used or not
969: Level: advanced
971: .seealso: EPSCISSGetUseST()
972: @*/
973: PetscErrorCode EPSCISSSetUseST(EPS eps,PetscBool usest)
974: {
977: PetscTryMethod(eps,"EPSCISSSetUseST_C",(EPS,PetscBool),(eps,usest));
978: return 0;
979: }
981: static PetscErrorCode EPSCISSGetUseST_CISS(EPS eps,PetscBool *usest)
982: {
983: EPS_CISS *ctx = (EPS_CISS*)eps->data;
985: *usest = ctx->usest;
986: return 0;
987: }
989: /*@
990: EPSCISSGetUseST - Gets the flag for using the ST object
991: in the CISS solver.
993: Not Collective
995: Input Parameter:
996: . eps - the eigenproblem solver context
998: Output Parameters:
999: . usest - boolean flag indicating if the ST object is being used
1001: Level: advanced
1003: .seealso: EPSCISSSetUseST()
1004: @*/
1005: PetscErrorCode EPSCISSGetUseST(EPS eps,PetscBool *usest)
1006: {
1009: PetscUseMethod(eps,"EPSCISSGetUseST_C",(EPS,PetscBool*),(eps,usest));
1010: return 0;
1011: }
1013: static PetscErrorCode EPSCISSSetQuadRule_CISS(EPS eps,EPSCISSQuadRule quad)
1014: {
1015: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1017: if (ctx->quad != quad) {
1018: ctx->quad = quad;
1019: eps->state = EPS_STATE_INITIAL;
1020: }
1021: return 0;
1022: }
1024: /*@
1025: EPSCISSSetQuadRule - Sets the quadrature rule used in the CISS solver.
1027: Logically Collective on eps
1029: Input Parameters:
1030: + eps - the eigenproblem solver context
1031: - quad - the quadrature rule
1033: Options Database Key:
1034: . -eps_ciss_quadrule - Sets the quadrature rule (either 'trapezoidal' or
1035: 'chebyshev')
1037: Notes:
1038: By default, the trapezoidal rule is used (EPS_CISS_QUADRULE_TRAPEZOIDAL).
1040: If the 'chebyshev' option is specified (EPS_CISS_QUADRULE_CHEBYSHEV), then
1041: Chebyshev points are used as quadrature points.
1043: Level: advanced
1045: .seealso: EPSCISSGetQuadRule(), EPSCISSQuadRule
1046: @*/
1047: PetscErrorCode EPSCISSSetQuadRule(EPS eps,EPSCISSQuadRule quad)
1048: {
1051: PetscTryMethod(eps,"EPSCISSSetQuadRule_C",(EPS,EPSCISSQuadRule),(eps,quad));
1052: return 0;
1053: }
1055: static PetscErrorCode EPSCISSGetQuadRule_CISS(EPS eps,EPSCISSQuadRule *quad)
1056: {
1057: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1059: *quad = ctx->quad;
1060: return 0;
1061: }
1063: /*@
1064: EPSCISSGetQuadRule - Gets the quadrature rule used in the CISS solver.
1066: Not Collective
1068: Input Parameter:
1069: . eps - the eigenproblem solver context
1071: Output Parameters:
1072: . quad - quadrature rule
1074: Level: advanced
1076: .seealso: EPSCISSSetQuadRule() EPSCISSQuadRule
1077: @*/
1078: PetscErrorCode EPSCISSGetQuadRule(EPS eps,EPSCISSQuadRule *quad)
1079: {
1082: PetscUseMethod(eps,"EPSCISSGetQuadRule_C",(EPS,EPSCISSQuadRule*),(eps,quad));
1083: return 0;
1084: }
1086: static PetscErrorCode EPSCISSSetExtraction_CISS(EPS eps,EPSCISSExtraction extraction)
1087: {
1088: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1090: if (ctx->extraction != extraction) {
1091: ctx->extraction = extraction;
1092: eps->state = EPS_STATE_INITIAL;
1093: }
1094: return 0;
1095: }
1097: /*@
1098: EPSCISSSetExtraction - Sets the extraction technique used in the CISS solver.
1100: Logically Collective on eps
1102: Input Parameters:
1103: + eps - the eigenproblem solver context
1104: - extraction - the extraction technique
1106: Options Database Key:
1107: . -eps_ciss_extraction - Sets the extraction technique (either 'ritz' or
1108: 'hankel')
1110: Notes:
1111: By default, the Rayleigh-Ritz extraction is used (EPS_CISS_EXTRACTION_RITZ).
1113: If the 'hankel' option is specified (EPS_CISS_EXTRACTION_HANKEL), then
1114: the Block Hankel method is used for extracting eigenpairs.
1116: Level: advanced
1118: .seealso: EPSCISSGetExtraction(), EPSCISSExtraction
1119: @*/
1120: PetscErrorCode EPSCISSSetExtraction(EPS eps,EPSCISSExtraction extraction)
1121: {
1124: PetscTryMethod(eps,"EPSCISSSetExtraction_C",(EPS,EPSCISSExtraction),(eps,extraction));
1125: return 0;
1126: }
1128: static PetscErrorCode EPSCISSGetExtraction_CISS(EPS eps,EPSCISSExtraction *extraction)
1129: {
1130: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1132: *extraction = ctx->extraction;
1133: return 0;
1134: }
1136: /*@
1137: EPSCISSGetExtraction - Gets the extraction technique used in the CISS solver.
1139: Not Collective
1141: Input Parameter:
1142: . eps - the eigenproblem solver context
1144: Output Parameters:
1145: . extraction - extraction technique
1147: Level: advanced
1149: .seealso: EPSCISSSetExtraction() EPSCISSExtraction
1150: @*/
1151: PetscErrorCode EPSCISSGetExtraction(EPS eps,EPSCISSExtraction *extraction)
1152: {
1155: PetscUseMethod(eps,"EPSCISSGetExtraction_C",(EPS,EPSCISSExtraction*),(eps,extraction));
1156: return 0;
1157: }
1159: static PetscErrorCode EPSCISSGetKSPs_CISS(EPS eps,PetscInt *nsolve,KSP **ksp)
1160: {
1161: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1162: SlepcContourData contour;
1163: PetscInt i,nsplit;
1164: PC pc;
1165: MPI_Comm child;
1167: if (!ctx->contour) { /* initialize contour data structure first */
1168: RGCanUseConjugates(eps->rg,ctx->isreal,&ctx->useconj);
1169: SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)eps,&ctx->contour);
1170: }
1171: contour = ctx->contour;
1172: if (!contour->ksp) {
1173: PetscMalloc1(contour->npoints,&contour->ksp);
1174: EPSGetST(eps,&eps->st);
1175: STGetSplitPreconditionerInfo(eps->st,&nsplit,NULL);
1176: PetscSubcommGetChild(contour->subcomm,&child);
1177: for (i=0;i<contour->npoints;i++) {
1178: KSPCreate(child,&contour->ksp[i]);
1179: PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)eps,1);
1180: KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)eps)->prefix);
1181: KSPAppendOptionsPrefix(contour->ksp[i],"eps_ciss_");
1182: PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)eps)->options);
1183: KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE);
1184: KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(eps->tol),PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
1185: KSPGetPC(contour->ksp[i],&pc);
1186: if (nsplit) {
1187: KSPSetType(contour->ksp[i],KSPBCGS);
1188: PCSetType(pc,PCBJACOBI);
1189: } else {
1190: KSPSetType(contour->ksp[i],KSPPREONLY);
1191: PCSetType(pc,PCLU);
1192: }
1193: }
1194: }
1195: if (nsolve) *nsolve = contour->npoints;
1196: if (ksp) *ksp = contour->ksp;
1197: return 0;
1198: }
1200: /*@C
1201: EPSCISSGetKSPs - Retrieve the array of linear solver objects associated with
1202: the CISS solver.
1204: Not Collective
1206: Input Parameter:
1207: . eps - the eigenproblem solver solver
1209: Output Parameters:
1210: + nsolve - number of solver objects
1211: - ksp - array of linear solver object
1213: Notes:
1214: The number of KSP solvers is equal to the number of integration points divided by
1215: the number of partitions. This value is halved in the case of real matrices with
1216: a region centered at the real axis.
1218: Level: advanced
1220: .seealso: EPSCISSSetSizes()
1221: @*/
1222: PetscErrorCode EPSCISSGetKSPs(EPS eps,PetscInt *nsolve,KSP **ksp)
1223: {
1225: PetscUseMethod(eps,"EPSCISSGetKSPs_C",(EPS,PetscInt*,KSP**),(eps,nsolve,ksp));
1226: return 0;
1227: }
1229: PetscErrorCode EPSReset_CISS(EPS eps)
1230: {
1231: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1233: BVDestroy(&ctx->S);
1234: BVDestroy(&ctx->V);
1235: BVDestroy(&ctx->Y);
1236: if (!ctx->usest) SlepcContourDataReset(ctx->contour);
1237: BVDestroy(&ctx->pV);
1238: return 0;
1239: }
1241: PetscErrorCode EPSSetFromOptions_CISS(EPS eps,PetscOptionItems *PetscOptionsObject)
1242: {
1243: PetscReal r3,r4;
1244: PetscInt i,i1,i2,i3,i4,i5,i6,i7;
1245: PetscBool b1,b2,flg,flg2,flg3,flg4,flg5,flg6;
1246: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1247: EPSCISSQuadRule quad;
1248: EPSCISSExtraction extraction;
1250: PetscOptionsHeadBegin(PetscOptionsObject,"EPS CISS Options");
1252: EPSCISSGetSizes(eps,&i1,&i2,&i3,&i4,&i5,&b1);
1253: PetscOptionsInt("-eps_ciss_integration_points","Number of integration points","EPSCISSSetSizes",i1,&i1,&flg);
1254: PetscOptionsInt("-eps_ciss_blocksize","Block size","EPSCISSSetSizes",i2,&i2,&flg2);
1255: PetscOptionsInt("-eps_ciss_moments","Moment size","EPSCISSSetSizes",i3,&i3,&flg3);
1256: PetscOptionsInt("-eps_ciss_partitions","Number of partitions","EPSCISSSetSizes",i4,&i4,&flg4);
1257: PetscOptionsInt("-eps_ciss_maxblocksize","Maximum block size","EPSCISSSetSizes",i5,&i5,&flg5);
1258: PetscOptionsBool("-eps_ciss_realmats","True if A and B are real","EPSCISSSetSizes",b1,&b1,&flg6);
1259: if (flg || flg2 || flg3 || flg4 || flg5 || flg6) EPSCISSSetSizes(eps,i1,i2,i3,i4,i5,b1);
1261: EPSCISSGetThreshold(eps,&r3,&r4);
1262: PetscOptionsReal("-eps_ciss_delta","Threshold for numerical rank","EPSCISSSetThreshold",r3,&r3,&flg);
1263: PetscOptionsReal("-eps_ciss_spurious_threshold","Threshold for the spurious eigenpairs","EPSCISSSetThreshold",r4,&r4,&flg2);
1264: if (flg || flg2) EPSCISSSetThreshold(eps,r3,r4);
1266: EPSCISSGetRefinement(eps,&i6,&i7);
1267: PetscOptionsInt("-eps_ciss_refine_inner","Number of inner iterative refinement iterations","EPSCISSSetRefinement",i6,&i6,&flg);
1268: PetscOptionsInt("-eps_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","EPSCISSSetRefinement",i7,&i7,&flg2);
1269: if (flg || flg2) EPSCISSSetRefinement(eps,i6,i7);
1271: EPSCISSGetUseST(eps,&b2);
1272: PetscOptionsBool("-eps_ciss_usest","Use ST for linear solves","EPSCISSSetUseST",b2,&b2,&flg);
1273: if (flg) EPSCISSSetUseST(eps,b2);
1275: PetscOptionsEnum("-eps_ciss_quadrule","Quadrature rule","EPSCISSSetQuadRule",EPSCISSQuadRules,(PetscEnum)ctx->quad,(PetscEnum*)&quad,&flg);
1276: if (flg) EPSCISSSetQuadRule(eps,quad);
1278: PetscOptionsEnum("-eps_ciss_extraction","Extraction technique","EPSCISSSetExtraction",EPSCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg);
1279: if (flg) EPSCISSSetExtraction(eps,extraction);
1281: PetscOptionsHeadEnd();
1283: if (!eps->rg) EPSGetRG(eps,&eps->rg);
1284: RGSetFromOptions(eps->rg); /* this is necessary here to set useconj */
1285: if (!ctx->contour || !ctx->contour->ksp) EPSCISSGetKSPs(eps,NULL,NULL);
1286: for (i=0;i<ctx->contour->npoints;i++) KSPSetFromOptions(ctx->contour->ksp[i]);
1287: PetscSubcommSetFromOptions(ctx->contour->subcomm);
1288: return 0;
1289: }
1291: PetscErrorCode EPSDestroy_CISS(EPS eps)
1292: {
1293: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1295: SlepcContourDataDestroy(&ctx->contour);
1296: PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
1297: PetscFree(eps->data);
1298: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",NULL);
1299: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",NULL);
1300: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",NULL);
1301: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",NULL);
1302: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",NULL);
1303: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",NULL);
1304: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetUseST_C",NULL);
1305: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetUseST_C",NULL);
1306: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetQuadRule_C",NULL);
1307: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetQuadRule_C",NULL);
1308: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetExtraction_C",NULL);
1309: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetExtraction_C",NULL);
1310: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetKSPs_C",NULL);
1311: return 0;
1312: }
1314: PetscErrorCode EPSView_CISS(EPS eps,PetscViewer viewer)
1315: {
1316: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1317: PetscBool isascii;
1318: PetscViewer sviewer;
1320: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1321: if (isascii) {
1322: PetscViewerASCIIPrintf(viewer," sizes { integration points: %" PetscInt_FMT ", block size: %" PetscInt_FMT ", moment size: %" PetscInt_FMT ", partitions: %" PetscInt_FMT ", maximum block size: %" PetscInt_FMT " }\n",ctx->N,ctx->L,ctx->M,ctx->npart,ctx->L_max);
1323: if (ctx->isreal) PetscViewerASCIIPrintf(viewer," exploiting symmetry of integration points\n");
1324: PetscViewerASCIIPrintf(viewer," threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold);
1325: PetscViewerASCIIPrintf(viewer," iterative refinement { inner: %" PetscInt_FMT ", blocksize: %" PetscInt_FMT " }\n",ctx->refine_inner, ctx->refine_blocksize);
1326: PetscViewerASCIIPrintf(viewer," extraction: %s\n",EPSCISSExtractions[ctx->extraction]);
1327: PetscViewerASCIIPrintf(viewer," quadrature rule: %s\n",EPSCISSQuadRules[ctx->quad]);
1328: if (ctx->usest) PetscViewerASCIIPrintf(viewer," using ST for linear solves\n");
1329: else {
1330: if (!ctx->contour || !ctx->contour->ksp) EPSCISSGetKSPs(eps,NULL,NULL);
1331: PetscViewerASCIIPushTab(viewer);
1332: if (ctx->npart>1 && ctx->contour->subcomm) {
1333: PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1334: if (!ctx->contour->subcomm->color) KSPView(ctx->contour->ksp[0],sviewer);
1335: PetscViewerFlush(sviewer);
1336: PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1337: PetscViewerFlush(viewer);
1338: /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1339: PetscViewerASCIIPopSynchronized(viewer);
1340: } else KSPView(ctx->contour->ksp[0],viewer);
1341: PetscViewerASCIIPopTab(viewer);
1342: }
1343: }
1344: return 0;
1345: }
1347: PetscErrorCode EPSSetDefaultST_CISS(EPS eps)
1348: {
1349: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1350: PetscBool usest = ctx->usest;
1351: KSP ksp;
1352: PC pc;
1354: if (!((PetscObject)eps->st)->type_name) {
1355: if (!ctx->usest_set) usest = (ctx->npart>1)? PETSC_FALSE: PETSC_TRUE;
1356: if (usest) STSetType(eps->st,STSINVERT);
1357: else {
1358: /* we are not going to use ST, so avoid factorizing the matrix */
1359: STSetType(eps->st,STSHIFT);
1360: if (eps->isgeneralized) {
1361: STGetKSP(eps->st,&ksp);
1362: KSPGetPC(ksp,&pc);
1363: PCSetType(pc,PCNONE);
1364: }
1365: }
1366: }
1367: return 0;
1368: }
1370: SLEPC_EXTERN PetscErrorCode EPSCreate_CISS(EPS eps)
1371: {
1372: EPS_CISS *ctx = (EPS_CISS*)eps->data;
1374: PetscNew(&ctx);
1375: eps->data = ctx;
1377: eps->useds = PETSC_TRUE;
1378: eps->categ = EPS_CATEGORY_CONTOUR;
1380: eps->ops->solve = EPSSolve_CISS;
1381: eps->ops->setup = EPSSetUp_CISS;
1382: eps->ops->setupsort = EPSSetUpSort_CISS;
1383: eps->ops->setfromoptions = EPSSetFromOptions_CISS;
1384: eps->ops->destroy = EPSDestroy_CISS;
1385: eps->ops->reset = EPSReset_CISS;
1386: eps->ops->view = EPSView_CISS;
1387: eps->ops->computevectors = EPSComputeVectors_CISS;
1388: eps->ops->setdefaultst = EPSSetDefaultST_CISS;
1390: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",EPSCISSSetSizes_CISS);
1391: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",EPSCISSGetSizes_CISS);
1392: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",EPSCISSSetThreshold_CISS);
1393: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",EPSCISSGetThreshold_CISS);
1394: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",EPSCISSSetRefinement_CISS);
1395: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",EPSCISSGetRefinement_CISS);
1396: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetUseST_C",EPSCISSSetUseST_CISS);
1397: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetUseST_C",EPSCISSGetUseST_CISS);
1398: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetQuadRule_C",EPSCISSSetQuadRule_CISS);
1399: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetQuadRule_C",EPSCISSGetQuadRule_CISS);
1400: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetExtraction_C",EPSCISSSetExtraction_CISS);
1401: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetExtraction_C",EPSCISSGetExtraction_CISS);
1402: PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetKSPs_C",EPSCISSGetKSPs_CISS);
1404: /* set default values of parameters */
1405: ctx->N = 32;
1406: ctx->L = 16;
1407: ctx->M = ctx->N/4;
1408: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
1409: ctx->L_max = 64;
1410: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1411: ctx->usest = PETSC_TRUE;
1412: ctx->usest_set = PETSC_FALSE;
1413: ctx->isreal = PETSC_FALSE;
1414: ctx->refine_inner = 0;
1415: ctx->refine_blocksize = 0;
1416: ctx->npart = 1;
1417: ctx->quad = (EPSCISSQuadRule)0;
1418: ctx->extraction = EPS_CISS_EXTRACTION_RITZ;
1419: return 0;
1420: }