Actual source code: pciss.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] 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: A = (contour->pA)?contour->pA:pep->A;
72: K = (contour->pP)?contour->pP:ctx->Psplit;
73: PetscMalloc1(pep->nmat,&coeff);
74: if (deriv) PEPEvaluateBasisDerivative(pep,lambda,0,coeff,NULL);
75: else PEPEvaluateBasis(pep,lambda,0,coeff,NULL);
76: STGetMatStructure(pep->st,&str);
77: MatZeroEntries(T);
78: if (!deriv && T != P) {
79: STGetSplitPreconditionerInfo(pep->st,NULL,&strp);
80: MatZeroEntries(P);
81: }
82: i = deriv?1:0;
83: for (;i<pep->nmat;i++) {
84: MatAXPY(T,coeff[i],A[i],str);
85: if (!deriv && T != P) MatAXPY(P,coeff[i],K[i],strp);
86: }
87: PetscFree(coeff);
88: return 0;
89: }
91: /*
92: Set up KSP solvers for every integration point
93: */
94: static PetscErrorCode PEPCISSSetUp(PEP pep,Mat T,Mat P)
95: {
96: PEP_CISS *ctx = (PEP_CISS*)pep->data;
97: SlepcContourData contour;
98: PetscInt i,p_id;
99: Mat Amat,Pmat;
101: if (!ctx->contour || !ctx->contour->ksp) PEPCISSGetKSPs(pep,NULL,NULL);
102: contour = ctx->contour;
103: for (i=0;i<contour->npoints;i++) {
104: p_id = i*contour->subcomm->n + contour->subcomm->color;
105: MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&Amat);
106: if (T != P) MatDuplicate(P,MAT_DO_NOT_COPY_VALUES,&Pmat); else Pmat = Amat;
107: PEPComputeFunction(pep,ctx->omega[p_id],Amat,Pmat,PETSC_FALSE);
108: PEP_KSPSetOperators(contour->ksp[i],Amat,Pmat);
109: MatDestroy(&Amat);
110: if (T != P) MatDestroy(&Pmat);
111: }
112: return 0;
113: }
115: /*
116: Y_i = F(z_i)^{-1}Fp(z_i)V for every integration point, Y=[Y_i] is in the context
117: */
118: static PetscErrorCode PEPCISSSolve(PEP pep,Mat dT,BV V,PetscInt L_start,PetscInt L_end)
119: {
120: PEP_CISS *ctx = (PEP_CISS*)pep->data;
121: SlepcContourData contour;
122: PetscInt i,p_id;
123: Mat MV,BMV=NULL,MC;
125: contour = ctx->contour;
126: BVSetActiveColumns(V,L_start,L_end);
127: BVGetMat(V,&MV);
128: for (i=0;i<contour->npoints;i++) {
129: p_id = i*contour->subcomm->n + contour->subcomm->color;
130: PEPComputeFunction(pep,ctx->omega[p_id],dT,NULL,PETSC_TRUE);
131: BVSetActiveColumns(ctx->Y,i*ctx->L+L_start,i*ctx->L+L_end);
132: BVGetMat(ctx->Y,&MC);
133: if (!i) {
134: MatProductCreate(dT,MV,NULL,&BMV);
135: MatProductSetType(BMV,MATPRODUCT_AB);
136: MatProductSetFromOptions(BMV);
137: MatProductSymbolic(BMV);
138: }
139: MatProductNumeric(BMV);
140: KSPMatSolve(contour->ksp[i],BMV,MC);
141: BVRestoreMat(ctx->Y,&MC);
142: }
143: MatDestroy(&BMV);
144: BVRestoreMat(V,&MV);
145: return 0;
146: }
148: PetscErrorCode PEPSetUp_CISS(PEP pep)
149: {
150: PEP_CISS *ctx = (PEP_CISS*)pep->data;
151: SlepcContourData contour;
152: PetscInt i,nwork,nsplit;
153: PetscBool istrivial,isellipse,flg;
154: PetscObjectId id;
155: PetscObjectState state;
156: Vec v0;
158: if (pep->ncv==PETSC_DEFAULT) pep->ncv = ctx->L_max*ctx->M;
159: else {
160: ctx->L_max = pep->ncv/ctx->M;
161: if (!ctx->L_max) {
162: ctx->L_max = 1;
163: pep->ncv = ctx->L_max*ctx->M;
164: }
165: }
166: ctx->L = PetscMin(ctx->L,ctx->L_max);
167: if (pep->max_it==PETSC_DEFAULT) pep->max_it = 5;
168: if (pep->mpd==PETSC_DEFAULT) pep->mpd = pep->ncv;
169: if (!pep->which) pep->which = PEP_ALL;
171: PEPCheckUnsupported(pep,PEP_FEATURE_STOPPING);
172: PEPCheckIgnored(pep,PEP_FEATURE_SCALE);
174: /* check region */
175: RGIsTrivial(pep->rg,&istrivial);
177: RGGetComplement(pep->rg,&flg);
179: PetscObjectTypeCompare((PetscObject)pep->rg,RGELLIPSE,&isellipse);
182: /* if the region has changed, then reset contour data */
183: PetscObjectGetId((PetscObject)pep->rg,&id);
184: PetscObjectStateGet((PetscObject)pep->rg,&state);
185: if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
186: SlepcContourDataDestroy(&ctx->contour);
187: PetscInfo(pep,"Resetting the contour data structure due to a change of region\n");
188: ctx->rgid = id; ctx->rgstate = state;
189: }
191: /* create contour data structure */
192: if (!ctx->contour) {
193: RGCanUseConjugates(pep->rg,ctx->isreal,&ctx->useconj);
194: SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)pep,&ctx->contour);
195: }
197: PEPAllocateSolution(pep,0);
198: if (ctx->weight) PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
199: PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma);
201: /* allocate basis vectors */
202: BVDestroy(&ctx->S);
203: BVDuplicateResize(pep->V,ctx->L*ctx->M,&ctx->S);
204: BVDestroy(&ctx->V);
205: BVDuplicateResize(pep->V,ctx->L,&ctx->V);
207: /* check if a user-defined split preconditioner has been set */
208: STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL);
209: if (nsplit) {
210: PetscFree(ctx->Psplit);
211: PetscMalloc1(nsplit,&ctx->Psplit);
212: for (i=0;i<nsplit;i++) STGetSplitPreconditionerTerm(pep->st,i,&ctx->Psplit[i]);
213: }
215: contour = ctx->contour;
216: SlepcContourRedundantMat(contour,pep->nmat,pep->A,ctx->Psplit);
217: if (!ctx->J) MatDuplicate(contour->pA?contour->pA[0]:pep->A[0],MAT_DO_NOT_COPY_VALUES,&ctx->J);
218: if (contour->pA) {
219: BVGetColumn(ctx->V,0,&v0);
220: SlepcContourScatterCreate(contour,v0);
221: BVRestoreColumn(ctx->V,0,&v0);
222: BVDestroy(&ctx->pV);
223: BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV);
224: BVSetSizesFromVec(ctx->pV,contour->xsub,pep->n);
225: BVSetFromOptions(ctx->pV);
226: BVResize(ctx->pV,ctx->L,PETSC_FALSE);
227: }
229: BVDestroy(&ctx->Y);
230: if (contour->pA) {
231: BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y);
232: BVSetSizesFromVec(ctx->Y,contour->xsub,pep->n);
233: BVSetFromOptions(ctx->Y);
234: BVResize(ctx->Y,contour->npoints*ctx->L,PETSC_FALSE);
235: } else BVDuplicateResize(pep->V,contour->npoints*ctx->L,&ctx->Y);
237: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) DSSetType(pep->ds,DSGNHEP);
238: else if (ctx->extraction == PEP_CISS_EXTRACTION_CAA) DSSetType(pep->ds,DSNHEP);
239: else {
240: DSSetType(pep->ds,DSPEP);
241: DSPEPSetDegree(pep->ds,pep->nmat-1);
242: DSPEPSetCoefficients(pep->ds,pep->pbc);
243: }
244: DSAllocate(pep->ds,pep->ncv);
245: nwork = 2;
246: PEPSetWorkVecs(pep,nwork);
247: return 0;
248: }
250: PetscErrorCode PEPSolve_CISS(PEP pep)
251: {
252: PEP_CISS *ctx = (PEP_CISS*)pep->data;
253: SlepcContourData contour = ctx->contour;
254: Mat X,M,E,T,P;
255: PetscInt i,j,ld,L_add=0,nv=0,L_base=ctx->L,inner,*inside,nsplit;
256: PetscScalar *Mu,*H0,*H1,*rr,*temp,center;
257: PetscReal error,max_error,radius,rgscale,est_eig,eta;
258: PetscBool isellipse,*fl1;
259: Vec si;
260: SlepcSC sc;
261: PetscRandom rand;
263: DSSetFromOptions(pep->ds);
264: DSGetSlepcSC(pep->ds,&sc);
265: sc->comparison = SlepcCompareLargestMagnitude;
266: sc->comparisonctx = NULL;
267: sc->map = NULL;
268: sc->mapobj = NULL;
269: DSGetLeadingDimension(pep->ds,&ld);
270: RGComputeQuadrature(pep->rg,RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight);
271: STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL);
272: if (contour->pA) {
273: T = contour->pA[0];
274: P = nsplit? contour->pP[0]: T;
275: } else {
276: T = pep->A[0];
277: P = nsplit? ctx->Psplit[0]: T;
278: }
279: PEPCISSSetUp(pep,T,P);
280: BVSetActiveColumns(ctx->V,0,ctx->L);
281: BVSetRandomSign(ctx->V);
282: BVGetRandomContext(ctx->V,&rand);
283: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
284: PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L);
285: PetscObjectTypeCompare((PetscObject)pep->rg,RGELLIPSE,&isellipse);
286: if (isellipse) {
287: BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig);
288: PetscInfo(pep,"Estimated eigenvalue count: %f\n",(double)est_eig);
289: eta = PetscPowReal(10.0,-PetscLog10Real(pep->tol)/ctx->N);
290: L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
291: if (L_add>ctx->L_max-ctx->L) {
292: PetscInfo(pep,"Number of eigenvalues inside the contour path may be too large\n");
293: L_add = ctx->L_max-ctx->L;
294: }
295: }
296: /* Updates L after estimate the number of eigenvalue */
297: if (L_add>0) {
298: PetscInfo(pep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by Estimate #Eig\n",ctx->L,ctx->L+L_add);
299: BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints);
300: BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
301: BVSetRandomSign(ctx->V);
302: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
303: ctx->L += L_add;
304: PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,ctx->L-L_add,ctx->L);
305: }
307: PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
308: for (i=0;i<ctx->refine_blocksize;i++) {
309: 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);
310: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
311: PetscLogEventBegin(PEP_CISS_SVD,pep,0,0,0);
312: SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
313: PetscLogEventEnd(PEP_CISS_SVD,pep,0,0,0);
314: if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
315: L_add = L_base;
316: if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
317: PetscInfo(pep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by SVD(H0)\n",ctx->L,ctx->L+L_add);
318: BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints);
319: BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
320: BVSetRandomSign(ctx->V);
321: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
322: ctx->L += L_add;
323: PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,ctx->L-L_add,ctx->L);
324: if (L_add) {
325: PetscFree2(Mu,H0);
326: PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
327: }
328: }
330: RGGetScale(pep->rg,&rgscale);
331: RGEllipseGetParameters(pep->rg,¢er,&radius,NULL);
333: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1);
335: while (pep->reason == PEP_CONVERGED_ITERATING) {
336: pep->its++;
337: for (inner=0;inner<=ctx->refine_inner;inner++) {
338: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
339: 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);
340: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
341: PetscLogEventBegin(PEP_CISS_SVD,pep,0,0,0);
342: SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
343: PetscLogEventEnd(PEP_CISS_SVD,pep,0,0,0);
344: } else {
345: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
346: /* compute SVD of S */
347: 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);
348: }
349: PetscInfo(pep,"Estimated rank: nv = %" PetscInt_FMT "\n",nv);
350: if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
351: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
352: BVSetActiveColumns(ctx->S,0,ctx->L);
353: BVSetActiveColumns(ctx->V,0,ctx->L);
354: BVCopy(ctx->S,ctx->V);
355: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
356: PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L);
357: } else break;
358: }
359: pep->nconv = 0;
360: if (nv == 0) { pep->reason = PEP_CONVERGED_TOL; break; }
361: else {
362: /* Extracting eigenpairs */
363: DSSetDimensions(pep->ds,nv,0,0);
364: DSSetState(pep->ds,DS_STATE_RAW);
365: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
366: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
367: CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1);
368: DSGetArray(pep->ds,DS_MAT_A,&temp);
369: for (j=0;j<nv;j++)
370: for (i=0;i<nv;i++)
371: temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
372: DSRestoreArray(pep->ds,DS_MAT_A,&temp);
373: DSGetArray(pep->ds,DS_MAT_B,&temp);
374: for (j=0;j<nv;j++)
375: for (i=0;i<nv;i++)
376: temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
377: DSRestoreArray(pep->ds,DS_MAT_B,&temp);
378: } else if (ctx->extraction == PEP_CISS_EXTRACTION_CAA) {
379: BVSetActiveColumns(ctx->S,0,nv);
380: DSGetArray(pep->ds,DS_MAT_A,&temp);
381: for (i=0;i<nv;i++) PetscArraycpy(temp+i*ld,H0+i*nv,nv);
382: DSRestoreArray(pep->ds,DS_MAT_A,&temp);
383: } else {
384: BVSetActiveColumns(ctx->S,0,nv);
385: for (i=0;i<pep->nmat;i++) {
386: DSGetMat(pep->ds,DSMatExtra[i],&E);
387: BVMatProject(ctx->S,pep->A[i],ctx->S,E);
388: DSRestoreMat(pep->ds,DSMatExtra[i],&E);
389: }
390: nv = (pep->nmat-1)*nv;
391: }
392: DSSolve(pep->ds,pep->eigr,pep->eigi);
393: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
394: if (ctx->extraction == PEP_CISS_EXTRACTION_CAA || ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
395: for (i=0;i<nv;i++) {
396: pep->eigr[i] = (pep->eigr[i]*radius+center)*rgscale;
397: }
398: }
399: PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr);
400: DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
401: DSGetMat(pep->ds,DS_MAT_X,&X);
402: SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1);
403: DSRestoreMat(pep->ds,DS_MAT_X,&X);
404: RGCheckInside(pep->rg,nv,pep->eigr,pep->eigi,inside);
405: for (i=0;i<nv;i++) {
406: if (fl1[i] && inside[i]>=0) {
407: rr[i] = 1.0;
408: pep->nconv++;
409: } else rr[i] = 0.0;
410: }
411: DSSort(pep->ds,pep->eigr,pep->eigi,rr,NULL,&pep->nconv);
412: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
413: if (ctx->extraction == PEP_CISS_EXTRACTION_CAA || ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
414: for (i=0;i<nv;i++) pep->eigr[i] = (pep->eigr[i]*radius+center)*rgscale;
415: }
416: PetscFree3(fl1,inside,rr);
417: BVSetActiveColumns(pep->V,0,nv);
418: DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
419: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
420: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
421: BVSetActiveColumns(ctx->S,0,nv);
422: BVCopy(ctx->S,pep->V);
423: DSGetMat(pep->ds,DS_MAT_X,&X);
424: BVMultInPlace(ctx->S,X,0,pep->nconv);
425: BVMultInPlace(pep->V,X,0,pep->nconv);
426: DSRestoreMat(pep->ds,DS_MAT_X,&X);
427: } else {
428: DSGetMat(pep->ds,DS_MAT_X,&X);
429: BVMultInPlace(ctx->S,X,0,pep->nconv);
430: DSRestoreMat(pep->ds,DS_MAT_X,&X);
431: BVSetActiveColumns(ctx->S,0,pep->nconv);
432: BVCopy(ctx->S,pep->V);
433: }
434: max_error = 0.0;
435: for (i=0;i<pep->nconv;i++) {
436: BVGetColumn(pep->V,i,&si);
437: VecNormalize(si,NULL);
438: PEPComputeResidualNorm_Private(pep,pep->eigr[i],0,si,NULL,pep->work,&error);
439: (*pep->converged)(pep,pep->eigr[i],0,error,&error,pep->convergedctx);
440: BVRestoreColumn(pep->V,i,&si);
441: max_error = PetscMax(max_error,error);
442: }
443: if (max_error <= pep->tol) pep->reason = PEP_CONVERGED_TOL;
444: else if (pep->its > pep->max_it) pep->reason = PEP_DIVERGED_ITS;
445: else {
446: if (pep->nconv > ctx->L) nv = pep->nconv;
447: else if (ctx->L > nv) nv = ctx->L;
448: nv = PetscMin(nv,ctx->L*ctx->M);
449: MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M);
450: MatSetRandom(M,rand);
451: BVSetActiveColumns(ctx->S,0,nv);
452: BVMultInPlace(ctx->S,M,0,ctx->L);
453: MatDestroy(&M);
454: BVSetActiveColumns(ctx->S,0,ctx->L);
455: BVSetActiveColumns(ctx->V,0,ctx->L);
456: BVCopy(ctx->S,ctx->V);
457: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
458: PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L);
459: }
460: }
461: }
462: PetscFree2(Mu,H0);
463: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) PetscFree(H1);
464: return 0;
465: }
467: static PetscErrorCode PEPCISSSetSizes_CISS(PEP pep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
468: {
469: PEP_CISS *ctx = (PEP_CISS*)pep->data;
470: PetscInt oN,oL,oM,oLmax,onpart;
471: PetscMPIInt size;
473: oN = ctx->N;
474: if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
475: if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
476: } else {
479: if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
480: }
481: oL = ctx->L;
482: if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
483: ctx->L = 16;
484: } else {
486: ctx->L = bs;
487: }
488: oM = ctx->M;
489: if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
490: ctx->M = ctx->N/4;
491: } else {
494: ctx->M = PetscMax(ms,2);
495: }
496: onpart = ctx->npart;
497: if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
498: ctx->npart = 1;
499: } else {
500: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size);
502: ctx->npart = npart;
503: }
504: oLmax = ctx->L_max;
505: if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
506: ctx->L_max = 64;
507: } else {
509: ctx->L_max = PetscMax(bsmax,ctx->L);
510: }
511: if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
512: SlepcContourDataDestroy(&ctx->contour);
513: PetscInfo(pep,"Resetting the contour data structure due to a change of parameters\n");
514: pep->state = PEP_STATE_INITIAL;
515: }
516: ctx->isreal = realmats;
517: if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) pep->state = PEP_STATE_INITIAL;
518: return 0;
519: }
521: /*@
522: PEPCISSSetSizes - Sets the values of various size parameters in the CISS solver.
524: Logically Collective on pep
526: Input Parameters:
527: + pep - the polynomial eigensolver context
528: . ip - number of integration points
529: . bs - block size
530: . ms - moment size
531: . npart - number of partitions when splitting the communicator
532: . bsmax - max block size
533: - realmats - all coefficient matrices of P(.) are real
535: Options Database Keys:
536: + -pep_ciss_integration_points - Sets the number of integration points
537: . -pep_ciss_blocksize - Sets the block size
538: . -pep_ciss_moments - Sets the moment size
539: . -pep_ciss_partitions - Sets the number of partitions
540: . -pep_ciss_maxblocksize - Sets the maximum block size
541: - -pep_ciss_realmats - all coefficient matrices of P(.) are real
543: Notes:
544: The default number of partitions is 1. This means the internal KSP object is shared
545: among all processes of the PEP communicator. Otherwise, the communicator is split
546: into npart communicators, so that npart KSP solves proceed simultaneously.
548: Level: advanced
550: .seealso: PEPCISSGetSizes()
551: @*/
552: PetscErrorCode PEPCISSSetSizes(PEP pep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
553: {
561: PetscTryMethod(pep,"PEPCISSSetSizes_C",(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(pep,ip,bs,ms,npart,bsmax,realmats));
562: return 0;
563: }
565: static PetscErrorCode PEPCISSGetSizes_CISS(PEP pep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
566: {
567: PEP_CISS *ctx = (PEP_CISS*)pep->data;
569: if (ip) *ip = ctx->N;
570: if (bs) *bs = ctx->L;
571: if (ms) *ms = ctx->M;
572: if (npart) *npart = ctx->npart;
573: if (bsmax) *bsmax = ctx->L_max;
574: if (realmats) *realmats = ctx->isreal;
575: return 0;
576: }
578: /*@
579: PEPCISSGetSizes - Gets the values of various size parameters in the CISS solver.
581: Not Collective
583: Input Parameter:
584: . pep - the polynomial eigensolver context
586: Output Parameters:
587: + ip - number of integration points
588: . bs - block size
589: . ms - moment size
590: . npart - number of partitions when splitting the communicator
591: . bsmax - max block size
592: - realmats - all coefficient matrices of P(.) are real
594: Level: advanced
596: .seealso: PEPCISSSetSizes()
597: @*/
598: PetscErrorCode PEPCISSGetSizes(PEP pep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
599: {
601: PetscUseMethod(pep,"PEPCISSGetSizes_C",(PEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(pep,ip,bs,ms,npart,bsmax,realmats));
602: return 0;
603: }
605: static PetscErrorCode PEPCISSSetThreshold_CISS(PEP pep,PetscReal delta,PetscReal spur)
606: {
607: PEP_CISS *ctx = (PEP_CISS*)pep->data;
609: if (delta == PETSC_DEFAULT) {
610: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
611: } else {
613: ctx->delta = delta;
614: }
615: if (spur == PETSC_DEFAULT) {
616: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
617: } else {
619: ctx->spurious_threshold = spur;
620: }
621: return 0;
622: }
624: /*@
625: PEPCISSSetThreshold - Sets the values of various threshold parameters in
626: the CISS solver.
628: Logically Collective on pep
630: Input Parameters:
631: + pep - the polynomial eigensolver context
632: . delta - threshold for numerical rank
633: - spur - spurious threshold (to discard spurious eigenpairs)
635: Options Database Keys:
636: + -pep_ciss_delta - Sets the delta
637: - -pep_ciss_spurious_threshold - Sets the spurious threshold
639: Level: advanced
641: .seealso: PEPCISSGetThreshold()
642: @*/
643: PetscErrorCode PEPCISSSetThreshold(PEP pep,PetscReal delta,PetscReal spur)
644: {
648: PetscTryMethod(pep,"PEPCISSSetThreshold_C",(PEP,PetscReal,PetscReal),(pep,delta,spur));
649: return 0;
650: }
652: static PetscErrorCode PEPCISSGetThreshold_CISS(PEP pep,PetscReal *delta,PetscReal *spur)
653: {
654: PEP_CISS *ctx = (PEP_CISS*)pep->data;
656: if (delta) *delta = ctx->delta;
657: if (spur) *spur = ctx->spurious_threshold;
658: return 0;
659: }
661: /*@
662: PEPCISSGetThreshold - Gets the values of various threshold parameters in
663: the CISS solver.
665: Not Collective
667: Input Parameter:
668: . pep - the polynomial eigensolver context
670: Output Parameters:
671: + delta - threshold for numerical rank
672: - spur - spurious threshold (to discard spurious eigenpairs)
674: Level: advanced
676: .seealso: PEPCISSSetThreshold()
677: @*/
678: PetscErrorCode PEPCISSGetThreshold(PEP pep,PetscReal *delta,PetscReal *spur)
679: {
681: PetscUseMethod(pep,"PEPCISSGetThreshold_C",(PEP,PetscReal*,PetscReal*),(pep,delta,spur));
682: return 0;
683: }
685: static PetscErrorCode PEPCISSSetRefinement_CISS(PEP pep,PetscInt inner,PetscInt blsize)
686: {
687: PEP_CISS *ctx = (PEP_CISS*)pep->data;
689: if (inner == PETSC_DEFAULT) {
690: ctx->refine_inner = 0;
691: } else {
693: ctx->refine_inner = inner;
694: }
695: if (blsize == PETSC_DEFAULT) {
696: ctx->refine_blocksize = 0;
697: } else {
699: ctx->refine_blocksize = blsize;
700: }
701: return 0;
702: }
704: /*@
705: PEPCISSSetRefinement - Sets the values of various refinement parameters
706: in the CISS solver.
708: Logically Collective on pep
710: Input Parameters:
711: + pep - the polynomial eigensolver context
712: . inner - number of iterative refinement iterations (inner loop)
713: - blsize - number of iterative refinement iterations (blocksize loop)
715: Options Database Keys:
716: + -pep_ciss_refine_inner - Sets number of inner iterations
717: - -pep_ciss_refine_blocksize - Sets number of blocksize iterations
719: Level: advanced
721: .seealso: PEPCISSGetRefinement()
722: @*/
723: PetscErrorCode PEPCISSSetRefinement(PEP pep,PetscInt inner,PetscInt blsize)
724: {
728: PetscTryMethod(pep,"PEPCISSSetRefinement_C",(PEP,PetscInt,PetscInt),(pep,inner,blsize));
729: return 0;
730: }
732: static PetscErrorCode PEPCISSGetRefinement_CISS(PEP pep,PetscInt *inner,PetscInt *blsize)
733: {
734: PEP_CISS *ctx = (PEP_CISS*)pep->data;
736: if (inner) *inner = ctx->refine_inner;
737: if (blsize) *blsize = ctx->refine_blocksize;
738: return 0;
739: }
741: /*@
742: PEPCISSGetRefinement - Gets the values of various refinement parameters
743: in the CISS solver.
745: Not Collective
747: Input Parameter:
748: . pep - the polynomial eigensolver context
750: Output Parameters:
751: + inner - number of iterative refinement iterations (inner loop)
752: - blsize - number of iterative refinement iterations (blocksize loop)
754: Level: advanced
756: .seealso: PEPCISSSetRefinement()
757: @*/
758: PetscErrorCode PEPCISSGetRefinement(PEP pep, PetscInt *inner, PetscInt *blsize)
759: {
761: PetscUseMethod(pep,"PEPCISSGetRefinement_C",(PEP,PetscInt*,PetscInt*),(pep,inner,blsize));
762: return 0;
763: }
765: static PetscErrorCode PEPCISSSetExtraction_CISS(PEP pep,PEPCISSExtraction extraction)
766: {
767: PEP_CISS *ctx = (PEP_CISS*)pep->data;
769: if (ctx->extraction != extraction) {
770: ctx->extraction = extraction;
771: pep->state = PEP_STATE_INITIAL;
772: }
773: return 0;
774: }
776: /*@
777: PEPCISSSetExtraction - Sets the extraction technique used in the CISS solver.
779: Logically Collective on pep
781: Input Parameters:
782: + pep - the polynomial eigensolver context
783: - extraction - the extraction technique
785: Options Database Key:
786: . -pep_ciss_extraction - Sets the extraction technique (either 'ritz', 'hankel' or 'caa')
788: Notes:
789: By default, the Rayleigh-Ritz extraction is used (PEP_CISS_EXTRACTION_RITZ).
791: If the 'hankel' or the 'caa' option is specified (PEP_CISS_EXTRACTION_HANKEL or
792: PEP_CISS_EXTRACTION_CAA), then the Block Hankel method, or the Communication-avoiding
793: Arnoldi method, respectively, is used for extracting eigenpairs.
795: Level: advanced
797: .seealso: PEPCISSGetExtraction(), PEPCISSExtraction
798: @*/
799: PetscErrorCode PEPCISSSetExtraction(PEP pep,PEPCISSExtraction extraction)
800: {
803: PetscTryMethod(pep,"PEPCISSSetExtraction_C",(PEP,PEPCISSExtraction),(pep,extraction));
804: return 0;
805: }
807: static PetscErrorCode PEPCISSGetExtraction_CISS(PEP pep,PEPCISSExtraction *extraction)
808: {
809: PEP_CISS *ctx = (PEP_CISS*)pep->data;
811: *extraction = ctx->extraction;
812: return 0;
813: }
815: /*@
816: PEPCISSGetExtraction - Gets the extraction technique used in the CISS solver.
818: Not Collective
820: Input Parameter:
821: . pep - the polynomial eigensolver context
823: Output Parameters:
824: . extraction - extraction technique
826: Level: advanced
828: .seealso: PEPCISSSetExtraction() PEPCISSExtraction
829: @*/
830: PetscErrorCode PEPCISSGetExtraction(PEP pep,PEPCISSExtraction *extraction)
831: {
834: PetscUseMethod(pep,"PEPCISSGetExtraction_C",(PEP,PEPCISSExtraction*),(pep,extraction));
835: return 0;
836: }
838: static PetscErrorCode PEPCISSGetKSPs_CISS(PEP pep,PetscInt *nsolve,KSP **ksp)
839: {
840: PEP_CISS *ctx = (PEP_CISS*)pep->data;
841: SlepcContourData contour;
842: PetscInt i,nsplit;
843: PC pc;
844: MPI_Comm child;
846: if (!ctx->contour) { /* initialize contour data structure first */
847: RGCanUseConjugates(pep->rg,ctx->isreal,&ctx->useconj);
848: SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)pep,&ctx->contour);
849: }
850: contour = ctx->contour;
851: if (!contour->ksp) {
852: PetscMalloc1(contour->npoints,&contour->ksp);
853: PEPGetST(pep,&pep->st);
854: STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL);
855: PetscSubcommGetChild(contour->subcomm,&child);
856: for (i=0;i<contour->npoints;i++) {
857: KSPCreate(child,&contour->ksp[i]);
858: PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)pep,1);
859: KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)pep)->prefix);
860: KSPAppendOptionsPrefix(contour->ksp[i],"pep_ciss_");
861: PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)pep)->options);
862: KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE);
863: KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(pep->tol),PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
864: KSPGetPC(contour->ksp[i],&pc);
865: if (nsplit) {
866: KSPSetType(contour->ksp[i],KSPBCGS);
867: PCSetType(pc,PCBJACOBI);
868: } else {
869: KSPSetType(contour->ksp[i],KSPPREONLY);
870: PCSetType(pc,PCLU);
871: }
872: }
873: }
874: if (nsolve) *nsolve = contour->npoints;
875: if (ksp) *ksp = contour->ksp;
876: return 0;
877: }
879: /*@C
880: PEPCISSGetKSPs - Retrieve the array of linear solver objects associated with
881: the CISS solver.
883: Not Collective
885: Input Parameter:
886: . pep - polynomial eigenvalue solver
888: Output Parameters:
889: + nsolve - number of solver objects
890: - ksp - array of linear solver object
892: Notes:
893: The number of KSP solvers is equal to the number of integration points divided by
894: the number of partitions. This value is halved in the case of real matrices with
895: a region centered at the real axis.
897: Level: advanced
899: .seealso: PEPCISSSetSizes()
900: @*/
901: PetscErrorCode PEPCISSGetKSPs(PEP pep,PetscInt *nsolve,KSP **ksp)
902: {
904: PetscUseMethod(pep,"PEPCISSGetKSPs_C",(PEP,PetscInt*,KSP**),(pep,nsolve,ksp));
905: return 0;
906: }
908: PetscErrorCode PEPReset_CISS(PEP pep)
909: {
910: PEP_CISS *ctx = (PEP_CISS*)pep->data;
912: BVDestroy(&ctx->S);
913: BVDestroy(&ctx->V);
914: BVDestroy(&ctx->Y);
915: SlepcContourDataReset(ctx->contour);
916: MatDestroy(&ctx->J);
917: BVDestroy(&ctx->pV);
918: PetscFree(ctx->Psplit);
919: return 0;
920: }
922: PetscErrorCode PEPSetFromOptions_CISS(PEP pep,PetscOptionItems *PetscOptionsObject)
923: {
924: PEP_CISS *ctx = (PEP_CISS*)pep->data;
925: PetscReal r1,r2;
926: PetscInt i,i1,i2,i3,i4,i5,i6,i7;
927: PetscBool b1,flg,flg2,flg3,flg4,flg5,flg6;
928: PEPCISSExtraction extraction;
930: PetscOptionsHeadBegin(PetscOptionsObject,"PEP CISS Options");
932: PEPCISSGetSizes(pep,&i1,&i2,&i3,&i4,&i5,&b1);
933: PetscOptionsInt("-pep_ciss_integration_points","Number of integration points","PEPCISSSetSizes",i1,&i1,&flg);
934: PetscOptionsInt("-pep_ciss_blocksize","Block size","PEPCISSSetSizes",i2,&i2,&flg2);
935: PetscOptionsInt("-pep_ciss_moments","Moment size","PEPCISSSetSizes",i3,&i3,&flg3);
936: PetscOptionsInt("-pep_ciss_partitions","Number of partitions","PEPCISSSetSizes",i4,&i4,&flg4);
937: PetscOptionsInt("-pep_ciss_maxblocksize","Maximum block size","PEPCISSSetSizes",i5,&i5,&flg5);
938: PetscOptionsBool("-pep_ciss_realmats","True if all coefficient matrices of P(.) are real","PEPCISSSetSizes",b1,&b1,&flg6);
939: if (flg || flg2 || flg3 || flg4 || flg5 || flg6) PEPCISSSetSizes(pep,i1,i2,i3,i4,i5,b1);
941: PEPCISSGetThreshold(pep,&r1,&r2);
942: PetscOptionsReal("-pep_ciss_delta","Threshold for numerical rank","PEPCISSSetThreshold",r1,&r1,&flg);
943: PetscOptionsReal("-pep_ciss_spurious_threshold","Threshold for the spurious eigenpairs","PEPCISSSetThreshold",r2,&r2,&flg2);
944: if (flg || flg2) PEPCISSSetThreshold(pep,r1,r2);
946: PEPCISSGetRefinement(pep,&i6,&i7);
947: PetscOptionsInt("-pep_ciss_refine_inner","Number of inner iterative refinement iterations","PEPCISSSetRefinement",i6,&i6,&flg);
948: PetscOptionsInt("-pep_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","PEPCISSSetRefinement",i7,&i7,&flg2);
949: if (flg || flg2) PEPCISSSetRefinement(pep,i6,i7);
951: PetscOptionsEnum("-pep_ciss_extraction","Extraction technique","PEPCISSSetExtraction",PEPCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg);
952: if (flg) PEPCISSSetExtraction(pep,extraction);
954: PetscOptionsHeadEnd();
956: if (!pep->rg) PEPGetRG(pep,&pep->rg);
957: RGSetFromOptions(pep->rg); /* this is necessary here to set useconj */
958: if (!ctx->contour || !ctx->contour->ksp) PEPCISSGetKSPs(pep,NULL,NULL);
959: for (i=0;i<ctx->contour->npoints;i++) KSPSetFromOptions(ctx->contour->ksp[i]);
960: PetscSubcommSetFromOptions(ctx->contour->subcomm);
961: return 0;
962: }
964: PetscErrorCode PEPDestroy_CISS(PEP pep)
965: {
966: PEP_CISS *ctx = (PEP_CISS*)pep->data;
968: SlepcContourDataDestroy(&ctx->contour);
969: PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
970: PetscFree(pep->data);
971: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetSizes_C",NULL);
972: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetSizes_C",NULL);
973: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetThreshold_C",NULL);
974: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetThreshold_C",NULL);
975: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetRefinement_C",NULL);
976: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetRefinement_C",NULL);
977: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetExtraction_C",NULL);
978: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetExtraction_C",NULL);
979: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetKSPs_C",NULL);
980: return 0;
981: }
983: PetscErrorCode PEPView_CISS(PEP pep,PetscViewer viewer)
984: {
985: PEP_CISS *ctx = (PEP_CISS*)pep->data;
986: PetscBool isascii;
987: PetscViewer sviewer;
989: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
990: if (isascii) {
991: 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);
992: if (ctx->isreal) PetscViewerASCIIPrintf(viewer," exploiting symmetry of integration points\n");
993: PetscViewerASCIIPrintf(viewer," threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold);
994: PetscViewerASCIIPrintf(viewer," iterative refinement { inner: %" PetscInt_FMT ", blocksize: %" PetscInt_FMT " }\n",ctx->refine_inner, ctx->refine_blocksize);
995: PetscViewerASCIIPrintf(viewer," extraction: %s\n",PEPCISSExtractions[ctx->extraction]);
996: if (!ctx->contour || !ctx->contour->ksp) PEPCISSGetKSPs(pep,NULL,NULL);
997: PetscViewerASCIIPushTab(viewer);
998: if (ctx->npart>1 && ctx->contour->subcomm) {
999: PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1000: if (!ctx->contour->subcomm->color) KSPView(ctx->contour->ksp[0],sviewer);
1001: PetscViewerFlush(sviewer);
1002: PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1003: PetscViewerFlush(viewer);
1004: /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1005: PetscViewerASCIIPopSynchronized(viewer);
1006: } else KSPView(ctx->contour->ksp[0],viewer);
1007: PetscViewerASCIIPopTab(viewer);
1008: }
1009: return 0;
1010: }
1012: SLEPC_EXTERN PetscErrorCode PEPCreate_CISS(PEP pep)
1013: {
1014: PEP_CISS *ctx = (PEP_CISS*)pep->data;
1016: PetscNew(&ctx);
1017: pep->data = ctx;
1018: /* set default values of parameters */
1019: ctx->N = 32;
1020: ctx->L = 16;
1021: ctx->M = ctx->N/4;
1022: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
1023: ctx->L_max = 64;
1024: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1025: ctx->isreal = PETSC_FALSE;
1026: ctx->npart = 1;
1028: pep->ops->solve = PEPSolve_CISS;
1029: pep->ops->setup = PEPSetUp_CISS;
1030: pep->ops->setfromoptions = PEPSetFromOptions_CISS;
1031: pep->ops->reset = PEPReset_CISS;
1032: pep->ops->destroy = PEPDestroy_CISS;
1033: pep->ops->view = PEPView_CISS;
1035: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetSizes_C",PEPCISSSetSizes_CISS);
1036: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetSizes_C",PEPCISSGetSizes_CISS);
1037: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetThreshold_C",PEPCISSSetThreshold_CISS);
1038: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetThreshold_C",PEPCISSGetThreshold_CISS);
1039: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetRefinement_C",PEPCISSSetRefinement_CISS);
1040: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetRefinement_C",PEPCISSGetRefinement_CISS);
1041: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetExtraction_C",PEPCISSSetExtraction_CISS);
1042: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetExtraction_C",PEPCISSGetExtraction_CISS);
1043: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetKSPs_C",PEPCISSGetKSPs_CISS);
1044: return 0;
1045: }