Actual source code: nciss.c
slepc-3.16.3 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: SLEPc eigensolver: "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 nonlinear eigenvalue problems using contour
27: integrals", JSIAM Lett. 1:52-55, 2009.
29: [2] S. Yokota and T. Sakurai, "A projection method for nonlinear
30: eigenvalue problems using contour integrals", JSIAM Lett.
31: 5:41-44, 2013.
32: */
34: #include <slepc/private/nepimpl.h>
35: #include <slepc/private/slepccontour.h>
37: typedef struct _n_nep_ciss_project *NEP_CISS_PROJECT;
38: typedef struct {
39: /* parameters */
40: PetscInt N; /* number of integration points (32) */
41: PetscInt L; /* block size (16) */
42: PetscInt M; /* moment degree (N/4 = 4) */
43: PetscReal delta; /* threshold of singular value (1e-12) */
44: PetscInt L_max; /* maximum number of columns of the source matrix V */
45: PetscReal spurious_threshold; /* discard spurious eigenpairs */
46: PetscBool isreal; /* T(z) is real for real z */
47: PetscInt npart; /* number of partitions */
48: PetscInt refine_inner;
49: PetscInt refine_blocksize;
50: NEPCISSExtraction extraction;
51: /* private data */
52: SlepcContourData contour;
53: PetscReal *sigma; /* threshold for numerical rank */
54: PetscScalar *weight;
55: PetscScalar *omega;
56: PetscScalar *pp;
57: BV V;
58: BV S;
59: BV Y;
60: PetscBool useconj;
61: Mat T,J; /* auxiliary matrices when using subcomm */
62: BV pV;
63: NEP_CISS_PROJECT dsctxf;
64: PetscObjectId rgid;
65: PetscObjectState rgstate;
66: } NEP_CISS;
68: struct _n_nep_ciss_project {
69: NEP nep;
70: BV Q;
71: };
73: static PetscErrorCode NEPContourDSComputeMatrix(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat,void *ctx)
74: {
75: NEP_CISS_PROJECT proj = (NEP_CISS_PROJECT)ctx;
76: PetscErrorCode ierr;
77: Mat M,fun;
80: if (!deriv) {
81: NEPComputeFunction(proj->nep,lambda,proj->nep->function,proj->nep->function);
82: fun = proj->nep->function;
83: } else {
84: NEPComputeJacobian(proj->nep,lambda,proj->nep->jacobian);
85: fun = proj->nep->jacobian;
86: }
87: DSGetMat(ds,mat,&M);
88: BVMatProject(proj->Q,fun,proj->Q,M);
89: DSRestoreMat(ds,mat,&M);
90: return(0);
91: }
93: static PetscErrorCode NEPComputeFunctionSubcomm(NEP nep,PetscScalar lambda,Mat T,PetscBool deriv)
94: {
96: PetscInt i;
97: PetscScalar alpha;
98: NEP_CISS *ctx = (NEP_CISS*)nep->data;
101: switch (nep->fui) {
102: case NEP_USER_INTERFACE_CALLBACK:
103: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Partitions are not supported with callbacks");
104: case NEP_USER_INTERFACE_SPLIT:
105: MatZeroEntries(T);
106: for (i=0;i<nep->nt;i++) {
107: if (!deriv) {
108: FNEvaluateFunction(nep->f[i],lambda,&alpha);
109: } else {
110: FNEvaluateDerivative(nep->f[i],lambda,&alpha);
111: }
112: MatAXPY(T,alpha,ctx->contour->pA[i],nep->mstr);
113: }
114: break;
115: }
116: return(0);
117: }
119: /*
120: Y_i = F(z_i)^{-1}Fp(z_i)V for every integration point, Y=[Y_i] is in the context
121: */
122: static PetscErrorCode NEPCISSSolveSystem(NEP nep,Mat T,Mat dT,BV V,PetscInt L_start,PetscInt L_end,PetscBool initksp)
123: {
124: PetscErrorCode ierr;
125: NEP_CISS *ctx = (NEP_CISS*)nep->data;
126: SlepcContourData contour = ctx->contour;
127: PetscInt i,p_id;
128: Mat kspMat,MV,BMV=NULL,MC;
131: if (!ctx->contour || !ctx->contour->ksp) { NEPCISSGetKSPs(nep,NULL,NULL); }
132: BVSetActiveColumns(V,L_start,L_end);
133: BVGetMat(V,&MV);
134: for (i=0;i<contour->npoints;i++) {
135: p_id = i*contour->subcomm->n + contour->subcomm->color;
136: if (initksp) {
137: if (contour->subcomm->n == 1) {
138: NEPComputeFunction(nep,ctx->omega[p_id],T,T);
139: } else {
140: NEPComputeFunctionSubcomm(nep,ctx->omega[p_id],T,PETSC_FALSE);
141: }
142: MatDuplicate(T,MAT_COPY_VALUES,&kspMat);
143: KSPSetOperators(contour->ksp[i],kspMat,kspMat);
144: MatDestroy(&kspMat);
145: }
146: if (contour->subcomm->n==1) {
147: NEPComputeJacobian(nep,ctx->omega[p_id],dT);
148: } else {
149: NEPComputeFunctionSubcomm(nep,ctx->omega[p_id],dT,PETSC_TRUE);
150: }
151: BVSetActiveColumns(ctx->Y,i*ctx->L_max+L_start,i*ctx->L_max+L_end);
152: BVGetMat(ctx->Y,&MC);
153: if (!i) {
154: MatProductCreate(dT,MV,NULL,&BMV);
155: MatProductSetType(BMV,MATPRODUCT_AB);
156: MatProductSetFromOptions(BMV);
157: MatProductSymbolic(BMV);
158: }
159: MatProductNumeric(BMV);
160: KSPMatSolve(contour->ksp[i],BMV,MC);
161: BVRestoreMat(ctx->Y,&MC);
162: }
163: MatDestroy(&BMV);
164: BVRestoreMat(V,&MV);
165: return(0);
166: }
168: PetscErrorCode NEPSetUp_CISS(NEP nep)
169: {
170: PetscErrorCode ierr;
171: NEP_CISS *ctx = (NEP_CISS*)nep->data;
172: SlepcContourData contour;
173: PetscInt nwork;
174: PetscBool istrivial,isellipse,flg;
175: NEP_CISS_PROJECT dsctxf;
176: PetscObjectId id;
177: PetscObjectState state;
178: Vec v0;
181: if (nep->ncv==PETSC_DEFAULT) nep->ncv = ctx->L_max*ctx->M;
182: else {
183: ctx->L_max = nep->ncv/ctx->M;
184: if (!ctx->L_max) {
185: ctx->L_max = 1;
186: nep->ncv = ctx->L_max*ctx->M;
187: }
188: }
189: ctx->L = PetscMin(ctx->L,ctx->L_max);
190: if (nep->max_it==PETSC_DEFAULT) nep->max_it = 5;
191: if (nep->mpd==PETSC_DEFAULT) nep->mpd = nep->ncv;
192: if (!nep->which) nep->which = NEP_ALL;
193: if (nep->which!=NEP_ALL) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver supports only computing all eigenvalues");
194: NEPCheckUnsupported(nep,NEP_FEATURE_STOPPING | NEP_FEATURE_TWOSIDED);
196: /* check region */
197: RGIsTrivial(nep->rg,&istrivial);
198: if (istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
199: RGGetComplement(nep->rg,&flg);
200: if (flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
201: PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse);
202: if (!isellipse) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Currently only implemented for elliptic regions");
204: /* if the region has changed, then reset contour data */
205: PetscObjectGetId((PetscObject)nep->rg,&id);
206: PetscObjectStateGet((PetscObject)nep->rg,&state);
207: if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
208: SlepcContourDataDestroy(&ctx->contour);
209: PetscInfo(nep,"Resetting the contour data structure due to a change of region\n");
210: ctx->rgid = id; ctx->rgstate = state;
211: }
213: /* create contour data structure */
214: if (!ctx->contour) {
215: RGCanUseConjugates(nep->rg,ctx->isreal,&ctx->useconj);
216: SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)nep,&ctx->contour);
217: }
219: NEPAllocateSolution(nep,0);
220: if (ctx->weight) { PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma); }
221: PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma);
222: PetscLogObjectMemory((PetscObject)nep,3*ctx->N*sizeof(PetscScalar)+ctx->L_max*ctx->N*sizeof(PetscReal));
224: /* allocate basis vectors */
225: BVDestroy(&ctx->S);
226: BVDuplicateResize(nep->V,ctx->L_max*ctx->M,&ctx->S);
227: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->S);
228: BVDestroy(&ctx->V);
229: BVDuplicateResize(nep->V,ctx->L_max,&ctx->V);
230: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->V);
232: contour = ctx->contour;
233: SlepcContourRedundantMat(contour,nep->nt,nep->A);
234: if (contour->pA) {
235: if (!ctx->T) {
236: MatDuplicate(contour->pA[0],MAT_DO_NOT_COPY_VALUES,&ctx->T);
237: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->T);
238: }
239: if (!ctx->J) {
240: MatDuplicate(contour->pA[0],MAT_DO_NOT_COPY_VALUES,&ctx->J);
241: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->J);
242: }
243: BVGetColumn(ctx->V,0,&v0);
244: SlepcContourScatterCreate(contour,v0);
245: BVRestoreColumn(ctx->V,0,&v0);
246: BVDestroy(&ctx->pV);
247: BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV);
248: BVSetSizesFromVec(ctx->pV,contour->xsub,nep->n);
249: BVSetFromOptions(ctx->pV);
250: BVResize(ctx->pV,ctx->L_max,PETSC_FALSE);
251: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pV);
252: }
254: BVDestroy(&ctx->Y);
255: if (contour->pA) {
256: BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y);
257: BVSetSizesFromVec(ctx->Y,contour->xsub,nep->n);
258: BVSetFromOptions(ctx->Y);
259: BVResize(ctx->Y,contour->npoints*ctx->L_max,PETSC_FALSE);
260: } else {
261: BVDuplicateResize(nep->V,contour->npoints*ctx->L_max,&ctx->Y);
262: }
264: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
265: DSSetType(nep->ds,DSGNHEP);
266: } else if (ctx->extraction == NEP_CISS_EXTRACTION_CAA) {
267: DSSetType(nep->ds,DSNHEP);
268: } else {
269: DSSetType(nep->ds,DSNEP);
270: DSSetMethod(nep->ds,1);
271: DSNEPSetRG(nep->ds,nep->rg);
272: if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
273: DSNEPSetFN(nep->ds,nep->nt,nep->f);
274: } else {
275: PetscNew(&dsctxf);
276: DSNEPSetComputeMatrixFunction(nep->ds,NEPContourDSComputeMatrix,dsctxf);
277: dsctxf->nep = nep;
278: ctx->dsctxf = dsctxf;
279: }
280: }
281: DSAllocate(nep->ds,nep->ncv);
282: nwork = (nep->fui==NEP_USER_INTERFACE_SPLIT)? 2: 1;
283: NEPSetWorkVecs(nep,nwork);
284: return(0);
285: }
287: PetscErrorCode NEPSolve_CISS(NEP nep)
288: {
289: PetscErrorCode ierr;
290: NEP_CISS *ctx = (NEP_CISS*)nep->data;
291: SlepcContourData contour = ctx->contour;
292: Mat X,M,E;
293: PetscInt i,j,ld,L_add=0,nv=0,L_base=ctx->L,inner,*inside;
294: PetscScalar *Mu,*H0,*H1,*rr,*temp,center;
295: PetscReal error,max_error,radius,rgscale,est_eig,eta;
296: PetscBool isellipse,*fl1;
297: Vec si;
298: SlepcSC sc;
299: PetscRandom rand;
302: DSSetFromOptions(nep->ds);
303: DSGetSlepcSC(nep->ds,&sc);
304: sc->comparison = SlepcCompareLargestMagnitude;
305: sc->comparisonctx = NULL;
306: sc->map = NULL;
307: sc->mapobj = NULL;
308: DSGetLeadingDimension(nep->ds,&ld);
309: RGComputeQuadrature(nep->rg,RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight);
310: BVSetActiveColumns(ctx->V,0,ctx->L);
311: BVSetRandomSign(ctx->V);
312: BVGetRandomContext(ctx->V,&rand);
313: if (contour->pA) {
314: BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
315: NEPCISSSolveSystem(nep,ctx->T,ctx->J,ctx->pV,0,ctx->L,PETSC_TRUE);
316: } else {
317: NEPCISSSolveSystem(nep,nep->function,nep->jacobian,ctx->V,0,ctx->L,PETSC_TRUE);
318: }
319: PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse);
320: if (isellipse) {
321: BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L_max,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig);
322: PetscInfo1(nep,"Estimated eigenvalue count: %f\n",(double)est_eig);
323: eta = PetscPowReal(10.0,-PetscLog10Real(nep->tol)/ctx->N);
324: L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
325: if (L_add>ctx->L_max-ctx->L) {
326: PetscInfo(nep,"Number of eigenvalues inside the contour path may be too large\n");
327: L_add = ctx->L_max-ctx->L;
328: }
329: }
330: /* Updates L after estimate the number of eigenvalue */
331: if (L_add>0) {
332: PetscInfo2(nep,"Changing L %D -> %D by Estimate #Eig\n",ctx->L,ctx->L+L_add);
333: BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
334: BVSetRandomSign(ctx->V);
335: if (contour->pA) {
336: BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
337: NEPCISSSolveSystem(nep,ctx->T,ctx->J,ctx->pV,ctx->L,ctx->L+L_add,PETSC_FALSE);
338: } else {
339: NEPCISSSolveSystem(nep,nep->function,nep->jacobian,ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
340: }
341: ctx->L += L_add;
342: }
344: PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
345: for (i=0;i<ctx->refine_blocksize;i++) {
346: BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj);
347: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
348: PetscLogEventBegin(NEP_CISS_SVD,nep,0,0,0);
349: SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
350: PetscLogEventEnd(NEP_CISS_SVD,nep,0,0,0);
351: if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
352: L_add = L_base;
353: if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
354: PetscInfo2(nep,"Changing L %D -> %D by SVD(H0)\n",ctx->L,ctx->L+L_add);
355: BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
356: BVSetRandomSign(ctx->V);
357: if (contour->pA) {
358: BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
359: NEPCISSSolveSystem(nep,ctx->T,ctx->J,ctx->pV,ctx->L,ctx->L+L_add,PETSC_FALSE);
360: } else {
361: NEPCISSSolveSystem(nep,nep->function,nep->jacobian,ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
362: }
363: ctx->L += L_add;
364: if (L_add) {
365: PetscFree2(Mu,H0);
366: PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
367: }
368: }
370: RGGetScale(nep->rg,&rgscale);
371: RGEllipseGetParameters(nep->rg,¢er,&radius,NULL);
373: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
374: PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1);
375: }
377: while (nep->reason == NEP_CONVERGED_ITERATING) {
378: nep->its++;
379: for (inner=0;inner<=ctx->refine_inner;inner++) {
380: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
381: BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj);
382: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
383: PetscLogEventBegin(NEP_CISS_SVD,nep,0,0,0);
384: SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
385: PetscLogEventEnd(NEP_CISS_SVD,nep,0,0,0);
386: } else {
387: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
388: /* compute SVD of S */
389: BVSVDAndRank(ctx->S,ctx->M,ctx->L,ctx->delta,(ctx->extraction==NEP_CISS_EXTRACTION_CAA)?BV_SVD_METHOD_QR_CAA:BV_SVD_METHOD_QR,H0,ctx->sigma,&nv);
390: }
391: if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
392: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
393: BVSetActiveColumns(ctx->S,0,ctx->L);
394: BVSetActiveColumns(ctx->V,0,ctx->L);
395: BVCopy(ctx->S,ctx->V);
396: if (contour->pA) {
397: BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
398: NEPCISSSolveSystem(nep,ctx->T,ctx->J,ctx->pV,0,ctx->L,PETSC_FALSE);
399: } else {
400: NEPCISSSolveSystem(nep,nep->function,nep->jacobian,ctx->V,0,ctx->L,PETSC_FALSE);
401: }
402: } else break;
403: }
404: nep->nconv = 0;
405: if (nv == 0) { nep->reason = NEP_CONVERGED_TOL; break; }
406: else {
407: /* Extracting eigenpairs */
408: DSSetDimensions(nep->ds,nv,0,0);
409: DSSetState(nep->ds,DS_STATE_RAW);
410: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
411: CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
412: CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1);
413: DSGetArray(nep->ds,DS_MAT_A,&temp);
414: for (j=0;j<nv;j++)
415: for (i=0;i<nv;i++)
416: temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
417: DSRestoreArray(nep->ds,DS_MAT_A,&temp);
418: DSGetArray(nep->ds,DS_MAT_B,&temp);
419: for (j=0;j<nv;j++)
420: for (i=0;i<nv;i++)
421: temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
422: DSRestoreArray(nep->ds,DS_MAT_B,&temp);
423: } else if (ctx->extraction == NEP_CISS_EXTRACTION_CAA) {
424: BVSetActiveColumns(ctx->S,0,nv);
425: DSGetArray(nep->ds,DS_MAT_A,&temp);
426: for (i=0;i<nv;i++) {
427: PetscArraycpy(temp+i*ld,H0+i*nv,nv);
428: }
429: DSRestoreArray(nep->ds,DS_MAT_A,&temp);
430: } else {
431: BVSetActiveColumns(ctx->S,0,nv);
432: if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
433: for (i=0;i<nep->nt;i++) {
434: DSGetMat(nep->ds,DSMatExtra[i],&E);
435: BVMatProject(ctx->S,nep->A[i],ctx->S,E);
436: DSRestoreMat(nep->ds,DSMatExtra[i],&E);
437: }
438: } else { ctx->dsctxf->Q = ctx->S; }
439: }
440: DSSolve(nep->ds,nep->eigr,nep->eigi);
441: DSSynchronize(nep->ds,nep->eigr,nep->eigi);
442: DSGetDimensions(nep->ds,NULL,NULL,NULL,&nv);
443: if (ctx->extraction == NEP_CISS_EXTRACTION_CAA || ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
444: for (i=0;i<nv;i++) {
445: nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
446: }
447: }
448: PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr);
449: DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
450: DSGetMat(nep->ds,DS_MAT_X,&X);
451: SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1);
452: MatDestroy(&X);
453: RGCheckInside(nep->rg,nv,nep->eigr,nep->eigi,inside);
454: for (i=0;i<nv;i++) {
455: if (fl1[i] && inside[i]>=0) {
456: rr[i] = 1.0;
457: nep->nconv++;
458: } else rr[i] = 0.0;
459: }
460: DSSort(nep->ds,nep->eigr,nep->eigi,rr,NULL,&nep->nconv);
461: DSSynchronize(nep->ds,nep->eigr,nep->eigi);
462: if (ctx->extraction == NEP_CISS_EXTRACTION_CAA || ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
463: for (i=0;i<nv;i++) nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
464: }
465: PetscFree3(fl1,inside,rr);
466: BVSetActiveColumns(nep->V,0,nv);
467: DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
468: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
469: BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
470: BVSetActiveColumns(ctx->S,0,nv);
471: BVCopy(ctx->S,nep->V);
472: DSGetMat(nep->ds,DS_MAT_X,&X);
473: BVMultInPlace(ctx->S,X,0,nep->nconv);
474: BVMultInPlace(nep->V,X,0,nep->nconv);
475: MatDestroy(&X);
476: } else {
477: DSGetMat(nep->ds,DS_MAT_X,&X);
478: BVMultInPlace(ctx->S,X,0,nep->nconv);
479: MatDestroy(&X);
480: BVSetActiveColumns(ctx->S,0,nv);
481: BVCopy(ctx->S,nep->V);
482: }
483: max_error = 0.0;
484: for (i=0;i<nep->nconv;i++) {
485: BVGetColumn(nep->V,i,&si);
486: VecNormalize(si,NULL);
487: NEPComputeResidualNorm_Private(nep,PETSC_FALSE,nep->eigr[i],si,nep->work,&error);
488: (*nep->converged)(nep,nep->eigr[i],0,error,&error,nep->convergedctx);
489: BVRestoreColumn(nep->V,i,&si);
490: max_error = PetscMax(max_error,error);
491: }
492: if (max_error <= nep->tol) nep->reason = NEP_CONVERGED_TOL;
493: else if (nep->its > nep->max_it) nep->reason = NEP_DIVERGED_ITS;
494: else {
495: if (nep->nconv > ctx->L) nv = nep->nconv;
496: else if (ctx->L > nv) nv = ctx->L;
497: MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M);
498: MatSetRandom(M,rand);
499: BVSetActiveColumns(ctx->S,0,nv);
500: BVMultInPlace(ctx->S,M,0,ctx->L);
501: MatDestroy(&M);
502: BVSetActiveColumns(ctx->S,0,ctx->L);
503: BVSetActiveColumns(ctx->V,0,ctx->L);
504: BVCopy(ctx->S,ctx->V);
505: if (contour->pA) {
506: BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
507: NEPCISSSolveSystem(nep,ctx->T,ctx->J,ctx->pV,0,ctx->L,PETSC_FALSE);
508: } else {
509: NEPCISSSolveSystem(nep,nep->function,nep->jacobian,ctx->V,0,ctx->L,PETSC_FALSE);
510: }
511: }
512: }
513: }
514: PetscFree2(Mu,H0);
515: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
516: PetscFree(H1);
517: }
518: return(0);
519: }
521: static PetscErrorCode NEPCISSSetSizes_CISS(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
522: {
524: NEP_CISS *ctx = (NEP_CISS*)nep->data;
525: PetscInt oN,oL,oM,oLmax,onpart;
528: oN = ctx->N;
529: if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
530: if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
531: } else {
532: if (ip<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
533: if (ip%2) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
534: if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
535: }
536: oL = ctx->L;
537: if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
538: ctx->L = 16;
539: } else {
540: if (bs<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
541: ctx->L = bs;
542: }
543: oM = ctx->M;
544: if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
545: ctx->M = ctx->N/4;
546: } else {
547: if (ms<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
548: if (ms>ctx->N) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be less than or equal to the number of integration points");
549: ctx->M = PetscMax(ms,2);
550: }
551: onpart = ctx->npart;
552: if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
553: ctx->npart = 1;
554: } else {
555: if (npart<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The npart argument must be > 0");
556: ctx->npart = npart;
557: }
558: oLmax = ctx->L_max;
559: if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
560: ctx->L_max = 64;
561: } else {
562: if (bsmax<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
563: ctx->L_max = PetscMax(bsmax,ctx->L);
564: }
565: if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
566: SlepcContourDataDestroy(&ctx->contour);
567: PetscInfo(nep,"Resetting the contour data structure due to a change of parameters\n");
568: nep->state = NEP_STATE_INITIAL;
569: }
570: ctx->isreal = realmats;
571: if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) nep->state = NEP_STATE_INITIAL;
572: return(0);
573: }
575: /*@
576: NEPCISSSetSizes - Sets the values of various size parameters in the CISS solver.
578: Logically Collective on nep
580: Input Parameters:
581: + nep - the nonlinear eigensolver context
582: . ip - number of integration points
583: . bs - block size
584: . ms - moment size
585: . npart - number of partitions when splitting the communicator
586: . bsmax - max block size
587: - realmats - T(z) is real for real z
589: Options Database Keys:
590: + -nep_ciss_integration_points - Sets the number of integration points
591: . -nep_ciss_blocksize - Sets the block size
592: . -nep_ciss_moments - Sets the moment size
593: . -nep_ciss_partitions - Sets the number of partitions
594: . -nep_ciss_maxblocksize - Sets the maximum block size
595: - -nep_ciss_realmats - T(z) is real for real z
597: Notes:
598: The default number of partitions is 1. This means the internal KSP object is shared
599: among all processes of the NEP communicator. Otherwise, the communicator is split
600: into npart communicators, so that npart KSP solves proceed simultaneously.
602: The realmats flag can be set to true when T(.) is guaranteed to be real
603: when the argument is a real value, for example, when all matrices in
604: the split form are real. When set to true, the solver avoids some computations.
606: Level: advanced
608: .seealso: NEPCISSGetSizes()
609: @*/
610: PetscErrorCode NEPCISSSetSizes(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
611: {
622: PetscTryMethod(nep,"NEPCISSSetSizes_C",(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(nep,ip,bs,ms,npart,bsmax,realmats));
623: return(0);
624: }
626: static PetscErrorCode NEPCISSGetSizes_CISS(NEP nep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
627: {
628: NEP_CISS *ctx = (NEP_CISS*)nep->data;
631: if (ip) *ip = ctx->N;
632: if (bs) *bs = ctx->L;
633: if (ms) *ms = ctx->M;
634: if (npart) *npart = ctx->npart;
635: if (bsmax) *bsmax = ctx->L_max;
636: if (realmats) *realmats = ctx->isreal;
637: return(0);
638: }
640: /*@
641: NEPCISSGetSizes - Gets the values of various size parameters in the CISS solver.
643: Not Collective
645: Input Parameter:
646: . nep - the nonlinear eigensolver context
648: Output Parameters:
649: + ip - number of integration points
650: . bs - block size
651: . ms - moment size
652: . npart - number of partitions when splitting the communicator
653: . bsmax - max block size
654: - realmats - T(z) is real for real z
656: Level: advanced
658: .seealso: NEPCISSSetSizes()
659: @*/
660: PetscErrorCode NEPCISSGetSizes(NEP nep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
661: {
666: PetscUseMethod(nep,"NEPCISSGetSizes_C",(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(nep,ip,bs,ms,npart,bsmax,realmats));
667: return(0);
668: }
670: static PetscErrorCode NEPCISSSetThreshold_CISS(NEP nep,PetscReal delta,PetscReal spur)
671: {
672: NEP_CISS *ctx = (NEP_CISS*)nep->data;
675: if (delta == PETSC_DEFAULT) {
676: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
677: } else {
678: if (delta<=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
679: ctx->delta = delta;
680: }
681: if (spur == PETSC_DEFAULT) {
682: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
683: } else {
684: if (spur<=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
685: ctx->spurious_threshold = spur;
686: }
687: return(0);
688: }
690: /*@
691: NEPCISSSetThreshold - Sets the values of various threshold parameters in
692: the CISS solver.
694: Logically Collective on nep
696: Input Parameters:
697: + nep - the nonlinear eigensolver context
698: . delta - threshold for numerical rank
699: - spur - spurious threshold (to discard spurious eigenpairs)
701: Options Database Keys:
702: + -nep_ciss_delta - Sets the delta
703: - -nep_ciss_spurious_threshold - Sets the spurious threshold
705: Level: advanced
707: .seealso: NEPCISSGetThreshold()
708: @*/
709: PetscErrorCode NEPCISSSetThreshold(NEP nep,PetscReal delta,PetscReal spur)
710: {
717: PetscTryMethod(nep,"NEPCISSSetThreshold_C",(NEP,PetscReal,PetscReal),(nep,delta,spur));
718: return(0);
719: }
721: static PetscErrorCode NEPCISSGetThreshold_CISS(NEP nep,PetscReal *delta,PetscReal *spur)
722: {
723: NEP_CISS *ctx = (NEP_CISS*)nep->data;
726: if (delta) *delta = ctx->delta;
727: if (spur) *spur = ctx->spurious_threshold;
728: return(0);
729: }
731: /*@
732: NEPCISSGetThreshold - Gets the values of various threshold parameters in
733: the CISS solver.
735: Not Collective
737: Input Parameter:
738: . nep - the nonlinear eigensolver context
740: Output Parameters:
741: + delta - threshold for numerical rank
742: - spur - spurious threshold (to discard spurious eigenpairs)
744: Level: advanced
746: .seealso: NEPCISSSetThreshold()
747: @*/
748: PetscErrorCode NEPCISSGetThreshold(NEP nep,PetscReal *delta,PetscReal *spur)
749: {
754: PetscUseMethod(nep,"NEPCISSGetThreshold_C",(NEP,PetscReal*,PetscReal*),(nep,delta,spur));
755: return(0);
756: }
758: static PetscErrorCode NEPCISSSetRefinement_CISS(NEP nep,PetscInt inner,PetscInt blsize)
759: {
760: NEP_CISS *ctx = (NEP_CISS*)nep->data;
763: if (inner == PETSC_DEFAULT) {
764: ctx->refine_inner = 0;
765: } else {
766: if (inner<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
767: ctx->refine_inner = inner;
768: }
769: if (blsize == PETSC_DEFAULT) {
770: ctx->refine_blocksize = 0;
771: } else {
772: if (blsize<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
773: ctx->refine_blocksize = blsize;
774: }
775: return(0);
776: }
778: /*@
779: NEPCISSSetRefinement - Sets the values of various refinement parameters
780: in the CISS solver.
782: Logically Collective on nep
784: Input Parameters:
785: + nep - the nonlinear eigensolver context
786: . inner - number of iterative refinement iterations (inner loop)
787: - blsize - number of iterative refinement iterations (blocksize loop)
789: Options Database Keys:
790: + -nep_ciss_refine_inner - Sets number of inner iterations
791: - -nep_ciss_refine_blocksize - Sets number of blocksize iterations
793: Level: advanced
795: .seealso: NEPCISSGetRefinement()
796: @*/
797: PetscErrorCode NEPCISSSetRefinement(NEP nep,PetscInt inner,PetscInt blsize)
798: {
805: PetscTryMethod(nep,"NEPCISSSetRefinement_C",(NEP,PetscInt,PetscInt),(nep,inner,blsize));
806: return(0);
807: }
809: static PetscErrorCode NEPCISSGetRefinement_CISS(NEP nep,PetscInt *inner,PetscInt *blsize)
810: {
811: NEP_CISS *ctx = (NEP_CISS*)nep->data;
814: if (inner) *inner = ctx->refine_inner;
815: if (blsize) *blsize = ctx->refine_blocksize;
816: return(0);
817: }
819: /*@
820: NEPCISSGetRefinement - Gets the values of various refinement parameters
821: in the CISS solver.
823: Not Collective
825: Input Parameter:
826: . nep - the nonlinear eigensolver context
828: Output Parameters:
829: + inner - number of iterative refinement iterations (inner loop)
830: - blsize - number of iterative refinement iterations (blocksize loop)
832: Level: advanced
834: .seealso: NEPCISSSetRefinement()
835: @*/
836: PetscErrorCode NEPCISSGetRefinement(NEP nep, PetscInt *inner, PetscInt *blsize)
837: {
842: PetscUseMethod(nep,"NEPCISSGetRefinement_C",(NEP,PetscInt*,PetscInt*),(nep,inner,blsize));
843: return(0);
844: }
846: static PetscErrorCode NEPCISSSetExtraction_CISS(NEP nep,NEPCISSExtraction extraction)
847: {
848: NEP_CISS *ctx = (NEP_CISS*)nep->data;
851: if (ctx->extraction != extraction) {
852: ctx->extraction = extraction;
853: nep->state = NEP_STATE_INITIAL;
854: }
855: return(0);
856: }
858: /*@
859: NEPCISSSetExtraction - Sets the extraction technique used in the CISS solver.
861: Logically Collective on nep
863: Input Parameters:
864: + nep - the nonlinear eigensolver context
865: - extraction - the extraction technique
867: Options Database Key:
868: . -nep_ciss_extraction - Sets the extraction technique (either 'ritz', 'hankel' or 'caa')
870: Notes:
871: By default, the Rayleigh-Ritz extraction is used (NEP_CISS_EXTRACTION_RITZ).
873: If the 'hankel' or the 'caa' option is specified (NEP_CISS_EXTRACTION_HANKEL or
874: NEP_CISS_EXTRACTION_CAA), then the Block Hankel method, or the Communication-avoiding
875: Arnoldi method, respectively, is used for extracting eigenpairs.
877: Level: advanced
879: .seealso: NEPCISSGetExtraction(), NEPCISSExtraction
880: @*/
881: PetscErrorCode NEPCISSSetExtraction(NEP nep,NEPCISSExtraction extraction)
882: {
888: PetscTryMethod(nep,"NEPCISSSetExtraction_C",(NEP,NEPCISSExtraction),(nep,extraction));
889: return(0);
890: }
892: static PetscErrorCode NEPCISSGetExtraction_CISS(NEP nep,NEPCISSExtraction *extraction)
893: {
894: NEP_CISS *ctx = (NEP_CISS*)nep->data;
897: *extraction = ctx->extraction;
898: return(0);
899: }
901: /*@
902: NEPCISSGetExtraction - Gets the extraction technique used in the CISS solver.
904: Not Collective
906: Input Parameter:
907: . nep - the nonlinear eigensolver context
909: Output Parameters:
910: . extraction - extraction technique
912: Level: advanced
914: .seealso: NEPCISSSetExtraction() NEPCISSExtraction
915: @*/
916: PetscErrorCode NEPCISSGetExtraction(NEP nep,NEPCISSExtraction *extraction)
917: {
923: PetscUseMethod(nep,"NEPCISSGetExtraction_C",(NEP,NEPCISSExtraction*),(nep,extraction));
924: return(0);
925: }
927: static PetscErrorCode NEPCISSGetKSPs_CISS(NEP nep,PetscInt *nsolve,KSP **ksp)
928: {
929: PetscErrorCode ierr;
930: NEP_CISS *ctx = (NEP_CISS*)nep->data;
931: SlepcContourData contour;
932: PetscInt i;
933: PC pc;
936: if (!ctx->contour) { /* initialize contour data structure first */
937: RGCanUseConjugates(nep->rg,ctx->isreal,&ctx->useconj);
938: SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)nep,&ctx->contour);
939: }
940: contour = ctx->contour;
941: if (!contour->ksp) {
942: PetscMalloc1(contour->npoints,&contour->ksp);
943: for (i=0;i<contour->npoints;i++) {
944: KSPCreate(PetscSubcommChild(contour->subcomm),&contour->ksp[i]);
945: PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)nep,1);
946: KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)nep)->prefix);
947: KSPAppendOptionsPrefix(contour->ksp[i],"nep_ciss_");
948: PetscLogObjectParent((PetscObject)nep,(PetscObject)contour->ksp[i]);
949: PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)nep)->options);
950: KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE);
951: KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(nep->tol),PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
952: KSPGetPC(contour->ksp[i],&pc);
953: KSPSetType(contour->ksp[i],KSPPREONLY);
954: PCSetType(pc,PCLU);
955: }
956: }
957: if (nsolve) *nsolve = contour->npoints;
958: if (ksp) *ksp = contour->ksp;
959: return(0);
960: }
962: /*@C
963: NEPCISSGetKSPs - Retrieve the array of linear solver objects associated with
964: the CISS solver.
966: Not Collective
968: Input Parameter:
969: . nep - nonlinear eigenvalue solver
971: Output Parameters:
972: + nsolve - number of solver objects
973: - ksp - array of linear solver object
975: Notes:
976: The number of KSP solvers is equal to the number of integration points divided by
977: the number of partitions. This value is halved in the case of real matrices with
978: a region centered at the real axis.
980: Level: advanced
982: .seealso: NEPCISSSetSizes()
983: @*/
984: PetscErrorCode NEPCISSGetKSPs(NEP nep,PetscInt *nsolve,KSP **ksp)
985: {
990: PetscUseMethod(nep,"NEPCISSGetKSPs_C",(NEP,PetscInt*,KSP**),(nep,nsolve,ksp));
991: return(0);
992: }
994: PetscErrorCode NEPReset_CISS(NEP nep)
995: {
997: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1000: BVDestroy(&ctx->S);
1001: BVDestroy(&ctx->V);
1002: BVDestroy(&ctx->Y);
1003: SlepcContourDataReset(ctx->contour);
1004: MatDestroy(&ctx->T);
1005: MatDestroy(&ctx->J);
1006: BVDestroy(&ctx->pV);
1007: if (ctx->extraction == NEP_CISS_EXTRACTION_RITZ && nep->fui==NEP_USER_INTERFACE_CALLBACK) {
1008: PetscFree(ctx->dsctxf);
1009: }
1010: return(0);
1011: }
1013: PetscErrorCode NEPSetFromOptions_CISS(PetscOptionItems *PetscOptionsObject,NEP nep)
1014: {
1015: PetscErrorCode ierr;
1016: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1017: PetscReal r1,r2;
1018: PetscInt i,i1,i2,i3,i4,i5,i6,i7;
1019: PetscBool b1,flg,flg2,flg3,flg4,flg5,flg6;
1020: NEPCISSExtraction extraction;
1023: PetscOptionsHead(PetscOptionsObject,"NEP CISS Options");
1025: NEPCISSGetSizes(nep,&i1,&i2,&i3,&i4,&i5,&b1);
1026: PetscOptionsInt("-nep_ciss_integration_points","Number of integration points","NEPCISSSetSizes",i1,&i1,&flg);
1027: PetscOptionsInt("-nep_ciss_blocksize","Block size","NEPCISSSetSizes",i2,&i2,&flg2);
1028: PetscOptionsInt("-nep_ciss_moments","Moment size","NEPCISSSetSizes",i3,&i3,&flg3);
1029: PetscOptionsInt("-nep_ciss_partitions","Number of partitions","NEPCISSSetSizes",i4,&i4,&flg4);
1030: PetscOptionsInt("-nep_ciss_maxblocksize","Maximum block size","NEPCISSSetSizes",i5,&i5,&flg5);
1031: PetscOptionsBool("-nep_ciss_realmats","True if T(z) is real for real z","NEPCISSSetSizes",b1,&b1,&flg6);
1032: if (flg || flg2 || flg3 || flg4 || flg5 || flg6) { NEPCISSSetSizes(nep,i1,i2,i3,i4,i5,b1); }
1034: NEPCISSGetThreshold(nep,&r1,&r2);
1035: PetscOptionsReal("-nep_ciss_delta","Threshold for numerical rank","NEPCISSSetThreshold",r1,&r1,&flg);
1036: PetscOptionsReal("-nep_ciss_spurious_threshold","Threshold for the spurious eigenpairs","NEPCISSSetThreshold",r2,&r2,&flg2);
1037: if (flg || flg2) { NEPCISSSetThreshold(nep,r1,r2); }
1039: NEPCISSGetRefinement(nep,&i6,&i7);
1040: PetscOptionsInt("-nep_ciss_refine_inner","Number of inner iterative refinement iterations","NEPCISSSetRefinement",i6,&i6,&flg);
1041: PetscOptionsInt("-nep_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","NEPCISSSetRefinement",i7,&i7,&flg2);
1042: if (flg || flg2) { NEPCISSSetRefinement(nep,i6,i7); }
1044: PetscOptionsEnum("-nep_ciss_extraction","Extraction technique","NEPCISSSetExtraction",NEPCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg);
1045: if (flg) { NEPCISSSetExtraction(nep,extraction); }
1047: PetscOptionsTail();
1049: if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
1050: RGSetFromOptions(nep->rg); /* this is necessary here to set useconj */
1051: if (!ctx->contour || !ctx->contour->ksp) { NEPCISSGetKSPs(nep,NULL,NULL); }
1052: for (i=0;i<ctx->contour->npoints;i++) {
1053: KSPSetFromOptions(ctx->contour->ksp[i]);
1054: }
1055: PetscSubcommSetFromOptions(ctx->contour->subcomm);
1056: return(0);
1057: }
1059: PetscErrorCode NEPDestroy_CISS(NEP nep)
1060: {
1062: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1065: SlepcContourDataDestroy(&ctx->contour);
1066: PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
1067: PetscFree(nep->data);
1068: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NULL);
1069: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NULL);
1070: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NULL);
1071: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NULL);
1072: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NULL);
1073: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NULL);
1074: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetExtraction_C",NULL);
1075: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetExtraction_C",NULL);
1076: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetKSPs_C",NULL);
1077: return(0);
1078: }
1080: PetscErrorCode NEPView_CISS(NEP nep,PetscViewer viewer)
1081: {
1083: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1084: PetscBool isascii;
1085: PetscViewer sviewer;
1088: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1089: if (isascii) {
1090: PetscViewerASCIIPrintf(viewer," sizes { integration points: %D, block size: %D, moment size: %D, partitions: %D, maximum block size: %D }\n",ctx->N,ctx->L,ctx->M,ctx->npart,ctx->L_max);
1091: if (ctx->isreal) {
1092: PetscViewerASCIIPrintf(viewer," exploiting symmetry of integration points\n");
1093: }
1094: PetscViewerASCIIPrintf(viewer," threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold);
1095: PetscViewerASCIIPrintf(viewer," iterative refinement { inner: %D, blocksize: %D }\n",ctx->refine_inner, ctx->refine_blocksize);
1096: PetscViewerASCIIPrintf(viewer," extraction: %s\n",NEPCISSExtractions[ctx->extraction]);
1097: if (!ctx->contour || !ctx->contour->ksp) { NEPCISSGetKSPs(nep,NULL,NULL); }
1098: PetscViewerASCIIPushTab(viewer);
1099: if (ctx->npart>1 && ctx->contour->subcomm) {
1100: PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1101: if (!ctx->contour->subcomm->color) {
1102: KSPView(ctx->contour->ksp[0],sviewer);
1103: }
1104: PetscViewerFlush(sviewer);
1105: PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1106: PetscViewerFlush(viewer);
1107: /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1108: PetscViewerASCIIPopSynchronized(viewer);
1109: } else {
1110: KSPView(ctx->contour->ksp[0],viewer);
1111: }
1112: PetscViewerASCIIPopTab(viewer);
1113: }
1114: return(0);
1115: }
1117: SLEPC_EXTERN PetscErrorCode NEPCreate_CISS(NEP nep)
1118: {
1120: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1123: PetscNewLog(nep,&ctx);
1124: nep->data = ctx;
1125: /* set default values of parameters */
1126: ctx->N = 32;
1127: ctx->L = 16;
1128: ctx->M = ctx->N/4;
1129: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
1130: ctx->L_max = 64;
1131: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1132: ctx->isreal = PETSC_FALSE;
1133: ctx->npart = 1;
1135: nep->useds = PETSC_TRUE;
1137: nep->ops->solve = NEPSolve_CISS;
1138: nep->ops->setup = NEPSetUp_CISS;
1139: nep->ops->setfromoptions = NEPSetFromOptions_CISS;
1140: nep->ops->reset = NEPReset_CISS;
1141: nep->ops->destroy = NEPDestroy_CISS;
1142: nep->ops->view = NEPView_CISS;
1144: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NEPCISSSetSizes_CISS);
1145: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NEPCISSGetSizes_CISS);
1146: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NEPCISSSetThreshold_CISS);
1147: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NEPCISSGetThreshold_CISS);
1148: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NEPCISSSetRefinement_CISS);
1149: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NEPCISSGetRefinement_CISS);
1150: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetExtraction_C",NEPCISSSetExtraction_CISS);
1151: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetExtraction_C",NEPCISSGetExtraction_CISS);
1152: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetKSPs_C",NEPCISSGetKSPs_CISS);
1153: return(0);
1154: }