Actual source code: pciss.c
slepc-3.22.2 2024-12-02
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] J. Asakura, T. Sakurai, H. Tadano, T. Ikegami, K. Kimura, "A
26: numerical method for polynomial eigenvalue problems using contour
27: integral", Japan J. Indust. Appl. Math. 27:73-90, 2010.
28: */
30: #include <slepc/private/pepimpl.h>
31: #include <slepc/private/slepccontour.h>
33: typedef struct {
34: /* parameters */
35: PetscInt N; /* number of integration points (32) */
36: PetscInt L; /* block size (16) */
37: PetscInt M; /* moment degree (N/4 = 4) */
38: PetscReal delta; /* threshold of singular value (1e-12) */
39: PetscInt L_max; /* maximum number of columns of the source matrix V */
40: PetscReal spurious_threshold; /* discard spurious eigenpairs */
41: PetscBool isreal; /* T(z) is real for real z */
42: PetscInt npart; /* number of partitions */
43: PetscInt refine_inner;
44: PetscInt refine_blocksize;
45: PEPCISSExtraction extraction;
46: /* private data */
47: SlepcContourData contour;
48: PetscReal *sigma; /* threshold for numerical rank */
49: PetscScalar *weight;
50: PetscScalar *omega;
51: PetscScalar *pp;
52: BV V;
53: BV S;
54: BV Y;
55: PetscBool useconj;
56: Mat J,*Psplit; /* auxiliary matrices */
57: BV pV;
58: PetscObjectId rgid;
59: PetscObjectState rgstate;
60: } PEP_CISS;
62: static PetscErrorCode PEPComputeFunction(PEP pep,PetscScalar lambda,Mat T,Mat P,PetscBool deriv)
63: {
64: PetscInt i;
65: PetscScalar *coeff;
66: Mat *A,*K;
67: MatStructure str,strp;
68: PEP_CISS *ctx = (PEP_CISS*)pep->data;
69: SlepcContourData contour = ctx->contour;
71: PetscFunctionBegin;
72: A = (contour->pA)?contour->pA:pep->A;
73: K = (contour->pP)?contour->pP:ctx->Psplit;
74: PetscCall(PetscMalloc1(pep->nmat,&coeff));
75: if (deriv) PetscCall(PEPEvaluateBasisDerivative(pep,lambda,0,coeff,NULL));
76: else PetscCall(PEPEvaluateBasis(pep,lambda,0,coeff,NULL));
77: PetscCall(STGetMatStructure(pep->st,&str));
78: PetscCall(MatZeroEntries(T));
79: if (!deriv && T != P) {
80: PetscCall(STGetSplitPreconditionerInfo(pep->st,NULL,&strp));
81: PetscCall(MatZeroEntries(P));
82: }
83: i = deriv?1:0;
84: for (;i<pep->nmat;i++) {
85: PetscCall(MatAXPY(T,coeff[i],A[i],str));
86: if (!deriv && T != P) PetscCall(MatAXPY(P,coeff[i],K[i],strp));
87: }
88: PetscCall(PetscFree(coeff));
89: PetscFunctionReturn(PETSC_SUCCESS);
90: }
92: /*
93: Set up KSP solvers for every integration point
94: */
95: static PetscErrorCode PEPCISSSetUp(PEP pep,Mat T,Mat P)
96: {
97: PEP_CISS *ctx = (PEP_CISS*)pep->data;
98: SlepcContourData contour;
99: PetscInt i,p_id;
100: Mat Amat,Pmat;
102: PetscFunctionBegin;
103: if (!ctx->contour || !ctx->contour->ksp) PetscCall(PEPCISSGetKSPs(pep,NULL,NULL));
104: contour = ctx->contour;
105: PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Something went wrong with PEPCISSGetKSPs()");
106: for (i=0;i<contour->npoints;i++) {
107: p_id = i*contour->subcomm->n + contour->subcomm->color;
108: PetscCall(MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&Amat));
109: if (T != P) PetscCall(MatDuplicate(P,MAT_DO_NOT_COPY_VALUES,&Pmat)); else Pmat = Amat;
110: PetscCall(PEPComputeFunction(pep,ctx->omega[p_id],Amat,Pmat,PETSC_FALSE));
111: PetscCall(PEP_KSPSetOperators(contour->ksp[i],Amat,Pmat));
112: PetscCall(MatDestroy(&Amat));
113: if (T != P) PetscCall(MatDestroy(&Pmat));
114: }
115: PetscFunctionReturn(PETSC_SUCCESS);
116: }
118: /*
119: Y_i = F(z_i)^{-1}Fp(z_i)V for every integration point, Y=[Y_i] is in the context
120: */
121: static PetscErrorCode PEPCISSSolve(PEP pep,Mat dT,BV V,PetscInt L_start,PetscInt L_end)
122: {
123: PEP_CISS *ctx = (PEP_CISS*)pep->data;
124: SlepcContourData contour;
125: PetscInt i,p_id;
126: Mat MV,BMV=NULL,MC;
128: PetscFunctionBegin;
129: contour = ctx->contour;
130: PetscCall(BVSetActiveColumns(V,L_start,L_end));
131: PetscCall(BVGetMat(V,&MV));
132: for (i=0;i<contour->npoints;i++) {
133: p_id = i*contour->subcomm->n + contour->subcomm->color;
134: PetscCall(PEPComputeFunction(pep,ctx->omega[p_id],dT,NULL,PETSC_TRUE));
135: PetscCall(BVSetActiveColumns(ctx->Y,i*ctx->L+L_start,i*ctx->L+L_end));
136: PetscCall(BVGetMat(ctx->Y,&MC));
137: if (!i) {
138: PetscCall(MatProductCreate(dT,MV,NULL,&BMV));
139: PetscCall(MatProductSetType(BMV,MATPRODUCT_AB));
140: PetscCall(MatProductSetFromOptions(BMV));
141: PetscCall(MatProductSymbolic(BMV));
142: }
143: PetscCall(MatProductNumeric(BMV));
144: PetscCall(KSPMatSolve(contour->ksp[i],BMV,MC));
145: PetscCall(BVRestoreMat(ctx->Y,&MC));
146: }
147: PetscCall(MatDestroy(&BMV));
148: PetscCall(BVRestoreMat(V,&MV));
149: PetscFunctionReturn(PETSC_SUCCESS);
150: }
152: static PetscErrorCode PEPSetUp_CISS(PEP pep)
153: {
154: PEP_CISS *ctx = (PEP_CISS*)pep->data;
155: SlepcContourData contour;
156: PetscInt i,nwork,nsplit;
157: PetscBool istrivial,isellipse,flg;
158: PetscObjectId id;
159: PetscObjectState state;
160: Vec v0;
162: PetscFunctionBegin;
163: if (pep->ncv==PETSC_DETERMINE) pep->ncv = ctx->L_max*ctx->M;
164: else {
165: ctx->L_max = pep->ncv/ctx->M;
166: if (!ctx->L_max) {
167: ctx->L_max = 1;
168: pep->ncv = ctx->L_max*ctx->M;
169: }
170: }
171: ctx->L = PetscMin(ctx->L,ctx->L_max);
172: if (pep->max_it==PETSC_DETERMINE) pep->max_it = 5;
173: if (pep->mpd==PETSC_DETERMINE) pep->mpd = pep->ncv;
174: if (!pep->which) pep->which = PEP_ALL;
175: PetscCheck(pep->which==PEP_ALL,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"This solver supports only computing all eigenvalues");
176: PEPCheckUnsupported(pep,PEP_FEATURE_STOPPING);
177: PEPCheckIgnored(pep,PEP_FEATURE_SCALE);
179: /* check region */
180: PetscCall(RGIsTrivial(pep->rg,&istrivial));
181: PetscCheck(!istrivial,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
182: PetscCall(RGGetComplement(pep->rg,&flg));
183: PetscCheck(!flg,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
184: PetscCall(PetscObjectTypeCompare((PetscObject)pep->rg,RGELLIPSE,&isellipse));
185: PetscCheck(isellipse,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Currently only implemented for elliptic regions");
187: /* if the region has changed, then reset contour data */
188: PetscCall(PetscObjectGetId((PetscObject)pep->rg,&id));
189: PetscCall(PetscObjectStateGet((PetscObject)pep->rg,&state));
190: if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
191: PetscCall(SlepcContourDataDestroy(&ctx->contour));
192: PetscCall(PetscInfo(pep,"Resetting the contour data structure due to a change of region\n"));
193: ctx->rgid = id; ctx->rgstate = state;
194: }
196: /* create contour data structure */
197: if (!ctx->contour) {
198: PetscCall(RGCanUseConjugates(pep->rg,ctx->isreal,&ctx->useconj));
199: PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)pep,&ctx->contour));
200: }
202: PetscCall(PEPAllocateSolution(pep,0));
203: if (ctx->weight) PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
204: PetscCall(PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma));
206: /* allocate basis vectors */
207: PetscCall(BVDestroy(&ctx->S));
208: PetscCall(BVDuplicateResize(pep->V,ctx->L*ctx->M,&ctx->S));
209: PetscCall(BVDestroy(&ctx->V));
210: PetscCall(BVDuplicateResize(pep->V,ctx->L,&ctx->V));
212: /* check if a user-defined split preconditioner has been set */
213: PetscCall(STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL));
214: if (nsplit) {
215: PetscCall(PetscFree(ctx->Psplit));
216: PetscCall(PetscMalloc1(nsplit,&ctx->Psplit));
217: for (i=0;i<nsplit;i++) PetscCall(STGetSplitPreconditionerTerm(pep->st,i,&ctx->Psplit[i]));
218: }
220: contour = ctx->contour;
221: PetscCall(SlepcContourRedundantMat(contour,pep->nmat,pep->A,ctx->Psplit));
222: if (!ctx->J) PetscCall(MatDuplicate(contour->pA?contour->pA[0]:pep->A[0],MAT_DO_NOT_COPY_VALUES,&ctx->J));
223: if (contour->pA) {
224: PetscCall(BVGetColumn(ctx->V,0,&v0));
225: PetscCall(SlepcContourScatterCreate(contour,v0));
226: PetscCall(BVRestoreColumn(ctx->V,0,&v0));
227: PetscCall(BVDestroy(&ctx->pV));
228: PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV));
229: PetscCall(BVSetSizesFromVec(ctx->pV,contour->xsub,pep->n));
230: PetscCall(BVSetFromOptions(ctx->pV));
231: PetscCall(BVResize(ctx->pV,ctx->L,PETSC_FALSE));
232: }
234: PetscCall(BVDestroy(&ctx->Y));
235: if (contour->pA) {
236: PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y));
237: PetscCall(BVSetSizesFromVec(ctx->Y,contour->xsub,pep->n));
238: PetscCall(BVSetFromOptions(ctx->Y));
239: PetscCall(BVResize(ctx->Y,contour->npoints*ctx->L,PETSC_FALSE));
240: } else PetscCall(BVDuplicateResize(pep->V,contour->npoints*ctx->L,&ctx->Y));
242: if (ctx->extraction == PEP_CISS_EXTRACTION_RITZ) {
243: PetscCall(DSPEPSetDegree(pep->ds,pep->nmat-1));
244: PetscCall(DSPEPSetCoefficients(pep->ds,pep->pbc));
245: }
246: PetscCall(DSAllocate(pep->ds,pep->ncv));
247: nwork = 2;
248: PetscCall(PEPSetWorkVecs(pep,nwork));
249: PetscFunctionReturn(PETSC_SUCCESS);
250: }
252: static PetscErrorCode PEPSolve_CISS(PEP pep)
253: {
254: PEP_CISS *ctx = (PEP_CISS*)pep->data;
255: SlepcContourData contour = ctx->contour;
256: Mat X,M,E,T,P;
257: PetscInt i,j,ld,L_add=0,nv=0,L_base=ctx->L,inner,*inside,nsplit;
258: PetscScalar *Mu,*H0,*H1,*rr,*temp,center;
259: PetscReal error,max_error,radius,rgscale,est_eig,eta;
260: PetscBool isellipse,*fl1;
261: Vec si;
262: SlepcSC sc;
263: PetscRandom rand;
265: PetscFunctionBegin;
266: PetscCall(DSGetSlepcSC(pep->ds,&sc));
267: sc->comparison = SlepcCompareLargestMagnitude;
268: sc->comparisonctx = NULL;
269: sc->map = NULL;
270: sc->mapobj = NULL;
271: PetscCall(DSGetLeadingDimension(pep->ds,&ld));
272: PetscCall(RGComputeQuadrature(pep->rg,RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight));
273: PetscCall(STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL));
274: if (contour->pA) {
275: T = contour->pA[0];
276: P = nsplit? contour->pP[0]: T;
277: } else {
278: T = pep->A[0];
279: P = nsplit? ctx->Psplit[0]: T;
280: }
281: PetscCall(PEPCISSSetUp(pep,T,P));
282: PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
283: PetscCall(BVSetRandomSign(ctx->V));
284: PetscCall(BVGetRandomContext(ctx->V,&rand));
285: if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
286: PetscCall(PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L));
287: PetscCall(PetscObjectTypeCompare((PetscObject)pep->rg,RGELLIPSE,&isellipse));
288: if (isellipse) {
289: PetscCall(BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig));
290: PetscCall(PetscInfo(pep,"Estimated eigenvalue count: %f\n",(double)est_eig));
291: eta = PetscPowReal(10.0,-PetscLog10Real(pep->tol)/ctx->N);
292: L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
293: if (L_add>ctx->L_max-ctx->L) {
294: PetscCall(PetscInfo(pep,"Number of eigenvalues inside the contour path may be too large\n"));
295: L_add = ctx->L_max-ctx->L;
296: }
297: }
298: /* Updates L after estimate the number of eigenvalue */
299: if (L_add>0) {
300: PetscCall(PetscInfo(pep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by Estimate #Eig\n",ctx->L,ctx->L+L_add));
301: PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
302: PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
303: PetscCall(BVSetRandomSign(ctx->V));
304: if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
305: ctx->L += L_add;
306: PetscCall(PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,ctx->L-L_add,ctx->L));
307: }
309: PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
310: for (i=0;i<ctx->refine_blocksize;i++) {
311: PetscCall(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));
312: PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
313: PetscCall(PetscLogEventBegin(PEP_CISS_SVD,pep,0,0,0));
314: PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
315: PetscCall(PetscLogEventEnd(PEP_CISS_SVD,pep,0,0,0));
316: if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
317: L_add = L_base;
318: if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
319: PetscCall(PetscInfo(pep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by SVD(H0)\n",ctx->L,ctx->L+L_add));
320: PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
321: PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
322: PetscCall(BVSetRandomSign(ctx->V));
323: if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
324: ctx->L += L_add;
325: PetscCall(PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,ctx->L-L_add,ctx->L));
326: if (L_add) {
327: PetscCall(PetscFree2(Mu,H0));
328: PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
329: }
330: }
332: PetscCall(RGGetScale(pep->rg,&rgscale));
333: PetscCall(RGEllipseGetParameters(pep->rg,¢er,&radius,NULL));
335: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) PetscCall(PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1));
337: while (pep->reason == PEP_CONVERGED_ITERATING) {
338: pep->its++;
339: for (inner=0;inner<=ctx->refine_inner;inner++) {
340: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
341: PetscCall(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));
342: PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
343: PetscCall(PetscLogEventBegin(PEP_CISS_SVD,pep,0,0,0));
344: PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
345: PetscCall(PetscLogEventEnd(PEP_CISS_SVD,pep,0,0,0));
346: } else {
347: PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
348: /* compute SVD of S */
349: PetscCall(BVSVDAndRank(ctx->S,ctx->M,ctx->L,ctx->delta,(ctx->extraction==PEP_CISS_EXTRACTION_CAA)?BV_SVD_METHOD_QR_CAA:BV_SVD_METHOD_QR,H0,ctx->sigma,&nv));
350: }
351: PetscCall(PetscInfo(pep,"Estimated rank: nv = %" PetscInt_FMT "\n",nv));
352: if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
353: PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
354: PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
355: PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
356: PetscCall(BVCopy(ctx->S,ctx->V));
357: if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
358: PetscCall(PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L));
359: } else break;
360: }
361: pep->nconv = 0;
362: if (nv == 0) { pep->reason = PEP_CONVERGED_TOL; break; }
363: else {
364: /* Extracting eigenpairs */
365: PetscCall(DSSetDimensions(pep->ds,nv,0,0));
366: PetscCall(DSSetState(pep->ds,DS_STATE_RAW));
367: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
368: PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
369: PetscCall(CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1));
370: PetscCall(DSGetArray(pep->ds,DS_MAT_A,&temp));
371: for (j=0;j<nv;j++)
372: for (i=0;i<nv;i++)
373: temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
374: PetscCall(DSRestoreArray(pep->ds,DS_MAT_A,&temp));
375: PetscCall(DSGetArray(pep->ds,DS_MAT_B,&temp));
376: for (j=0;j<nv;j++)
377: for (i=0;i<nv;i++)
378: temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
379: PetscCall(DSRestoreArray(pep->ds,DS_MAT_B,&temp));
380: } else if (ctx->extraction == PEP_CISS_EXTRACTION_CAA) {
381: PetscCall(BVSetActiveColumns(ctx->S,0,nv));
382: PetscCall(DSGetArray(pep->ds,DS_MAT_A,&temp));
383: for (i=0;i<nv;i++) PetscCall(PetscArraycpy(temp+i*ld,H0+i*nv,nv));
384: PetscCall(DSRestoreArray(pep->ds,DS_MAT_A,&temp));
385: } else {
386: PetscCall(BVSetActiveColumns(ctx->S,0,nv));
387: for (i=0;i<pep->nmat;i++) {
388: PetscCall(DSGetMat(pep->ds,DSMatExtra[i],&E));
389: PetscCall(BVMatProject(ctx->S,pep->A[i],ctx->S,E));
390: PetscCall(DSRestoreMat(pep->ds,DSMatExtra[i],&E));
391: }
392: nv = (pep->nmat-1)*nv;
393: }
394: PetscCall(DSSolve(pep->ds,pep->eigr,pep->eigi));
395: PetscCall(DSSynchronize(pep->ds,pep->eigr,pep->eigi));
396: if (ctx->extraction == PEP_CISS_EXTRACTION_CAA || ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
397: for (i=0;i<nv;i++) {
398: pep->eigr[i] = (pep->eigr[i]*radius+center)*rgscale;
399: }
400: }
401: PetscCall(PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr));
402: PetscCall(DSVectors(pep->ds,DS_MAT_X,NULL,NULL));
403: PetscCall(DSGetMat(pep->ds,DS_MAT_X,&X));
404: PetscCall(SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1));
405: PetscCall(DSRestoreMat(pep->ds,DS_MAT_X,&X));
406: PetscCall(RGCheckInside(pep->rg,nv,pep->eigr,pep->eigi,inside));
407: for (i=0;i<nv;i++) {
408: if (fl1[i] && inside[i]>=0) {
409: rr[i] = 1.0;
410: pep->nconv++;
411: } else rr[i] = 0.0;
412: }
413: PetscCall(DSSort(pep->ds,pep->eigr,pep->eigi,rr,NULL,&pep->nconv));
414: PetscCall(DSSynchronize(pep->ds,pep->eigr,pep->eigi));
415: if (ctx->extraction == PEP_CISS_EXTRACTION_CAA || ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
416: for (i=0;i<nv;i++) pep->eigr[i] = (pep->eigr[i]*radius+center)*rgscale;
417: }
418: PetscCall(PetscFree3(fl1,inside,rr));
419: PetscCall(BVSetActiveColumns(pep->V,0,nv));
420: PetscCall(DSVectors(pep->ds,DS_MAT_X,NULL,NULL));
421: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
422: PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
423: PetscCall(BVSetActiveColumns(ctx->S,0,nv));
424: PetscCall(BVCopy(ctx->S,pep->V));
425: PetscCall(DSGetMat(pep->ds,DS_MAT_X,&X));
426: PetscCall(BVMultInPlace(ctx->S,X,0,pep->nconv));
427: PetscCall(BVMultInPlace(pep->V,X,0,pep->nconv));
428: PetscCall(DSRestoreMat(pep->ds,DS_MAT_X,&X));
429: } else {
430: PetscCall(DSGetMat(pep->ds,DS_MAT_X,&X));
431: PetscCall(BVMultInPlace(ctx->S,X,0,pep->nconv));
432: PetscCall(DSRestoreMat(pep->ds,DS_MAT_X,&X));
433: PetscCall(BVSetActiveColumns(ctx->S,0,pep->nconv));
434: PetscCall(BVCopy(ctx->S,pep->V));
435: }
436: max_error = 0.0;
437: for (i=0;i<pep->nconv;i++) {
438: PetscCall(BVGetColumn(pep->V,i,&si));
439: PetscCall(VecNormalize(si,NULL));
440: PetscCall(PEPComputeResidualNorm_Private(pep,pep->eigr[i],0,si,NULL,pep->work,&error));
441: PetscCall((*pep->converged)(pep,pep->eigr[i],0,error,&error,pep->convergedctx));
442: PetscCall(BVRestoreColumn(pep->V,i,&si));
443: max_error = PetscMax(max_error,error);
444: }
445: if (max_error <= pep->tol) pep->reason = PEP_CONVERGED_TOL;
446: else if (pep->its > pep->max_it) pep->reason = PEP_DIVERGED_ITS;
447: else {
448: if (pep->nconv > ctx->L) nv = pep->nconv;
449: else if (ctx->L > nv) nv = ctx->L;
450: nv = PetscMin(nv,ctx->L*ctx->M);
451: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M));
452: PetscCall(MatSetRandom(M,rand));
453: PetscCall(BVSetActiveColumns(ctx->S,0,nv));
454: PetscCall(BVMultInPlace(ctx->S,M,0,ctx->L));
455: PetscCall(MatDestroy(&M));
456: PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
457: PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
458: PetscCall(BVCopy(ctx->S,ctx->V));
459: if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
460: PetscCall(PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L));
461: }
462: }
463: }
464: PetscCall(PetscFree2(Mu,H0));
465: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) PetscCall(PetscFree(H1));
466: PetscFunctionReturn(PETSC_SUCCESS);
467: }
469: static PetscErrorCode PEPCISSSetSizes_CISS(PEP pep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
470: {
471: PEP_CISS *ctx = (PEP_CISS*)pep->data;
472: PetscInt oN,oL,oM,oLmax,onpart;
473: PetscMPIInt size;
475: PetscFunctionBegin;
476: oN = ctx->N;
477: if (ip == PETSC_DETERMINE) {
478: if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
479: } else if (ip != PETSC_CURRENT) {
480: PetscCheck(ip>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
481: PetscCheck(ip%2==0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
482: if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
483: }
484: oL = ctx->L;
485: if (bs == PETSC_DETERMINE) {
486: ctx->L = 16;
487: } else if (bs != PETSC_CURRENT) {
488: PetscCheck(bs>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
489: ctx->L = bs;
490: }
491: oM = ctx->M;
492: if (ms == PETSC_DETERMINE) {
493: ctx->M = ctx->N/4;
494: } else if (ms != PETSC_CURRENT) {
495: PetscCheck(ms>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
496: PetscCheck(ms<=ctx->N,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be less than or equal to the number of integration points");
497: ctx->M = ms;
498: }
499: onpart = ctx->npart;
500: if (npart == PETSC_DETERMINE) {
501: ctx->npart = 1;
502: } else if (npart != PETSC_CURRENT) {
503: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size));
504: PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
505: ctx->npart = npart;
506: }
507: oLmax = ctx->L_max;
508: if (bsmax == PETSC_DETERMINE) {
509: ctx->L_max = 64;
510: } else if (bsmax != PETSC_CURRENT) {
511: PetscCheck(bsmax>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
512: ctx->L_max = PetscMax(bsmax,ctx->L);
513: }
514: if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
515: PetscCall(SlepcContourDataDestroy(&ctx->contour));
516: PetscCall(PetscInfo(pep,"Resetting the contour data structure due to a change of parameters\n"));
517: pep->state = PEP_STATE_INITIAL;
518: }
519: ctx->isreal = realmats;
520: if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) pep->state = PEP_STATE_INITIAL;
521: PetscFunctionReturn(PETSC_SUCCESS);
522: }
524: /*@
525: PEPCISSSetSizes - Sets the values of various size parameters in the CISS solver.
527: Logically Collective
529: Input Parameters:
530: + pep - the polynomial eigensolver context
531: . ip - number of integration points
532: . bs - block size
533: . ms - moment size
534: . npart - number of partitions when splitting the communicator
535: . bsmax - max block size
536: - realmats - all coefficient matrices of P(.) are real
538: Options Database Keys:
539: + -pep_ciss_integration_points - Sets the number of integration points
540: . -pep_ciss_blocksize - Sets the block size
541: . -pep_ciss_moments - Sets the moment size
542: . -pep_ciss_partitions - Sets the number of partitions
543: . -pep_ciss_maxblocksize - Sets the maximum block size
544: - -pep_ciss_realmats - all coefficient matrices of P(.) are real
546: Notes:
547: For all integer arguments, you can use PETSC_CURRENT to keep the current value, and
548: PETSC_DETERMINE to set them to a default value.
550: The default number of partitions is 1. This means the internal KSP object is shared
551: among all processes of the PEP communicator. Otherwise, the communicator is split
552: into npart communicators, so that npart KSP solves proceed simultaneously.
554: Level: advanced
556: .seealso: PEPCISSGetSizes()
557: @*/
558: PetscErrorCode PEPCISSSetSizes(PEP pep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
559: {
560: PetscFunctionBegin;
568: PetscTryMethod(pep,"PEPCISSSetSizes_C",(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(pep,ip,bs,ms,npart,bsmax,realmats));
569: PetscFunctionReturn(PETSC_SUCCESS);
570: }
572: static PetscErrorCode PEPCISSGetSizes_CISS(PEP pep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
573: {
574: PEP_CISS *ctx = (PEP_CISS*)pep->data;
576: PetscFunctionBegin;
577: if (ip) *ip = ctx->N;
578: if (bs) *bs = ctx->L;
579: if (ms) *ms = ctx->M;
580: if (npart) *npart = ctx->npart;
581: if (bsmax) *bsmax = ctx->L_max;
582: if (realmats) *realmats = ctx->isreal;
583: PetscFunctionReturn(PETSC_SUCCESS);
584: }
586: /*@
587: PEPCISSGetSizes - Gets the values of various size parameters in the CISS solver.
589: Not Collective
591: Input Parameter:
592: . pep - the polynomial eigensolver context
594: Output Parameters:
595: + ip - number of integration points
596: . bs - block size
597: . ms - moment size
598: . npart - number of partitions when splitting the communicator
599: . bsmax - max block size
600: - realmats - all coefficient matrices of P(.) are real
602: Level: advanced
604: .seealso: PEPCISSSetSizes()
605: @*/
606: PetscErrorCode PEPCISSGetSizes(PEP pep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
607: {
608: PetscFunctionBegin;
610: PetscUseMethod(pep,"PEPCISSGetSizes_C",(PEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(pep,ip,bs,ms,npart,bsmax,realmats));
611: PetscFunctionReturn(PETSC_SUCCESS);
612: }
614: static PetscErrorCode PEPCISSSetThreshold_CISS(PEP pep,PetscReal delta,PetscReal spur)
615: {
616: PEP_CISS *ctx = (PEP_CISS*)pep->data;
618: PetscFunctionBegin;
619: if (delta == (PetscReal)PETSC_DETERMINE) {
620: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
621: } else if (delta != (PetscReal)PETSC_CURRENT) {
622: PetscCheck(delta>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
623: ctx->delta = delta;
624: }
625: if (spur == (PetscReal)PETSC_DETERMINE) {
626: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
627: } else if (spur != (PetscReal)PETSC_CURRENT) {
628: PetscCheck(spur>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
629: ctx->spurious_threshold = spur;
630: }
631: PetscFunctionReturn(PETSC_SUCCESS);
632: }
634: /*@
635: PEPCISSSetThreshold - Sets the values of various threshold parameters in
636: the CISS solver.
638: Logically Collective
640: Input Parameters:
641: + pep - the polynomial eigensolver context
642: . delta - threshold for numerical rank
643: - spur - spurious threshold (to discard spurious eigenpairs)
645: Options Database Keys:
646: + -pep_ciss_delta - Sets the delta
647: - -pep_ciss_spurious_threshold - Sets the spurious threshold
649: Note:
650: PETSC_CURRENT can be used to preserve the current value of any of the
651: arguments, and PETSC_DETERMINE to set them to a default value.
653: Level: advanced
655: .seealso: PEPCISSGetThreshold()
656: @*/
657: PetscErrorCode PEPCISSSetThreshold(PEP pep,PetscReal delta,PetscReal spur)
658: {
659: PetscFunctionBegin;
663: PetscTryMethod(pep,"PEPCISSSetThreshold_C",(PEP,PetscReal,PetscReal),(pep,delta,spur));
664: PetscFunctionReturn(PETSC_SUCCESS);
665: }
667: static PetscErrorCode PEPCISSGetThreshold_CISS(PEP pep,PetscReal *delta,PetscReal *spur)
668: {
669: PEP_CISS *ctx = (PEP_CISS*)pep->data;
671: PetscFunctionBegin;
672: if (delta) *delta = ctx->delta;
673: if (spur) *spur = ctx->spurious_threshold;
674: PetscFunctionReturn(PETSC_SUCCESS);
675: }
677: /*@
678: PEPCISSGetThreshold - Gets the values of various threshold parameters in
679: the CISS solver.
681: Not Collective
683: Input Parameter:
684: . pep - the polynomial eigensolver context
686: Output Parameters:
687: + delta - threshold for numerical rank
688: - spur - spurious threshold (to discard spurious eigenpairs)
690: Level: advanced
692: .seealso: PEPCISSSetThreshold()
693: @*/
694: PetscErrorCode PEPCISSGetThreshold(PEP pep,PetscReal *delta,PetscReal *spur)
695: {
696: PetscFunctionBegin;
698: PetscUseMethod(pep,"PEPCISSGetThreshold_C",(PEP,PetscReal*,PetscReal*),(pep,delta,spur));
699: PetscFunctionReturn(PETSC_SUCCESS);
700: }
702: static PetscErrorCode PEPCISSSetRefinement_CISS(PEP pep,PetscInt inner,PetscInt blsize)
703: {
704: PEP_CISS *ctx = (PEP_CISS*)pep->data;
706: PetscFunctionBegin;
707: if (inner == PETSC_DETERMINE) {
708: ctx->refine_inner = 0;
709: } else if (inner != PETSC_CURRENT) {
710: PetscCheck(inner>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
711: ctx->refine_inner = inner;
712: }
713: if (blsize == PETSC_DETERMINE) {
714: ctx->refine_blocksize = 0;
715: } else if (blsize != PETSC_CURRENT) {
716: PetscCheck(blsize>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
717: ctx->refine_blocksize = blsize;
718: }
719: PetscFunctionReturn(PETSC_SUCCESS);
720: }
722: /*@
723: PEPCISSSetRefinement - Sets the values of various refinement parameters
724: in the CISS solver.
726: Logically Collective
728: Input Parameters:
729: + pep - the polynomial eigensolver context
730: . inner - number of iterative refinement iterations (inner loop)
731: - blsize - number of iterative refinement iterations (blocksize loop)
733: Options Database Keys:
734: + -pep_ciss_refine_inner - Sets number of inner iterations
735: - -pep_ciss_refine_blocksize - Sets number of blocksize iterations
737: Note:
738: PETSC_CURRENT can be used to preserve the current value of any of the
739: arguments, and PETSC_DETERMINE to set them to a default value.
741: Level: advanced
743: .seealso: PEPCISSGetRefinement()
744: @*/
745: PetscErrorCode PEPCISSSetRefinement(PEP pep,PetscInt inner,PetscInt blsize)
746: {
747: PetscFunctionBegin;
751: PetscTryMethod(pep,"PEPCISSSetRefinement_C",(PEP,PetscInt,PetscInt),(pep,inner,blsize));
752: PetscFunctionReturn(PETSC_SUCCESS);
753: }
755: static PetscErrorCode PEPCISSGetRefinement_CISS(PEP pep,PetscInt *inner,PetscInt *blsize)
756: {
757: PEP_CISS *ctx = (PEP_CISS*)pep->data;
759: PetscFunctionBegin;
760: if (inner) *inner = ctx->refine_inner;
761: if (blsize) *blsize = ctx->refine_blocksize;
762: PetscFunctionReturn(PETSC_SUCCESS);
763: }
765: /*@
766: PEPCISSGetRefinement - Gets the values of various refinement parameters
767: in the CISS solver.
769: Not Collective
771: Input Parameter:
772: . pep - the polynomial eigensolver context
774: Output Parameters:
775: + inner - number of iterative refinement iterations (inner loop)
776: - blsize - number of iterative refinement iterations (blocksize loop)
778: Level: advanced
780: .seealso: PEPCISSSetRefinement()
781: @*/
782: PetscErrorCode PEPCISSGetRefinement(PEP pep, PetscInt *inner, PetscInt *blsize)
783: {
784: PetscFunctionBegin;
786: PetscUseMethod(pep,"PEPCISSGetRefinement_C",(PEP,PetscInt*,PetscInt*),(pep,inner,blsize));
787: PetscFunctionReturn(PETSC_SUCCESS);
788: }
790: static PetscErrorCode PEPCISSSetExtraction_CISS(PEP pep,PEPCISSExtraction extraction)
791: {
792: PEP_CISS *ctx = (PEP_CISS*)pep->data;
794: PetscFunctionBegin;
795: if (ctx->extraction != extraction) {
796: ctx->extraction = extraction;
797: pep->state = PEP_STATE_INITIAL;
798: }
799: PetscFunctionReturn(PETSC_SUCCESS);
800: }
802: /*@
803: PEPCISSSetExtraction - Sets the extraction technique used in the CISS solver.
805: Logically Collective
807: Input Parameters:
808: + pep - the polynomial eigensolver context
809: - extraction - the extraction technique
811: Options Database Key:
812: . -pep_ciss_extraction - Sets the extraction technique (either 'ritz', 'hankel' or 'caa')
814: Notes:
815: By default, the Rayleigh-Ritz extraction is used (PEP_CISS_EXTRACTION_RITZ).
817: If the 'hankel' or the 'caa' option is specified (PEP_CISS_EXTRACTION_HANKEL or
818: PEP_CISS_EXTRACTION_CAA), then the Block Hankel method, or the Communication-avoiding
819: Arnoldi method, respectively, is used for extracting eigenpairs.
821: Level: advanced
823: .seealso: PEPCISSGetExtraction(), PEPCISSExtraction
824: @*/
825: PetscErrorCode PEPCISSSetExtraction(PEP pep,PEPCISSExtraction extraction)
826: {
827: PetscFunctionBegin;
830: PetscTryMethod(pep,"PEPCISSSetExtraction_C",(PEP,PEPCISSExtraction),(pep,extraction));
831: PetscFunctionReturn(PETSC_SUCCESS);
832: }
834: static PetscErrorCode PEPCISSGetExtraction_CISS(PEP pep,PEPCISSExtraction *extraction)
835: {
836: PEP_CISS *ctx = (PEP_CISS*)pep->data;
838: PetscFunctionBegin;
839: *extraction = ctx->extraction;
840: PetscFunctionReturn(PETSC_SUCCESS);
841: }
843: /*@
844: PEPCISSGetExtraction - Gets the extraction technique used in the CISS solver.
846: Not Collective
848: Input Parameter:
849: . pep - the polynomial eigensolver context
851: Output Parameters:
852: . extraction - extraction technique
854: Level: advanced
856: .seealso: PEPCISSSetExtraction() PEPCISSExtraction
857: @*/
858: PetscErrorCode PEPCISSGetExtraction(PEP pep,PEPCISSExtraction *extraction)
859: {
860: PetscFunctionBegin;
862: PetscAssertPointer(extraction,2);
863: PetscUseMethod(pep,"PEPCISSGetExtraction_C",(PEP,PEPCISSExtraction*),(pep,extraction));
864: PetscFunctionReturn(PETSC_SUCCESS);
865: }
867: static PetscErrorCode PEPCISSGetKSPs_CISS(PEP pep,PetscInt *nsolve,KSP **ksp)
868: {
869: PEP_CISS *ctx = (PEP_CISS*)pep->data;
870: SlepcContourData contour;
871: PetscInt i,nsplit;
872: PC pc;
873: MPI_Comm child;
875: PetscFunctionBegin;
876: if (!ctx->contour) { /* initialize contour data structure first */
877: PetscCall(RGCanUseConjugates(pep->rg,ctx->isreal,&ctx->useconj));
878: PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)pep,&ctx->contour));
879: }
880: contour = ctx->contour;
881: if (!contour->ksp) {
882: PetscCall(PetscMalloc1(contour->npoints,&contour->ksp));
883: PetscCall(PEPGetST(pep,&pep->st));
884: PetscCall(STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL));
885: PetscCall(PetscSubcommGetChild(contour->subcomm,&child));
886: for (i=0;i<contour->npoints;i++) {
887: PetscCall(KSPCreate(child,&contour->ksp[i]));
888: PetscCall(PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)pep,1));
889: PetscCall(KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)pep)->prefix));
890: PetscCall(KSPAppendOptionsPrefix(contour->ksp[i],"pep_ciss_"));
891: PetscCall(PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)pep)->options));
892: PetscCall(KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE));
893: PetscCall(KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(pep->tol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
894: PetscCall(KSPGetPC(contour->ksp[i],&pc));
895: if (nsplit) {
896: PetscCall(KSPSetType(contour->ksp[i],KSPBCGS));
897: PetscCall(PCSetType(pc,PCBJACOBI));
898: } else {
899: PetscCall(KSPSetType(contour->ksp[i],KSPPREONLY));
900: PetscCall(PCSetType(pc,PCLU));
901: }
902: }
903: }
904: if (nsolve) *nsolve = contour->npoints;
905: if (ksp) *ksp = contour->ksp;
906: PetscFunctionReturn(PETSC_SUCCESS);
907: }
909: /*@C
910: PEPCISSGetKSPs - Retrieve the array of linear solver objects associated with
911: the CISS solver.
913: Collective
915: Input Parameter:
916: . pep - polynomial eigenvalue solver
918: Output Parameters:
919: + nsolve - number of solver objects
920: - ksp - array of linear solver object
922: Notes:
923: The number of KSP solvers is equal to the number of integration points divided by
924: the number of partitions. This value is halved in the case of real matrices with
925: a region centered at the real axis.
927: Level: advanced
929: .seealso: PEPCISSSetSizes()
930: @*/
931: PetscErrorCode PEPCISSGetKSPs(PEP pep,PetscInt *nsolve,KSP **ksp)
932: {
933: PetscFunctionBegin;
935: PetscUseMethod(pep,"PEPCISSGetKSPs_C",(PEP,PetscInt*,KSP**),(pep,nsolve,ksp));
936: PetscFunctionReturn(PETSC_SUCCESS);
937: }
939: static PetscErrorCode PEPReset_CISS(PEP pep)
940: {
941: PEP_CISS *ctx = (PEP_CISS*)pep->data;
943: PetscFunctionBegin;
944: PetscCall(BVDestroy(&ctx->S));
945: PetscCall(BVDestroy(&ctx->V));
946: PetscCall(BVDestroy(&ctx->Y));
947: PetscCall(SlepcContourDataReset(ctx->contour));
948: PetscCall(MatDestroy(&ctx->J));
949: PetscCall(BVDestroy(&ctx->pV));
950: PetscCall(PetscFree(ctx->Psplit));
951: PetscFunctionReturn(PETSC_SUCCESS);
952: }
954: static PetscErrorCode PEPSetFromOptions_CISS(PEP pep,PetscOptionItems *PetscOptionsObject)
955: {
956: PEP_CISS *ctx = (PEP_CISS*)pep->data;
957: PetscReal r1,r2;
958: PetscInt i,i1,i2,i3,i4,i5,i6,i7;
959: PetscBool b1,flg,flg2,flg3,flg4,flg5,flg6;
960: PEPCISSExtraction extraction;
962: PetscFunctionBegin;
963: PetscOptionsHeadBegin(PetscOptionsObject,"PEP CISS Options");
965: PetscCall(PEPCISSGetSizes(pep,&i1,&i2,&i3,&i4,&i5,&b1));
966: PetscCall(PetscOptionsInt("-pep_ciss_integration_points","Number of integration points","PEPCISSSetSizes",i1,&i1,&flg));
967: PetscCall(PetscOptionsInt("-pep_ciss_blocksize","Block size","PEPCISSSetSizes",i2,&i2,&flg2));
968: PetscCall(PetscOptionsInt("-pep_ciss_moments","Moment size","PEPCISSSetSizes",i3,&i3,&flg3));
969: PetscCall(PetscOptionsInt("-pep_ciss_partitions","Number of partitions","PEPCISSSetSizes",i4,&i4,&flg4));
970: PetscCall(PetscOptionsInt("-pep_ciss_maxblocksize","Maximum block size","PEPCISSSetSizes",i5,&i5,&flg5));
971: PetscCall(PetscOptionsBool("-pep_ciss_realmats","True if all coefficient matrices of P(.) are real","PEPCISSSetSizes",b1,&b1,&flg6));
972: if (flg || flg2 || flg3 || flg4 || flg5 || flg6) PetscCall(PEPCISSSetSizes(pep,i1,i2,i3,i4,i5,b1));
974: PetscCall(PEPCISSGetThreshold(pep,&r1,&r2));
975: PetscCall(PetscOptionsReal("-pep_ciss_delta","Threshold for numerical rank","PEPCISSSetThreshold",r1,&r1,&flg));
976: PetscCall(PetscOptionsReal("-pep_ciss_spurious_threshold","Threshold for the spurious eigenpairs","PEPCISSSetThreshold",r2,&r2,&flg2));
977: if (flg || flg2) PetscCall(PEPCISSSetThreshold(pep,r1,r2));
979: PetscCall(PEPCISSGetRefinement(pep,&i6,&i7));
980: PetscCall(PetscOptionsInt("-pep_ciss_refine_inner","Number of inner iterative refinement iterations","PEPCISSSetRefinement",i6,&i6,&flg));
981: PetscCall(PetscOptionsInt("-pep_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","PEPCISSSetRefinement",i7,&i7,&flg2));
982: if (flg || flg2) PetscCall(PEPCISSSetRefinement(pep,i6,i7));
984: PetscCall(PetscOptionsEnum("-pep_ciss_extraction","Extraction technique","PEPCISSSetExtraction",PEPCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg));
985: if (flg) PetscCall(PEPCISSSetExtraction(pep,extraction));
987: PetscOptionsHeadEnd();
989: if (!pep->rg) PetscCall(PEPGetRG(pep,&pep->rg));
990: PetscCall(RGSetFromOptions(pep->rg)); /* this is necessary here to set useconj */
991: if (!ctx->contour || !ctx->contour->ksp) PetscCall(PEPCISSGetKSPs(pep,NULL,NULL));
992: PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Something went wrong with PEPCISSGetKSPs()");
993: for (i=0;i<ctx->contour->npoints;i++) PetscCall(KSPSetFromOptions(ctx->contour->ksp[i]));
994: PetscCall(PetscSubcommSetFromOptions(ctx->contour->subcomm));
995: PetscFunctionReturn(PETSC_SUCCESS);
996: }
998: static PetscErrorCode PEPDestroy_CISS(PEP pep)
999: {
1000: PEP_CISS *ctx = (PEP_CISS*)pep->data;
1002: PetscFunctionBegin;
1003: PetscCall(SlepcContourDataDestroy(&ctx->contour));
1004: PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
1005: PetscCall(PetscFree(pep->data));
1006: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetSizes_C",NULL));
1007: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetSizes_C",NULL));
1008: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetThreshold_C",NULL));
1009: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetThreshold_C",NULL));
1010: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetRefinement_C",NULL));
1011: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetRefinement_C",NULL));
1012: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetExtraction_C",NULL));
1013: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetExtraction_C",NULL));
1014: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetKSPs_C",NULL));
1015: PetscFunctionReturn(PETSC_SUCCESS);
1016: }
1018: static PetscErrorCode PEPView_CISS(PEP pep,PetscViewer viewer)
1019: {
1020: PEP_CISS *ctx = (PEP_CISS*)pep->data;
1021: PetscBool isascii;
1022: PetscViewer sviewer;
1024: PetscFunctionBegin;
1025: PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
1026: if (isascii) {
1027: PetscCall(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));
1028: if (ctx->isreal) PetscCall(PetscViewerASCIIPrintf(viewer," exploiting symmetry of integration points\n"));
1029: PetscCall(PetscViewerASCIIPrintf(viewer," threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold));
1030: PetscCall(PetscViewerASCIIPrintf(viewer," iterative refinement { inner: %" PetscInt_FMT ", blocksize: %" PetscInt_FMT " }\n",ctx->refine_inner, ctx->refine_blocksize));
1031: PetscCall(PetscViewerASCIIPrintf(viewer," extraction: %s\n",PEPCISSExtractions[ctx->extraction]));
1032: if (!ctx->contour || !ctx->contour->ksp) PetscCall(PEPCISSGetKSPs(pep,NULL,NULL));
1033: PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Something went wrong with PEPCISSGetKSPs()");
1034: PetscCall(PetscViewerASCIIPushTab(viewer));
1035: if (ctx->npart>1 && ctx->contour->subcomm) {
1036: PetscCall(PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
1037: if (!ctx->contour->subcomm->color) PetscCall(KSPView(ctx->contour->ksp[0],sviewer));
1038: PetscCall(PetscViewerFlush(sviewer));
1039: PetscCall(PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
1040: /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1041: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1042: } else PetscCall(KSPView(ctx->contour->ksp[0],viewer));
1043: PetscCall(PetscViewerASCIIPopTab(viewer));
1044: }
1045: PetscFunctionReturn(PETSC_SUCCESS);
1046: }
1048: static PetscErrorCode PEPSetDSType_CISS(PEP pep)
1049: {
1050: PEP_CISS *ctx = (PEP_CISS*)pep->data;
1052: PetscFunctionBegin;
1053: switch (ctx->extraction) {
1054: case PEP_CISS_EXTRACTION_RITZ: PetscCall(DSSetType(pep->ds,DSPEP)); break;
1055: case PEP_CISS_EXTRACTION_HANKEL: PetscCall(DSSetType(pep->ds,DSGNHEP)); break;
1056: case PEP_CISS_EXTRACTION_CAA: PetscCall(DSSetType(pep->ds,DSNHEP)); break;
1057: }
1058: PetscFunctionReturn(PETSC_SUCCESS);
1059: }
1061: SLEPC_EXTERN PetscErrorCode PEPCreate_CISS(PEP pep)
1062: {
1063: PEP_CISS *ctx = (PEP_CISS*)pep->data;
1065: PetscFunctionBegin;
1066: PetscCall(PetscNew(&ctx));
1067: pep->data = ctx;
1068: /* set default values of parameters */
1069: ctx->N = 32;
1070: ctx->L = 16;
1071: ctx->M = ctx->N/4;
1072: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
1073: ctx->L_max = 64;
1074: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1075: ctx->isreal = PETSC_FALSE;
1076: ctx->npart = 1;
1078: pep->ops->solve = PEPSolve_CISS;
1079: pep->ops->setup = PEPSetUp_CISS;
1080: pep->ops->setfromoptions = PEPSetFromOptions_CISS;
1081: pep->ops->reset = PEPReset_CISS;
1082: pep->ops->destroy = PEPDestroy_CISS;
1083: pep->ops->view = PEPView_CISS;
1084: pep->ops->setdstype = PEPSetDSType_CISS;
1086: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetSizes_C",PEPCISSSetSizes_CISS));
1087: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetSizes_C",PEPCISSGetSizes_CISS));
1088: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetThreshold_C",PEPCISSSetThreshold_CISS));
1089: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetThreshold_C",PEPCISSGetThreshold_CISS));
1090: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetRefinement_C",PEPCISSSetRefinement_CISS));
1091: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetRefinement_C",PEPCISSGetRefinement_CISS));
1092: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetExtraction_C",PEPCISSSetExtraction_CISS));
1093: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetExtraction_C",PEPCISSGetExtraction_CISS));
1094: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetKSPs_C",PEPCISSGetKSPs_CISS));
1095: PetscFunctionReturn(PETSC_SUCCESS);
1096: }