Actual source code: pciss.c
slepc-3.18.2 2023-01-26
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: PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Something went wrong with PEPCISSGetKSPs()");
104: for (i=0;i<contour->npoints;i++) {
105: p_id = i*contour->subcomm->n + contour->subcomm->color;
106: MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&Amat);
107: if (T != P) MatDuplicate(P,MAT_DO_NOT_COPY_VALUES,&Pmat); else Pmat = Amat;
108: PEPComputeFunction(pep,ctx->omega[p_id],Amat,Pmat,PETSC_FALSE);
109: PEP_KSPSetOperators(contour->ksp[i],Amat,Pmat);
110: MatDestroy(&Amat);
111: if (T != P) MatDestroy(&Pmat);
112: }
113: return 0;
114: }
116: /*
117: Y_i = F(z_i)^{-1}Fp(z_i)V for every integration point, Y=[Y_i] is in the context
118: */
119: static PetscErrorCode PEPCISSSolve(PEP pep,Mat dT,BV V,PetscInt L_start,PetscInt L_end)
120: {
121: PEP_CISS *ctx = (PEP_CISS*)pep->data;
122: SlepcContourData contour;
123: PetscInt i,p_id;
124: Mat MV,BMV=NULL,MC;
126: contour = ctx->contour;
127: BVSetActiveColumns(V,L_start,L_end);
128: BVGetMat(V,&MV);
129: for (i=0;i<contour->npoints;i++) {
130: p_id = i*contour->subcomm->n + contour->subcomm->color;
131: PEPComputeFunction(pep,ctx->omega[p_id],dT,NULL,PETSC_TRUE);
132: BVSetActiveColumns(ctx->Y,i*ctx->L+L_start,i*ctx->L+L_end);
133: BVGetMat(ctx->Y,&MC);
134: if (!i) {
135: MatProductCreate(dT,MV,NULL,&BMV);
136: MatProductSetType(BMV,MATPRODUCT_AB);
137: MatProductSetFromOptions(BMV);
138: MatProductSymbolic(BMV);
139: }
140: MatProductNumeric(BMV);
141: KSPMatSolve(contour->ksp[i],BMV,MC);
142: BVRestoreMat(ctx->Y,&MC);
143: }
144: MatDestroy(&BMV);
145: BVRestoreMat(V,&MV);
146: return 0;
147: }
149: PetscErrorCode PEPSetUp_CISS(PEP pep)
150: {
151: PEP_CISS *ctx = (PEP_CISS*)pep->data;
152: SlepcContourData contour;
153: PetscInt i,nwork,nsplit;
154: PetscBool istrivial,isellipse,flg;
155: PetscObjectId id;
156: PetscObjectState state;
157: Vec v0;
159: if (pep->ncv==PETSC_DEFAULT) pep->ncv = ctx->L_max*ctx->M;
160: else {
161: ctx->L_max = pep->ncv/ctx->M;
162: if (!ctx->L_max) {
163: ctx->L_max = 1;
164: pep->ncv = ctx->L_max*ctx->M;
165: }
166: }
167: ctx->L = PetscMin(ctx->L,ctx->L_max);
168: if (pep->max_it==PETSC_DEFAULT) pep->max_it = 5;
169: if (pep->mpd==PETSC_DEFAULT) pep->mpd = pep->ncv;
170: if (!pep->which) pep->which = PEP_ALL;
172: PEPCheckUnsupported(pep,PEP_FEATURE_STOPPING);
173: PEPCheckIgnored(pep,PEP_FEATURE_SCALE);
175: /* check region */
176: RGIsTrivial(pep->rg,&istrivial);
178: RGGetComplement(pep->rg,&flg);
180: PetscObjectTypeCompare((PetscObject)pep->rg,RGELLIPSE,&isellipse);
183: /* if the region has changed, then reset contour data */
184: PetscObjectGetId((PetscObject)pep->rg,&id);
185: PetscObjectStateGet((PetscObject)pep->rg,&state);
186: if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
187: SlepcContourDataDestroy(&ctx->contour);
188: PetscInfo(pep,"Resetting the contour data structure due to a change of region\n");
189: ctx->rgid = id; ctx->rgstate = state;
190: }
192: /* create contour data structure */
193: if (!ctx->contour) {
194: RGCanUseConjugates(pep->rg,ctx->isreal,&ctx->useconj);
195: SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)pep,&ctx->contour);
196: }
198: PEPAllocateSolution(pep,0);
199: if (ctx->weight) PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
200: PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma);
202: /* allocate basis vectors */
203: BVDestroy(&ctx->S);
204: BVDuplicateResize(pep->V,ctx->L*ctx->M,&ctx->S);
205: BVDestroy(&ctx->V);
206: BVDuplicateResize(pep->V,ctx->L,&ctx->V);
208: /* check if a user-defined split preconditioner has been set */
209: STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL);
210: if (nsplit) {
211: PetscFree(ctx->Psplit);
212: PetscMalloc1(nsplit,&ctx->Psplit);
213: for (i=0;i<nsplit;i++) STGetSplitPreconditionerTerm(pep->st,i,&ctx->Psplit[i]);
214: }
216: contour = ctx->contour;
217: SlepcContourRedundantMat(contour,pep->nmat,pep->A,ctx->Psplit);
218: if (!ctx->J) MatDuplicate(contour->pA?contour->pA[0]:pep->A[0],MAT_DO_NOT_COPY_VALUES,&ctx->J);
219: if (contour->pA) {
220: BVGetColumn(ctx->V,0,&v0);
221: SlepcContourScatterCreate(contour,v0);
222: BVRestoreColumn(ctx->V,0,&v0);
223: BVDestroy(&ctx->pV);
224: BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV);
225: BVSetSizesFromVec(ctx->pV,contour->xsub,pep->n);
226: BVSetFromOptions(ctx->pV);
227: BVResize(ctx->pV,ctx->L,PETSC_FALSE);
228: }
230: BVDestroy(&ctx->Y);
231: if (contour->pA) {
232: BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y);
233: BVSetSizesFromVec(ctx->Y,contour->xsub,pep->n);
234: BVSetFromOptions(ctx->Y);
235: BVResize(ctx->Y,contour->npoints*ctx->L,PETSC_FALSE);
236: } else BVDuplicateResize(pep->V,contour->npoints*ctx->L,&ctx->Y);
238: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) DSSetType(pep->ds,DSGNHEP);
239: else if (ctx->extraction == PEP_CISS_EXTRACTION_CAA) DSSetType(pep->ds,DSNHEP);
240: else {
241: DSSetType(pep->ds,DSPEP);
242: DSPEPSetDegree(pep->ds,pep->nmat-1);
243: DSPEPSetCoefficients(pep->ds,pep->pbc);
244: }
245: DSAllocate(pep->ds,pep->ncv);
246: nwork = 2;
247: PEPSetWorkVecs(pep,nwork);
248: return 0;
249: }
251: PetscErrorCode PEPSolve_CISS(PEP pep)
252: {
253: PEP_CISS *ctx = (PEP_CISS*)pep->data;
254: SlepcContourData contour = ctx->contour;
255: Mat X,M,E,T,P;
256: PetscInt i,j,ld,L_add=0,nv=0,L_base=ctx->L,inner,*inside,nsplit;
257: PetscScalar *Mu,*H0,*H1,*rr,*temp,center;
258: PetscReal error,max_error,radius,rgscale,est_eig,eta;
259: PetscBool isellipse,*fl1;
260: Vec si;
261: SlepcSC sc;
262: PetscRandom rand;
264: DSSetFromOptions(pep->ds);
265: DSGetSlepcSC(pep->ds,&sc);
266: sc->comparison = SlepcCompareLargestMagnitude;
267: sc->comparisonctx = NULL;
268: sc->map = NULL;
269: sc->mapobj = NULL;
270: DSGetLeadingDimension(pep->ds,&ld);
271: RGComputeQuadrature(pep->rg,RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight);
272: STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL);
273: if (contour->pA) {
274: T = contour->pA[0];
275: P = nsplit? contour->pP[0]: T;
276: } else {
277: T = pep->A[0];
278: P = nsplit? ctx->Psplit[0]: T;
279: }
280: PEPCISSSetUp(pep,T,P);
281: BVSetActiveColumns(ctx->V,0,ctx->L);
282: BVSetRandomSign(ctx->V);
283: BVGetRandomContext(ctx->V,&rand);
284: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
285: PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L);
286: PetscObjectTypeCompare((PetscObject)pep->rg,RGELLIPSE,&isellipse);
287: if (isellipse) {
288: BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig);
289: PetscInfo(pep,"Estimated eigenvalue count: %f\n",(double)est_eig);
290: eta = PetscPowReal(10.0,-PetscLog10Real(pep->tol)/ctx->N);
291: L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
292: if (L_add>ctx->L_max-ctx->L) {
293: PetscInfo(pep,"Number of eigenvalues inside the contour path may be too large\n");
294: L_add = ctx->L_max-ctx->L;
295: }
296: }
297: /* Updates L after estimate the number of eigenvalue */
298: if (L_add>0) {
299: PetscInfo(pep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by Estimate #Eig\n",ctx->L,ctx->L+L_add);
300: BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints);
301: BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
302: BVSetRandomSign(ctx->V);
303: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
304: ctx->L += L_add;
305: PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,ctx->L-L_add,ctx->L);
306: }
308: PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
309: for (i=0;i<ctx->refine_blocksize;i++) {
310: 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);
311: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
312: PetscLogEventBegin(PEP_CISS_SVD,pep,0,0,0);
313: SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
314: PetscLogEventEnd(PEP_CISS_SVD,pep,0,0,0);
315: if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
316: L_add = L_base;
317: if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
318: PetscInfo(pep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by SVD(H0)\n",ctx->L,ctx->L+L_add);
319: BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints);
320: BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
321: BVSetRandomSign(ctx->V);
322: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
323: ctx->L += L_add;
324: PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,ctx->L-L_add,ctx->L);
325: if (L_add) {
326: PetscFree2(Mu,H0);
327: PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
328: }
329: }
331: RGGetScale(pep->rg,&rgscale);
332: RGEllipseGetParameters(pep->rg,¢er,&radius,NULL);
334: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1);
336: while (pep->reason == PEP_CONVERGED_ITERATING) {
337: pep->its++;
338: for (inner=0;inner<=ctx->refine_inner;inner++) {
339: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
340: 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);
341: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
342: PetscLogEventBegin(PEP_CISS_SVD,pep,0,0,0);
343: SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
344: PetscLogEventEnd(PEP_CISS_SVD,pep,0,0,0);
345: } else {
346: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
347: /* compute SVD of S */
348: 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);
349: }
350: PetscInfo(pep,"Estimated rank: nv = %" PetscInt_FMT "\n",nv);
351: if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
352: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
353: BVSetActiveColumns(ctx->S,0,ctx->L);
354: BVSetActiveColumns(ctx->V,0,ctx->L);
355: BVCopy(ctx->S,ctx->V);
356: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
357: PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L);
358: } else break;
359: }
360: pep->nconv = 0;
361: if (nv == 0) { pep->reason = PEP_CONVERGED_TOL; break; }
362: else {
363: /* Extracting eigenpairs */
364: DSSetDimensions(pep->ds,nv,0,0);
365: DSSetState(pep->ds,DS_STATE_RAW);
366: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
367: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
368: CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1);
369: DSGetArray(pep->ds,DS_MAT_A,&temp);
370: for (j=0;j<nv;j++)
371: for (i=0;i<nv;i++)
372: temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
373: DSRestoreArray(pep->ds,DS_MAT_A,&temp);
374: DSGetArray(pep->ds,DS_MAT_B,&temp);
375: for (j=0;j<nv;j++)
376: for (i=0;i<nv;i++)
377: temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
378: DSRestoreArray(pep->ds,DS_MAT_B,&temp);
379: } else if (ctx->extraction == PEP_CISS_EXTRACTION_CAA) {
380: BVSetActiveColumns(ctx->S,0,nv);
381: DSGetArray(pep->ds,DS_MAT_A,&temp);
382: for (i=0;i<nv;i++) PetscArraycpy(temp+i*ld,H0+i*nv,nv);
383: DSRestoreArray(pep->ds,DS_MAT_A,&temp);
384: } else {
385: BVSetActiveColumns(ctx->S,0,nv);
386: for (i=0;i<pep->nmat;i++) {
387: DSGetMat(pep->ds,DSMatExtra[i],&E);
388: BVMatProject(ctx->S,pep->A[i],ctx->S,E);
389: DSRestoreMat(pep->ds,DSMatExtra[i],&E);
390: }
391: nv = (pep->nmat-1)*nv;
392: }
393: DSSolve(pep->ds,pep->eigr,pep->eigi);
394: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
395: if (ctx->extraction == PEP_CISS_EXTRACTION_CAA || ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
396: for (i=0;i<nv;i++) {
397: pep->eigr[i] = (pep->eigr[i]*radius+center)*rgscale;
398: }
399: }
400: PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr);
401: DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
402: DSGetMat(pep->ds,DS_MAT_X,&X);
403: SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1);
404: DSRestoreMat(pep->ds,DS_MAT_X,&X);
405: RGCheckInside(pep->rg,nv,pep->eigr,pep->eigi,inside);
406: for (i=0;i<nv;i++) {
407: if (fl1[i] && inside[i]>=0) {
408: rr[i] = 1.0;
409: pep->nconv++;
410: } else rr[i] = 0.0;
411: }
412: DSSort(pep->ds,pep->eigr,pep->eigi,rr,NULL,&pep->nconv);
413: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
414: if (ctx->extraction == PEP_CISS_EXTRACTION_CAA || ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
415: for (i=0;i<nv;i++) pep->eigr[i] = (pep->eigr[i]*radius+center)*rgscale;
416: }
417: PetscFree3(fl1,inside,rr);
418: BVSetActiveColumns(pep->V,0,nv);
419: DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
420: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
421: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
422: BVSetActiveColumns(ctx->S,0,nv);
423: BVCopy(ctx->S,pep->V);
424: DSGetMat(pep->ds,DS_MAT_X,&X);
425: BVMultInPlace(ctx->S,X,0,pep->nconv);
426: BVMultInPlace(pep->V,X,0,pep->nconv);
427: DSRestoreMat(pep->ds,DS_MAT_X,&X);
428: } else {
429: DSGetMat(pep->ds,DS_MAT_X,&X);
430: BVMultInPlace(ctx->S,X,0,pep->nconv);
431: DSRestoreMat(pep->ds,DS_MAT_X,&X);
432: BVSetActiveColumns(ctx->S,0,pep->nconv);
433: BVCopy(ctx->S,pep->V);
434: }
435: max_error = 0.0;
436: for (i=0;i<pep->nconv;i++) {
437: BVGetColumn(pep->V,i,&si);
438: VecNormalize(si,NULL);
439: PEPComputeResidualNorm_Private(pep,pep->eigr[i],0,si,NULL,pep->work,&error);
440: (*pep->converged)(pep,pep->eigr[i],0,error,&error,pep->convergedctx);
441: BVRestoreColumn(pep->V,i,&si);
442: max_error = PetscMax(max_error,error);
443: }
444: if (max_error <= pep->tol) pep->reason = PEP_CONVERGED_TOL;
445: else if (pep->its > pep->max_it) pep->reason = PEP_DIVERGED_ITS;
446: else {
447: if (pep->nconv > ctx->L) nv = pep->nconv;
448: else if (ctx->L > nv) nv = ctx->L;
449: nv = PetscMin(nv,ctx->L*ctx->M);
450: MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M);
451: MatSetRandom(M,rand);
452: BVSetActiveColumns(ctx->S,0,nv);
453: BVMultInPlace(ctx->S,M,0,ctx->L);
454: MatDestroy(&M);
455: BVSetActiveColumns(ctx->S,0,ctx->L);
456: BVSetActiveColumns(ctx->V,0,ctx->L);
457: BVCopy(ctx->S,ctx->V);
458: if (contour->pA) BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
459: PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L);
460: }
461: }
462: }
463: PetscFree2(Mu,H0);
464: if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) PetscFree(H1);
465: return 0;
466: }
468: static PetscErrorCode PEPCISSSetSizes_CISS(PEP pep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
469: {
470: PEP_CISS *ctx = (PEP_CISS*)pep->data;
471: PetscInt oN,oL,oM,oLmax,onpart;
472: PetscMPIInt size;
474: oN = ctx->N;
475: if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
476: if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
477: } else {
480: if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
481: }
482: oL = ctx->L;
483: if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
484: ctx->L = 16;
485: } else {
487: ctx->L = bs;
488: }
489: oM = ctx->M;
490: if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
491: ctx->M = ctx->N/4;
492: } else {
495: ctx->M = PetscMax(ms,2);
496: }
497: onpart = ctx->npart;
498: if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
499: ctx->npart = 1;
500: } else {
501: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size);
503: ctx->npart = npart;
504: }
505: oLmax = ctx->L_max;
506: if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
507: ctx->L_max = 64;
508: } else {
510: ctx->L_max = PetscMax(bsmax,ctx->L);
511: }
512: if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
513: SlepcContourDataDestroy(&ctx->contour);
514: PetscInfo(pep,"Resetting the contour data structure due to a change of parameters\n");
515: pep->state = PEP_STATE_INITIAL;
516: }
517: ctx->isreal = realmats;
518: if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) pep->state = PEP_STATE_INITIAL;
519: return 0;
520: }
522: /*@
523: PEPCISSSetSizes - Sets the values of various size parameters in the CISS solver.
525: Logically Collective on pep
527: Input Parameters:
528: + pep - the polynomial eigensolver context
529: . ip - number of integration points
530: . bs - block size
531: . ms - moment size
532: . npart - number of partitions when splitting the communicator
533: . bsmax - max block size
534: - realmats - all coefficient matrices of P(.) are real
536: Options Database Keys:
537: + -pep_ciss_integration_points - Sets the number of integration points
538: . -pep_ciss_blocksize - Sets the block size
539: . -pep_ciss_moments - Sets the moment size
540: . -pep_ciss_partitions - Sets the number of partitions
541: . -pep_ciss_maxblocksize - Sets the maximum block size
542: - -pep_ciss_realmats - all coefficient matrices of P(.) are real
544: Notes:
545: The default number of partitions is 1. This means the internal KSP object is shared
546: among all processes of the PEP communicator. Otherwise, the communicator is split
547: into npart communicators, so that npart KSP solves proceed simultaneously.
549: Level: advanced
551: .seealso: PEPCISSGetSizes()
552: @*/
553: PetscErrorCode PEPCISSSetSizes(PEP pep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
554: {
562: PetscTryMethod(pep,"PEPCISSSetSizes_C",(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(pep,ip,bs,ms,npart,bsmax,realmats));
563: return 0;
564: }
566: static PetscErrorCode PEPCISSGetSizes_CISS(PEP pep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
567: {
568: PEP_CISS *ctx = (PEP_CISS*)pep->data;
570: if (ip) *ip = ctx->N;
571: if (bs) *bs = ctx->L;
572: if (ms) *ms = ctx->M;
573: if (npart) *npart = ctx->npart;
574: if (bsmax) *bsmax = ctx->L_max;
575: if (realmats) *realmats = ctx->isreal;
576: return 0;
577: }
579: /*@
580: PEPCISSGetSizes - Gets the values of various size parameters in the CISS solver.
582: Not Collective
584: Input Parameter:
585: . pep - the polynomial eigensolver context
587: Output Parameters:
588: + ip - number of integration points
589: . bs - block size
590: . ms - moment size
591: . npart - number of partitions when splitting the communicator
592: . bsmax - max block size
593: - realmats - all coefficient matrices of P(.) are real
595: Level: advanced
597: .seealso: PEPCISSSetSizes()
598: @*/
599: PetscErrorCode PEPCISSGetSizes(PEP pep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
600: {
602: PetscUseMethod(pep,"PEPCISSGetSizes_C",(PEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(pep,ip,bs,ms,npart,bsmax,realmats));
603: return 0;
604: }
606: static PetscErrorCode PEPCISSSetThreshold_CISS(PEP pep,PetscReal delta,PetscReal spur)
607: {
608: PEP_CISS *ctx = (PEP_CISS*)pep->data;
610: if (delta == PETSC_DEFAULT) {
611: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
612: } else {
614: ctx->delta = delta;
615: }
616: if (spur == PETSC_DEFAULT) {
617: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
618: } else {
620: ctx->spurious_threshold = spur;
621: }
622: return 0;
623: }
625: /*@
626: PEPCISSSetThreshold - Sets the values of various threshold parameters in
627: the CISS solver.
629: Logically Collective on pep
631: Input Parameters:
632: + pep - the polynomial eigensolver context
633: . delta - threshold for numerical rank
634: - spur - spurious threshold (to discard spurious eigenpairs)
636: Options Database Keys:
637: + -pep_ciss_delta - Sets the delta
638: - -pep_ciss_spurious_threshold - Sets the spurious threshold
640: Level: advanced
642: .seealso: PEPCISSGetThreshold()
643: @*/
644: PetscErrorCode PEPCISSSetThreshold(PEP pep,PetscReal delta,PetscReal spur)
645: {
649: PetscTryMethod(pep,"PEPCISSSetThreshold_C",(PEP,PetscReal,PetscReal),(pep,delta,spur));
650: return 0;
651: }
653: static PetscErrorCode PEPCISSGetThreshold_CISS(PEP pep,PetscReal *delta,PetscReal *spur)
654: {
655: PEP_CISS *ctx = (PEP_CISS*)pep->data;
657: if (delta) *delta = ctx->delta;
658: if (spur) *spur = ctx->spurious_threshold;
659: return 0;
660: }
662: /*@
663: PEPCISSGetThreshold - Gets the values of various threshold parameters in
664: the CISS solver.
666: Not Collective
668: Input Parameter:
669: . pep - the polynomial eigensolver context
671: Output Parameters:
672: + delta - threshold for numerical rank
673: - spur - spurious threshold (to discard spurious eigenpairs)
675: Level: advanced
677: .seealso: PEPCISSSetThreshold()
678: @*/
679: PetscErrorCode PEPCISSGetThreshold(PEP pep,PetscReal *delta,PetscReal *spur)
680: {
682: PetscUseMethod(pep,"PEPCISSGetThreshold_C",(PEP,PetscReal*,PetscReal*),(pep,delta,spur));
683: return 0;
684: }
686: static PetscErrorCode PEPCISSSetRefinement_CISS(PEP pep,PetscInt inner,PetscInt blsize)
687: {
688: PEP_CISS *ctx = (PEP_CISS*)pep->data;
690: if (inner == PETSC_DEFAULT) {
691: ctx->refine_inner = 0;
692: } else {
694: ctx->refine_inner = inner;
695: }
696: if (blsize == PETSC_DEFAULT) {
697: ctx->refine_blocksize = 0;
698: } else {
700: ctx->refine_blocksize = blsize;
701: }
702: return 0;
703: }
705: /*@
706: PEPCISSSetRefinement - Sets the values of various refinement parameters
707: in the CISS solver.
709: Logically Collective on pep
711: Input Parameters:
712: + pep - the polynomial eigensolver context
713: . inner - number of iterative refinement iterations (inner loop)
714: - blsize - number of iterative refinement iterations (blocksize loop)
716: Options Database Keys:
717: + -pep_ciss_refine_inner - Sets number of inner iterations
718: - -pep_ciss_refine_blocksize - Sets number of blocksize iterations
720: Level: advanced
722: .seealso: PEPCISSGetRefinement()
723: @*/
724: PetscErrorCode PEPCISSSetRefinement(PEP pep,PetscInt inner,PetscInt blsize)
725: {
729: PetscTryMethod(pep,"PEPCISSSetRefinement_C",(PEP,PetscInt,PetscInt),(pep,inner,blsize));
730: return 0;
731: }
733: static PetscErrorCode PEPCISSGetRefinement_CISS(PEP pep,PetscInt *inner,PetscInt *blsize)
734: {
735: PEP_CISS *ctx = (PEP_CISS*)pep->data;
737: if (inner) *inner = ctx->refine_inner;
738: if (blsize) *blsize = ctx->refine_blocksize;
739: return 0;
740: }
742: /*@
743: PEPCISSGetRefinement - Gets the values of various refinement parameters
744: in the CISS solver.
746: Not Collective
748: Input Parameter:
749: . pep - the polynomial eigensolver context
751: Output Parameters:
752: + inner - number of iterative refinement iterations (inner loop)
753: - blsize - number of iterative refinement iterations (blocksize loop)
755: Level: advanced
757: .seealso: PEPCISSSetRefinement()
758: @*/
759: PetscErrorCode PEPCISSGetRefinement(PEP pep, PetscInt *inner, PetscInt *blsize)
760: {
762: PetscUseMethod(pep,"PEPCISSGetRefinement_C",(PEP,PetscInt*,PetscInt*),(pep,inner,blsize));
763: return 0;
764: }
766: static PetscErrorCode PEPCISSSetExtraction_CISS(PEP pep,PEPCISSExtraction extraction)
767: {
768: PEP_CISS *ctx = (PEP_CISS*)pep->data;
770: if (ctx->extraction != extraction) {
771: ctx->extraction = extraction;
772: pep->state = PEP_STATE_INITIAL;
773: }
774: return 0;
775: }
777: /*@
778: PEPCISSSetExtraction - Sets the extraction technique used in the CISS solver.
780: Logically Collective on pep
782: Input Parameters:
783: + pep - the polynomial eigensolver context
784: - extraction - the extraction technique
786: Options Database Key:
787: . -pep_ciss_extraction - Sets the extraction technique (either 'ritz', 'hankel' or 'caa')
789: Notes:
790: By default, the Rayleigh-Ritz extraction is used (PEP_CISS_EXTRACTION_RITZ).
792: If the 'hankel' or the 'caa' option is specified (PEP_CISS_EXTRACTION_HANKEL or
793: PEP_CISS_EXTRACTION_CAA), then the Block Hankel method, or the Communication-avoiding
794: Arnoldi method, respectively, is used for extracting eigenpairs.
796: Level: advanced
798: .seealso: PEPCISSGetExtraction(), PEPCISSExtraction
799: @*/
800: PetscErrorCode PEPCISSSetExtraction(PEP pep,PEPCISSExtraction extraction)
801: {
804: PetscTryMethod(pep,"PEPCISSSetExtraction_C",(PEP,PEPCISSExtraction),(pep,extraction));
805: return 0;
806: }
808: static PetscErrorCode PEPCISSGetExtraction_CISS(PEP pep,PEPCISSExtraction *extraction)
809: {
810: PEP_CISS *ctx = (PEP_CISS*)pep->data;
812: *extraction = ctx->extraction;
813: return 0;
814: }
816: /*@
817: PEPCISSGetExtraction - Gets the extraction technique used in the CISS solver.
819: Not Collective
821: Input Parameter:
822: . pep - the polynomial eigensolver context
824: Output Parameters:
825: . extraction - extraction technique
827: Level: advanced
829: .seealso: PEPCISSSetExtraction() PEPCISSExtraction
830: @*/
831: PetscErrorCode PEPCISSGetExtraction(PEP pep,PEPCISSExtraction *extraction)
832: {
835: PetscUseMethod(pep,"PEPCISSGetExtraction_C",(PEP,PEPCISSExtraction*),(pep,extraction));
836: return 0;
837: }
839: static PetscErrorCode PEPCISSGetKSPs_CISS(PEP pep,PetscInt *nsolve,KSP **ksp)
840: {
841: PEP_CISS *ctx = (PEP_CISS*)pep->data;
842: SlepcContourData contour;
843: PetscInt i,nsplit;
844: PC pc;
845: MPI_Comm child;
847: if (!ctx->contour) { /* initialize contour data structure first */
848: RGCanUseConjugates(pep->rg,ctx->isreal,&ctx->useconj);
849: SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)pep,&ctx->contour);
850: }
851: contour = ctx->contour;
852: if (!contour->ksp) {
853: PetscMalloc1(contour->npoints,&contour->ksp);
854: PEPGetST(pep,&pep->st);
855: STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL);
856: PetscSubcommGetChild(contour->subcomm,&child);
857: for (i=0;i<contour->npoints;i++) {
858: KSPCreate(child,&contour->ksp[i]);
859: PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)pep,1);
860: KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)pep)->prefix);
861: KSPAppendOptionsPrefix(contour->ksp[i],"pep_ciss_");
862: PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)pep)->options);
863: KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE);
864: KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(pep->tol),PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
865: KSPGetPC(contour->ksp[i],&pc);
866: if (nsplit) {
867: KSPSetType(contour->ksp[i],KSPBCGS);
868: PCSetType(pc,PCBJACOBI);
869: } else {
870: KSPSetType(contour->ksp[i],KSPPREONLY);
871: PCSetType(pc,PCLU);
872: }
873: }
874: }
875: if (nsolve) *nsolve = contour->npoints;
876: if (ksp) *ksp = contour->ksp;
877: return 0;
878: }
880: /*@C
881: PEPCISSGetKSPs - Retrieve the array of linear solver objects associated with
882: the CISS solver.
884: Not Collective
886: Input Parameter:
887: . pep - polynomial eigenvalue solver
889: Output Parameters:
890: + nsolve - number of solver objects
891: - ksp - array of linear solver object
893: Notes:
894: The number of KSP solvers is equal to the number of integration points divided by
895: the number of partitions. This value is halved in the case of real matrices with
896: a region centered at the real axis.
898: Level: advanced
900: .seealso: PEPCISSSetSizes()
901: @*/
902: PetscErrorCode PEPCISSGetKSPs(PEP pep,PetscInt *nsolve,KSP **ksp)
903: {
905: PetscUseMethod(pep,"PEPCISSGetKSPs_C",(PEP,PetscInt*,KSP**),(pep,nsolve,ksp));
906: return 0;
907: }
909: PetscErrorCode PEPReset_CISS(PEP pep)
910: {
911: PEP_CISS *ctx = (PEP_CISS*)pep->data;
913: BVDestroy(&ctx->S);
914: BVDestroy(&ctx->V);
915: BVDestroy(&ctx->Y);
916: SlepcContourDataReset(ctx->contour);
917: MatDestroy(&ctx->J);
918: BVDestroy(&ctx->pV);
919: PetscFree(ctx->Psplit);
920: return 0;
921: }
923: PetscErrorCode PEPSetFromOptions_CISS(PEP pep,PetscOptionItems *PetscOptionsObject)
924: {
925: PEP_CISS *ctx = (PEP_CISS*)pep->data;
926: PetscReal r1,r2;
927: PetscInt i,i1,i2,i3,i4,i5,i6,i7;
928: PetscBool b1,flg,flg2,flg3,flg4,flg5,flg6;
929: PEPCISSExtraction extraction;
931: PetscOptionsHeadBegin(PetscOptionsObject,"PEP CISS Options");
933: PEPCISSGetSizes(pep,&i1,&i2,&i3,&i4,&i5,&b1);
934: PetscOptionsInt("-pep_ciss_integration_points","Number of integration points","PEPCISSSetSizes",i1,&i1,&flg);
935: PetscOptionsInt("-pep_ciss_blocksize","Block size","PEPCISSSetSizes",i2,&i2,&flg2);
936: PetscOptionsInt("-pep_ciss_moments","Moment size","PEPCISSSetSizes",i3,&i3,&flg3);
937: PetscOptionsInt("-pep_ciss_partitions","Number of partitions","PEPCISSSetSizes",i4,&i4,&flg4);
938: PetscOptionsInt("-pep_ciss_maxblocksize","Maximum block size","PEPCISSSetSizes",i5,&i5,&flg5);
939: PetscOptionsBool("-pep_ciss_realmats","True if all coefficient matrices of P(.) are real","PEPCISSSetSizes",b1,&b1,&flg6);
940: if (flg || flg2 || flg3 || flg4 || flg5 || flg6) PEPCISSSetSizes(pep,i1,i2,i3,i4,i5,b1);
942: PEPCISSGetThreshold(pep,&r1,&r2);
943: PetscOptionsReal("-pep_ciss_delta","Threshold for numerical rank","PEPCISSSetThreshold",r1,&r1,&flg);
944: PetscOptionsReal("-pep_ciss_spurious_threshold","Threshold for the spurious eigenpairs","PEPCISSSetThreshold",r2,&r2,&flg2);
945: if (flg || flg2) PEPCISSSetThreshold(pep,r1,r2);
947: PEPCISSGetRefinement(pep,&i6,&i7);
948: PetscOptionsInt("-pep_ciss_refine_inner","Number of inner iterative refinement iterations","PEPCISSSetRefinement",i6,&i6,&flg);
949: PetscOptionsInt("-pep_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","PEPCISSSetRefinement",i7,&i7,&flg2);
950: if (flg || flg2) PEPCISSSetRefinement(pep,i6,i7);
952: PetscOptionsEnum("-pep_ciss_extraction","Extraction technique","PEPCISSSetExtraction",PEPCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg);
953: if (flg) PEPCISSSetExtraction(pep,extraction);
955: PetscOptionsHeadEnd();
957: if (!pep->rg) PEPGetRG(pep,&pep->rg);
958: RGSetFromOptions(pep->rg); /* this is necessary here to set useconj */
959: if (!ctx->contour || !ctx->contour->ksp) PEPCISSGetKSPs(pep,NULL,NULL);
960: PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Something went wrong with PEPCISSGetKSPs()");
961: for (i=0;i<ctx->contour->npoints;i++) KSPSetFromOptions(ctx->contour->ksp[i]);
962: PetscSubcommSetFromOptions(ctx->contour->subcomm);
963: return 0;
964: }
966: PetscErrorCode PEPDestroy_CISS(PEP pep)
967: {
968: PEP_CISS *ctx = (PEP_CISS*)pep->data;
970: SlepcContourDataDestroy(&ctx->contour);
971: PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
972: PetscFree(pep->data);
973: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetSizes_C",NULL);
974: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetSizes_C",NULL);
975: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetThreshold_C",NULL);
976: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetThreshold_C",NULL);
977: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetRefinement_C",NULL);
978: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetRefinement_C",NULL);
979: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetExtraction_C",NULL);
980: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetExtraction_C",NULL);
981: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetKSPs_C",NULL);
982: return 0;
983: }
985: PetscErrorCode PEPView_CISS(PEP pep,PetscViewer viewer)
986: {
987: PEP_CISS *ctx = (PEP_CISS*)pep->data;
988: PetscBool isascii;
989: PetscViewer sviewer;
991: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
992: if (isascii) {
993: 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);
994: if (ctx->isreal) PetscViewerASCIIPrintf(viewer," exploiting symmetry of integration points\n");
995: PetscViewerASCIIPrintf(viewer," threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold);
996: PetscViewerASCIIPrintf(viewer," iterative refinement { inner: %" PetscInt_FMT ", blocksize: %" PetscInt_FMT " }\n",ctx->refine_inner, ctx->refine_blocksize);
997: PetscViewerASCIIPrintf(viewer," extraction: %s\n",PEPCISSExtractions[ctx->extraction]);
998: if (!ctx->contour || !ctx->contour->ksp) PEPCISSGetKSPs(pep,NULL,NULL);
999: PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Something went wrong with PEPCISSGetKSPs()");
1000: PetscViewerASCIIPushTab(viewer);
1001: if (ctx->npart>1 && ctx->contour->subcomm) {
1002: PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1003: if (!ctx->contour->subcomm->color) KSPView(ctx->contour->ksp[0],sviewer);
1004: PetscViewerFlush(sviewer);
1005: PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1006: PetscViewerFlush(viewer);
1007: /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1008: PetscViewerASCIIPopSynchronized(viewer);
1009: } else KSPView(ctx->contour->ksp[0],viewer);
1010: PetscViewerASCIIPopTab(viewer);
1011: }
1012: return 0;
1013: }
1015: SLEPC_EXTERN PetscErrorCode PEPCreate_CISS(PEP pep)
1016: {
1017: PEP_CISS *ctx = (PEP_CISS*)pep->data;
1019: PetscNew(&ctx);
1020: pep->data = ctx;
1021: /* set default values of parameters */
1022: ctx->N = 32;
1023: ctx->L = 16;
1024: ctx->M = ctx->N/4;
1025: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
1026: ctx->L_max = 64;
1027: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1028: ctx->isreal = PETSC_FALSE;
1029: ctx->npart = 1;
1031: pep->ops->solve = PEPSolve_CISS;
1032: pep->ops->setup = PEPSetUp_CISS;
1033: pep->ops->setfromoptions = PEPSetFromOptions_CISS;
1034: pep->ops->reset = PEPReset_CISS;
1035: pep->ops->destroy = PEPDestroy_CISS;
1036: pep->ops->view = PEPView_CISS;
1038: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetSizes_C",PEPCISSSetSizes_CISS);
1039: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetSizes_C",PEPCISSGetSizes_CISS);
1040: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetThreshold_C",PEPCISSSetThreshold_CISS);
1041: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetThreshold_C",PEPCISSGetThreshold_CISS);
1042: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetRefinement_C",PEPCISSSetRefinement_CISS);
1043: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetRefinement_C",PEPCISSGetRefinement_CISS);
1044: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetExtraction_C",PEPCISSSetExtraction_CISS);
1045: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetExtraction_C",PEPCISSGetExtraction_CISS);
1046: PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetKSPs_C",PEPCISSGetKSPs_CISS);
1047: return 0;
1048: }