Actual source code: nciss.c
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 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 J; /* auxiliary matrix when using subcomm */
62: BV pV;
63: PetscObjectId rgid;
64: PetscObjectState rgstate;
65: } NEP_CISS;
67: struct _n_nep_ciss_project {
68: NEP nep;
69: BV Q;
70: };
72: static PetscErrorCode NEPContourDSComputeMatrix(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat,void *ctx)
73: {
74: NEP_CISS_PROJECT proj = (NEP_CISS_PROJECT)ctx;
75: Mat M,fun;
77: PetscFunctionBegin;
78: if (!deriv) {
79: PetscCall(NEPComputeFunction(proj->nep,lambda,proj->nep->function,proj->nep->function));
80: fun = proj->nep->function;
81: } else {
82: PetscCall(NEPComputeJacobian(proj->nep,lambda,proj->nep->jacobian));
83: fun = proj->nep->jacobian;
84: }
85: PetscCall(DSGetMat(ds,mat,&M));
86: PetscCall(BVMatProject(proj->Q,fun,proj->Q,M));
87: PetscCall(DSRestoreMat(ds,mat,&M));
88: PetscFunctionReturn(PETSC_SUCCESS);
89: }
91: static PetscErrorCode NEPComputeFunctionSubcomm(NEP nep,PetscScalar lambda,Mat T,Mat P,PetscBool deriv)
92: {
93: PetscInt i;
94: PetscScalar alpha;
95: NEP_CISS *ctx = (NEP_CISS*)nep->data;
97: PetscFunctionBegin;
98: PetscAssert(nep->fui!=NEP_USER_INTERFACE_CALLBACK,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONGSTATE,"Should not arrive here with callbacks");
99: PetscCall(MatZeroEntries(T));
100: if (!deriv && T != P) PetscCall(MatZeroEntries(P));
101: for (i=0;i<nep->nt;i++) {
102: if (!deriv) PetscCall(FNEvaluateFunction(nep->f[i],lambda,&alpha));
103: else PetscCall(FNEvaluateDerivative(nep->f[i],lambda,&alpha));
104: PetscCall(MatAXPY(T,alpha,ctx->contour->pA[i],nep->mstr));
105: if (!deriv && T != P) PetscCall(MatAXPY(P,alpha,ctx->contour->pP[i],nep->mstrp));
106: }
107: PetscFunctionReturn(PETSC_SUCCESS);
108: }
110: /*
111: Set up KSP solvers for every integration point
112: */
113: static PetscErrorCode NEPCISSSetUp(NEP nep,Mat T,Mat P)
114: {
115: NEP_CISS *ctx = (NEP_CISS*)nep->data;
116: SlepcContourData contour;
117: PetscInt i,p_id;
118: Mat Amat,Pmat;
120: PetscFunctionBegin;
121: if (!ctx->contour || !ctx->contour->ksp) PetscCall(NEPCISSGetKSPs(nep,NULL,NULL));
122: contour = ctx->contour;
123: PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)nep),PETSC_ERR_PLIB,"Something went wrong with NEPCISSGetKSPs()");
124: for (i=0;i<contour->npoints;i++) {
125: p_id = i*contour->subcomm->n + contour->subcomm->color;
126: PetscCall(MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&Amat));
127: if (T != P) PetscCall(MatDuplicate(P,MAT_DO_NOT_COPY_VALUES,&Pmat)); else Pmat = Amat;
128: if (contour->subcomm->n == 1 || nep->fui==NEP_USER_INTERFACE_CALLBACK) PetscCall(NEPComputeFunction(nep,ctx->omega[p_id],Amat,Pmat));
129: else PetscCall(NEPComputeFunctionSubcomm(nep,ctx->omega[p_id],Amat,Pmat,PETSC_FALSE));
130: PetscCall(NEP_KSPSetOperators(contour->ksp[i],Amat,Pmat));
131: PetscCall(MatDestroy(&Amat));
132: if (T != P) PetscCall(MatDestroy(&Pmat));
133: }
134: PetscFunctionReturn(PETSC_SUCCESS);
135: }
137: /*
138: Y_i = F(z_i)^{-1}Fp(z_i)V for every integration point, Y=[Y_i] is in the context
139: */
140: static PetscErrorCode NEPCISSSolve(NEP nep,Mat dT,BV V,PetscInt L_start,PetscInt L_end)
141: {
142: NEP_CISS *ctx = (NEP_CISS*)nep->data;
143: SlepcContourData contour = ctx->contour;
144: PetscInt i,p_id;
145: Mat MV,BMV=NULL,MC;
147: PetscFunctionBegin;
148: PetscCall(BVSetActiveColumns(V,L_start,L_end));
149: PetscCall(BVGetMat(V,&MV));
150: for (i=0;i<contour->npoints;i++) {
151: p_id = i*contour->subcomm->n + contour->subcomm->color;
152: if (contour->subcomm->n==1 || nep->fui==NEP_USER_INTERFACE_CALLBACK) PetscCall(NEPComputeJacobian(nep,ctx->omega[p_id],dT));
153: else PetscCall(NEPComputeFunctionSubcomm(nep,ctx->omega[p_id],dT,NULL,PETSC_TRUE));
154: PetscCall(BVSetActiveColumns(ctx->Y,i*ctx->L+L_start,i*ctx->L+L_end));
155: PetscCall(BVGetMat(ctx->Y,&MC));
156: if (!i) {
157: PetscCall(MatProductCreate(dT,MV,NULL,&BMV));
158: PetscCall(MatProductSetType(BMV,MATPRODUCT_AB));
159: PetscCall(MatProductSetFromOptions(BMV));
160: PetscCall(MatProductSymbolic(BMV));
161: }
162: PetscCall(MatProductNumeric(BMV));
163: PetscCall(KSPMatSolve(contour->ksp[i],BMV,MC));
164: PetscCall(BVRestoreMat(ctx->Y,&MC));
165: }
166: PetscCall(MatDestroy(&BMV));
167: PetscCall(BVRestoreMat(V,&MV));
168: PetscFunctionReturn(PETSC_SUCCESS);
169: }
171: static PetscErrorCode NEPSetUp_CISS(NEP nep)
172: {
173: NEP_CISS *ctx = (NEP_CISS*)nep->data;
174: SlepcContourData contour;
175: PetscInt nwork;
176: PetscBool istrivial,isellipse,flg;
177: NEP_CISS_PROJECT dsctxf;
178: PetscObjectId id;
179: PetscObjectState state;
180: Vec v0;
182: PetscFunctionBegin;
183: if (nep->ncv==PETSC_DETERMINE) nep->ncv = ctx->L_max*ctx->M;
184: else {
185: ctx->L_max = nep->ncv/ctx->M;
186: if (!ctx->L_max) {
187: ctx->L_max = 1;
188: nep->ncv = ctx->L_max*ctx->M;
189: }
190: }
191: ctx->L = PetscMin(ctx->L,ctx->L_max);
192: if (nep->max_it==PETSC_DETERMINE) nep->max_it = 5;
193: if (nep->mpd==PETSC_DETERMINE) nep->mpd = nep->ncv;
194: if (!nep->which) nep->which = NEP_ALL;
195: PetscCheck(nep->which==NEP_ALL,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver supports only computing all eigenvalues");
196: NEPCheckUnsupported(nep,NEP_FEATURE_STOPPING | NEP_FEATURE_TWOSIDED);
198: /* check region */
199: PetscCall(RGIsTrivial(nep->rg,&istrivial));
200: PetscCheck(!istrivial,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
201: PetscCall(RGGetComplement(nep->rg,&flg));
202: PetscCheck(!flg,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
203: PetscCall(PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse));
204: PetscCheck(isellipse,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Currently only implemented for elliptic regions");
206: /* if the region has changed, then reset contour data */
207: PetscCall(PetscObjectGetId((PetscObject)nep->rg,&id));
208: PetscCall(PetscObjectStateGet((PetscObject)nep->rg,&state));
209: if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
210: PetscCall(SlepcContourDataDestroy(&ctx->contour));
211: PetscCall(PetscInfo(nep,"Resetting the contour data structure due to a change of region\n"));
212: ctx->rgid = id; ctx->rgstate = state;
213: }
215: /* create contour data structure */
216: if (!ctx->contour) {
217: PetscCall(RGCanUseConjugates(nep->rg,ctx->isreal,&ctx->useconj));
218: PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)nep,&ctx->contour));
219: }
221: PetscCall(NEPAllocateSolution(nep,0));
222: if (ctx->weight) PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
223: PetscCall(PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma));
225: /* allocate basis vectors */
226: PetscCall(BVDestroy(&ctx->S));
227: PetscCall(BVDuplicateResize(nep->V,ctx->L*ctx->M,&ctx->S));
228: PetscCall(BVDestroy(&ctx->V));
229: PetscCall(BVDuplicateResize(nep->V,ctx->L,&ctx->V));
231: contour = ctx->contour;
232: if (contour->subcomm && contour->subcomm->n != 1 && nep->fui==NEP_USER_INTERFACE_CALLBACK) {
233: PetscCall(NEPComputeFunction(nep,0,nep->function,nep->function_pre));
234: PetscCall(SlepcContourRedundantMat(contour,1,&nep->function,(nep->function!=nep->function_pre)?&nep->function_pre:NULL));
235: } else PetscCall(SlepcContourRedundantMat(contour,nep->nt,nep->A,nep->P));
236: if (contour->pA) {
237: if (!ctx->J) PetscCall(MatDuplicate(contour->pA[0],MAT_DO_NOT_COPY_VALUES,&ctx->J));
238: PetscCall(BVGetColumn(ctx->V,0,&v0));
239: PetscCall(SlepcContourScatterCreate(contour,v0));
240: PetscCall(BVRestoreColumn(ctx->V,0,&v0));
241: PetscCall(BVDestroy(&ctx->pV));
242: PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV));
243: PetscCall(BVSetSizesFromVec(ctx->pV,contour->xsub,nep->n));
244: PetscCall(BVSetFromOptions(ctx->pV));
245: PetscCall(BVResize(ctx->pV,ctx->L,PETSC_FALSE));
246: }
248: PetscCall(BVDestroy(&ctx->Y));
249: if (contour->pA) {
250: PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y));
251: PetscCall(BVSetSizesFromVec(ctx->Y,contour->xsub,nep->n));
252: PetscCall(BVSetFromOptions(ctx->Y));
253: PetscCall(BVResize(ctx->Y,contour->npoints*ctx->L,PETSC_FALSE));
254: } else PetscCall(BVDuplicateResize(nep->V,contour->npoints*ctx->L,&ctx->Y));
256: if (ctx->extraction == NEP_CISS_EXTRACTION_RITZ) {
257: PetscCall(DSSetMethod(nep->ds,1));
258: PetscCall(DSNEPSetRG(nep->ds,nep->rg));
259: if (nep->fui==NEP_USER_INTERFACE_SPLIT) PetscCall(DSNEPSetFN(nep->ds,nep->nt,nep->f));
260: else {
261: PetscCall(PetscNew(&dsctxf));
262: dsctxf->nep = nep;
263: PetscCall(DSNEPSetComputeMatrixFunction(nep->ds,NEPContourDSComputeMatrix,dsctxf,PetscCtxDestroyDefault));
264: }
265: }
266: PetscCall(DSAllocate(nep->ds,nep->ncv));
267: nwork = (nep->fui==NEP_USER_INTERFACE_SPLIT)? 2: 1;
268: PetscCall(NEPSetWorkVecs(nep,nwork));
269: PetscFunctionReturn(PETSC_SUCCESS);
270: }
272: static PetscErrorCode NEPSolve_CISS(NEP nep)
273: {
274: NEP_CISS *ctx = (NEP_CISS*)nep->data;
275: SlepcContourData contour = ctx->contour;
276: Mat X,M,E,T,P,J;
277: BV V;
278: PetscInt i,j,ld,L_add=0,nv=0,L_base=ctx->L,inner,*inside;
279: PetscScalar *Mu,*H0,*H1,*rr,*temp,center;
280: PetscReal error,max_error,radius,rgscale,est_eig,eta;
281: PetscBool isellipse,*fl1;
282: Vec si;
283: SlepcSC sc;
284: PetscRandom rand;
285: NEP_CISS_PROJECT dsctxf;
287: PetscFunctionBegin;
288: PetscCall(DSGetSlepcSC(nep->ds,&sc));
289: sc->comparison = SlepcCompareLargestMagnitude;
290: sc->comparisonctx = NULL;
291: sc->map = NULL;
292: sc->mapobj = NULL;
293: PetscCall(DSGetLeadingDimension(nep->ds,&ld));
294: PetscCall(RGComputeQuadrature(nep->rg,RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight));
295: if (contour->pA) {
296: T = contour->pA[0];
297: P = ((nep->fui==NEP_USER_INTERFACE_SPLIT && nep->P) || (nep->fui==NEP_USER_INTERFACE_CALLBACK && contour->pP))? contour->pP[0]: T;
298: } else {
299: T = nep->function;
300: P = nep->function_pre? nep->function_pre: nep->function;
301: }
302: PetscCall(NEPCISSSetUp(nep,T,P));
303: PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
304: PetscCall(BVSetRandomSign(ctx->V));
305: PetscCall(BVGetRandomContext(ctx->V,&rand));
306: if (contour->pA) {
307: J = ctx->J;
308: V = ctx->pV;
309: PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
310: } else {
311: J = nep->jacobian;
312: V = ctx->V;
313: }
314: PetscCall(NEPCISSSolve(nep,J,V,0,ctx->L));
315: PetscCall(PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse));
316: if (isellipse) {
317: PetscCall(BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig));
318: PetscCall(PetscInfo(nep,"Estimated eigenvalue count: %f\n",(double)est_eig));
319: eta = PetscPowReal(10.0,-PetscLog10Real(nep->tol)/ctx->N);
320: L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
321: if (L_add>ctx->L_max-ctx->L) {
322: PetscCall(PetscInfo(nep,"Number of eigenvalues inside the contour path may be too large\n"));
323: L_add = ctx->L_max-ctx->L;
324: }
325: }
326: /* Updates L after estimate the number of eigenvalue */
327: if (L_add>0) {
328: PetscCall(PetscInfo(nep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by Estimate #Eig\n",ctx->L,ctx->L+L_add));
329: PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
330: PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
331: PetscCall(BVSetRandomSign(ctx->V));
332: if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
333: ctx->L += L_add;
334: PetscCall(NEPCISSSolve(nep,J,V,ctx->L-L_add,ctx->L));
335: }
337: PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
338: for (i=0;i<ctx->refine_blocksize;i++) {
339: PetscCall(BVDotQuadrature(ctx->Y,V,Mu,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj));
340: PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
341: PetscCall(PetscLogEventBegin(NEP_CISS_SVD,nep,0,0,0));
342: PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
343: PetscCall(PetscLogEventEnd(NEP_CISS_SVD,nep,0,0,0));
344: if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
345: L_add = L_base;
346: if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
347: PetscCall(PetscInfo(nep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by SVD(H0)\n",ctx->L,ctx->L+L_add));
348: PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
349: PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
350: PetscCall(BVSetRandomSign(ctx->V));
351: if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
352: ctx->L += L_add;
353: PetscCall(NEPCISSSolve(nep,J,V,ctx->L-L_add,ctx->L));
354: if (L_add) {
355: PetscCall(PetscFree2(Mu,H0));
356: PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
357: }
358: }
360: PetscCall(RGGetScale(nep->rg,&rgscale));
361: PetscCall(RGEllipseGetParameters(nep->rg,¢er,&radius,NULL));
363: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) PetscCall(PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1));
365: while (nep->reason == NEP_CONVERGED_ITERATING) {
366: nep->its++;
367: for (inner=0;inner<=ctx->refine_inner;inner++) {
368: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
369: PetscCall(BVDotQuadrature(ctx->Y,V,Mu,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj));
370: PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
371: PetscCall(PetscLogEventBegin(NEP_CISS_SVD,nep,0,0,0));
372: PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
373: PetscCall(PetscLogEventEnd(NEP_CISS_SVD,nep,0,0,0));
374: } else {
375: PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
376: /* compute SVD of S */
377: PetscCall(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));
378: }
379: PetscCall(PetscInfo(nep,"Estimated rank: nv = %" PetscInt_FMT "\n",nv));
380: if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
381: PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
382: PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
383: PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
384: PetscCall(BVCopy(ctx->S,ctx->V));
385: if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
386: PetscCall(NEPCISSSolve(nep,J,V,0,ctx->L));
387: } else break;
388: }
389: nep->nconv = 0;
390: if (nv == 0) { nep->reason = NEP_CONVERGED_TOL; break; }
391: else {
392: /* Extracting eigenpairs */
393: PetscCall(DSSetDimensions(nep->ds,nv,0,0));
394: PetscCall(DSSetState(nep->ds,DS_STATE_RAW));
395: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
396: PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
397: PetscCall(CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1));
398: PetscCall(DSGetArray(nep->ds,DS_MAT_A,&temp));
399: for (j=0;j<nv;j++)
400: for (i=0;i<nv;i++)
401: temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
402: PetscCall(DSRestoreArray(nep->ds,DS_MAT_A,&temp));
403: PetscCall(DSGetArray(nep->ds,DS_MAT_B,&temp));
404: for (j=0;j<nv;j++)
405: for (i=0;i<nv;i++)
406: temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
407: PetscCall(DSRestoreArray(nep->ds,DS_MAT_B,&temp));
408: } else if (ctx->extraction == NEP_CISS_EXTRACTION_CAA) {
409: PetscCall(BVSetActiveColumns(ctx->S,0,nv));
410: PetscCall(DSGetArray(nep->ds,DS_MAT_A,&temp));
411: for (i=0;i<nv;i++) PetscCall(PetscArraycpy(temp+i*ld,H0+i*nv,nv));
412: PetscCall(DSRestoreArray(nep->ds,DS_MAT_A,&temp));
413: } else {
414: PetscCall(BVSetActiveColumns(ctx->S,0,nv));
415: if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
416: for (i=0;i<nep->nt;i++) {
417: PetscCall(DSGetMat(nep->ds,DSMatExtra[i],&E));
418: PetscCall(BVMatProject(ctx->S,nep->A[i],ctx->S,E));
419: PetscCall(DSRestoreMat(nep->ds,DSMatExtra[i],&E));
420: }
421: } else {
422: PetscCall(DSNEPGetComputeMatrixFunction(nep->ds,NULL,(PetscCtxRt)&dsctxf,NULL));
423: dsctxf->Q = ctx->S;
424: }
425: }
426: PetscCall(DSSolve(nep->ds,nep->eigr,nep->eigi));
427: PetscCall(DSSynchronize(nep->ds,nep->eigr,nep->eigi));
428: PetscCall(DSGetDimensions(nep->ds,NULL,NULL,NULL,&nv));
429: if (ctx->extraction == NEP_CISS_EXTRACTION_CAA || ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
430: for (i=0;i<nv;i++) {
431: nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
432: }
433: }
434: PetscCall(PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr));
435: PetscCall(DSVectors(nep->ds,DS_MAT_X,NULL,NULL));
436: PetscCall(DSGetMat(nep->ds,DS_MAT_X,&X));
437: PetscCall(SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1));
438: PetscCall(DSRestoreMat(nep->ds,DS_MAT_X,&X));
439: PetscCall(RGCheckInside(nep->rg,nv,nep->eigr,nep->eigi,inside));
440: for (i=0;i<nv;i++) {
441: if (fl1[i] && inside[i]>=0) {
442: rr[i] = 1.0;
443: nep->nconv++;
444: } else rr[i] = 0.0;
445: }
446: PetscCall(DSSort(nep->ds,nep->eigr,nep->eigi,rr,NULL,&nep->nconv));
447: PetscCall(DSSynchronize(nep->ds,nep->eigr,nep->eigi));
448: if (ctx->extraction == NEP_CISS_EXTRACTION_CAA || ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
449: for (i=0;i<nv;i++) nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
450: }
451: PetscCall(PetscFree3(fl1,inside,rr));
452: PetscCall(BVSetActiveColumns(nep->V,0,nv));
453: PetscCall(DSVectors(nep->ds,DS_MAT_X,NULL,NULL));
454: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
455: PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
456: PetscCall(BVSetActiveColumns(ctx->S,0,nv));
457: PetscCall(BVCopy(ctx->S,nep->V));
458: PetscCall(DSGetMat(nep->ds,DS_MAT_X,&X));
459: PetscCall(BVMultInPlace(ctx->S,X,0,nep->nconv));
460: PetscCall(BVMultInPlace(nep->V,X,0,nep->nconv));
461: PetscCall(DSRestoreMat(nep->ds,DS_MAT_X,&X));
462: } else {
463: PetscCall(DSGetMat(nep->ds,DS_MAT_X,&X));
464: PetscCall(BVMultInPlace(ctx->S,X,0,nep->nconv));
465: PetscCall(DSRestoreMat(nep->ds,DS_MAT_X,&X));
466: PetscCall(BVSetActiveColumns(ctx->S,0,nep->nconv));
467: PetscCall(BVCopy(ctx->S,nep->V));
468: }
469: max_error = 0.0;
470: for (i=0;i<nep->nconv;i++) {
471: PetscCall(BVGetColumn(nep->V,i,&si));
472: PetscCall(VecNormalize(si,NULL));
473: PetscCall(NEPComputeResidualNorm_Private(nep,PETSC_FALSE,nep->eigr[i],si,nep->work,&error));
474: PetscCall((*nep->converged)(nep,nep->eigr[i],0,error,&error,nep->convergedctx));
475: PetscCall(BVRestoreColumn(nep->V,i,&si));
476: max_error = PetscMax(max_error,error);
477: }
478: if (max_error <= nep->tol) nep->reason = NEP_CONVERGED_TOL;
479: else if (nep->its > nep->max_it) nep->reason = NEP_DIVERGED_ITS;
480: else {
481: if (nep->nconv > ctx->L) nv = nep->nconv;
482: else if (ctx->L > nv) nv = ctx->L;
483: nv = PetscMin(nv,ctx->L*ctx->M);
484: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M));
485: PetscCall(MatSetRandom(M,rand));
486: PetscCall(BVSetActiveColumns(ctx->S,0,nv));
487: PetscCall(BVMultInPlace(ctx->S,M,0,ctx->L));
488: PetscCall(MatDestroy(&M));
489: PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
490: PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
491: PetscCall(BVCopy(ctx->S,ctx->V));
492: if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
493: PetscCall(NEPCISSSolve(nep,J,V,0,ctx->L));
494: }
495: }
496: }
497: PetscCall(PetscFree2(Mu,H0));
498: if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) PetscCall(PetscFree(H1));
499: PetscFunctionReturn(PETSC_SUCCESS);
500: }
502: static PetscErrorCode NEPCISSSetSizes_CISS(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
503: {
504: NEP_CISS *ctx = (NEP_CISS*)nep->data;
505: PetscInt oN,oL,oM,oLmax,onpart;
506: PetscMPIInt size;
508: PetscFunctionBegin;
509: oN = ctx->N;
510: if (ip == PETSC_DETERMINE) {
511: if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
512: } else if (ip != PETSC_CURRENT) {
513: PetscCheck(ip>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
514: PetscCheck(ip%2==0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
515: if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
516: }
517: oL = ctx->L;
518: if (bs == PETSC_DETERMINE) {
519: ctx->L = 16;
520: } else if (bs != PETSC_CURRENT) {
521: PetscCheck(bs>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
522: ctx->L = bs;
523: }
524: oM = ctx->M;
525: if (ms == PETSC_DETERMINE) {
526: ctx->M = ctx->N/4;
527: } else if (ms != PETSC_CURRENT) {
528: PetscCheck(ms>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
529: PetscCheck(ms<=ctx->N,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be less than or equal to the number of integration points");
530: ctx->M = ms;
531: }
532: onpart = ctx->npart;
533: if (npart == PETSC_DETERMINE) {
534: ctx->npart = 1;
535: } else if (npart != PETSC_CURRENT) {
536: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)nep),&size));
537: PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
538: ctx->npart = npart;
539: }
540: oLmax = ctx->L_max;
541: if (bsmax == PETSC_DETERMINE) {
542: ctx->L_max = 64;
543: } else if (bsmax != PETSC_CURRENT) {
544: PetscCheck(bsmax>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
545: ctx->L_max = PetscMax(bsmax,ctx->L);
546: }
547: if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
548: PetscCall(SlepcContourDataDestroy(&ctx->contour));
549: PetscCall(PetscInfo(nep,"Resetting the contour data structure due to a change of parameters\n"));
550: nep->state = NEP_STATE_INITIAL;
551: }
552: ctx->isreal = realmats;
553: if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) nep->state = NEP_STATE_INITIAL;
554: PetscFunctionReturn(PETSC_SUCCESS);
555: }
557: /*@
558: NEPCISSSetSizes - Sets the values of various size parameters in the CISS solver.
560: Logically Collective
562: Input Parameters:
563: + nep - the nonlinear eigensolver context
564: . ip - number of integration points
565: . bs - block size
566: . ms - moment size
567: . npart - number of partitions when splitting the communicator
568: . bsmax - max block size
569: - realmats - $T(z)$ is real for real $z$
571: Options Database Keys:
572: + -nep_ciss_integration_points ip - sets the number of integration points
573: . -nep_ciss_blocksize bs - sets the block size
574: . -nep_ciss_moments ms - sets the moment size
575: . -nep_ciss_partitions npart - sets the number of partitions
576: . -nep_ciss_maxblocksize bsmax - sets the maximum block size
577: - -nep_ciss_realmats (true|false) - $T(z)$ is real for real $z$
579: Notes:
580: For all integer arguments, you can use `PETSC_CURRENT` to keep the current value, and
581: `PETSC_DETERMINE` to set them to a default value.
583: The default number of partitions is 1. This means the internal `KSP` object is shared
584: among all processes of the `NEP` communicator. Otherwise, the communicator is split
585: into `npart` communicators, so that `npart` `KSP` solves proceed simultaneously.
587: The `realmats` flag can be set to `PETSC_TRUE` when $T(\cdot)$ is guaranteed to be real
588: when the argument is a real value, for example, when all matrices in
589: the split form are real. When set to `PETSC_TRUE`, the solver avoids some computations.
591: For a detailed description of the parameters see {cite:p}`Mae16`.
593: Level: advanced
595: .seealso: [](ch:nep), `NEPCISS`, `NEPCISSGetSizes()`
596: @*/
597: PetscErrorCode NEPCISSSetSizes(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
598: {
599: PetscFunctionBegin;
607: PetscTryMethod(nep,"NEPCISSSetSizes_C",(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(nep,ip,bs,ms,npart,bsmax,realmats));
608: PetscFunctionReturn(PETSC_SUCCESS);
609: }
611: static PetscErrorCode NEPCISSGetSizes_CISS(NEP nep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
612: {
613: NEP_CISS *ctx = (NEP_CISS*)nep->data;
615: PetscFunctionBegin;
616: if (ip) *ip = ctx->N;
617: if (bs) *bs = ctx->L;
618: if (ms) *ms = ctx->M;
619: if (npart) *npart = ctx->npart;
620: if (bsmax) *bsmax = ctx->L_max;
621: if (realmats) *realmats = ctx->isreal;
622: PetscFunctionReturn(PETSC_SUCCESS);
623: }
625: /*@
626: NEPCISSGetSizes - Gets the values of various size parameters in the CISS solver.
628: Not Collective
630: Input Parameter:
631: . nep - the nonlinear eigensolver context
633: Output Parameters:
634: + ip - number of integration points
635: . bs - block size
636: . ms - moment size
637: . npart - number of partitions when splitting the communicator
638: . bsmax - max block size
639: - realmats - $T(z)$ is real for real $z$
641: Level: advanced
643: .seealso: [](ch:nep), `NEPCISS`, `NEPCISSSetSizes()`
644: @*/
645: PetscErrorCode NEPCISSGetSizes(NEP nep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
646: {
647: PetscFunctionBegin;
649: PetscUseMethod(nep,"NEPCISSGetSizes_C",(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(nep,ip,bs,ms,npart,bsmax,realmats));
650: PetscFunctionReturn(PETSC_SUCCESS);
651: }
653: static PetscErrorCode NEPCISSSetThreshold_CISS(NEP nep,PetscReal delta,PetscReal spur)
654: {
655: NEP_CISS *ctx = (NEP_CISS*)nep->data;
657: PetscFunctionBegin;
658: if (delta == (PetscReal)PETSC_DETERMINE) {
659: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
660: } else if (delta != (PetscReal)PETSC_CURRENT) {
661: PetscCheck(delta>0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
662: ctx->delta = delta;
663: }
664: if (spur == (PetscReal)PETSC_DETERMINE) {
665: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
666: } else if (spur != (PetscReal)PETSC_CURRENT) {
667: PetscCheck(spur>0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
668: ctx->spurious_threshold = spur;
669: }
670: PetscFunctionReturn(PETSC_SUCCESS);
671: }
673: /*@
674: NEPCISSSetThreshold - Sets the values of various threshold parameters in
675: the CISS solver.
677: Logically Collective
679: Input Parameters:
680: + nep - the nonlinear eigensolver context
681: . delta - threshold for numerical rank
682: - spur - spurious threshold (to discard spurious eigenpairs)
684: Options Database Keys:
685: + -nep_ciss_delta delta - sets the delta
686: - -nep_ciss_spurious_threshold spur - sets the spurious threshold
688: Notes:
689: `PETSC_CURRENT` can be used to preserve the current value of any of the
690: arguments, and `PETSC_DETERMINE` to set them to a default value.
692: For a detailed description of the parameters see {cite:p}`Mae16`.
694: Level: advanced
696: .seealso: [](ch:nep), `NEPCISS`, `NEPCISSGetThreshold()`
697: @*/
698: PetscErrorCode NEPCISSSetThreshold(NEP nep,PetscReal delta,PetscReal spur)
699: {
700: PetscFunctionBegin;
704: PetscTryMethod(nep,"NEPCISSSetThreshold_C",(NEP,PetscReal,PetscReal),(nep,delta,spur));
705: PetscFunctionReturn(PETSC_SUCCESS);
706: }
708: static PetscErrorCode NEPCISSGetThreshold_CISS(NEP nep,PetscReal *delta,PetscReal *spur)
709: {
710: NEP_CISS *ctx = (NEP_CISS*)nep->data;
712: PetscFunctionBegin;
713: if (delta) *delta = ctx->delta;
714: if (spur) *spur = ctx->spurious_threshold;
715: PetscFunctionReturn(PETSC_SUCCESS);
716: }
718: /*@
719: NEPCISSGetThreshold - Gets the values of various threshold parameters in
720: the CISS solver.
722: Not Collective
724: Input Parameter:
725: . nep - the nonlinear eigensolver context
727: Output Parameters:
728: + delta - threshold for numerical rank
729: - spur - spurious threshold (to discard spurious eigenpairs)
731: Level: advanced
733: .seealso: [](ch:nep), `NEPCISS`, `NEPCISSSetThreshold()`
734: @*/
735: PetscErrorCode NEPCISSGetThreshold(NEP nep,PetscReal *delta,PetscReal *spur)
736: {
737: PetscFunctionBegin;
739: PetscUseMethod(nep,"NEPCISSGetThreshold_C",(NEP,PetscReal*,PetscReal*),(nep,delta,spur));
740: PetscFunctionReturn(PETSC_SUCCESS);
741: }
743: static PetscErrorCode NEPCISSSetRefinement_CISS(NEP nep,PetscInt inner,PetscInt blsize)
744: {
745: NEP_CISS *ctx = (NEP_CISS*)nep->data;
747: PetscFunctionBegin;
748: if (inner == PETSC_DETERMINE) {
749: ctx->refine_inner = 0;
750: } else if (inner != PETSC_CURRENT) {
751: PetscCheck(inner>=0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
752: ctx->refine_inner = inner;
753: }
754: if (blsize == PETSC_DETERMINE) {
755: ctx->refine_blocksize = 0;
756: } else if (blsize != PETSC_CURRENT) {
757: PetscCheck(blsize>=0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
758: ctx->refine_blocksize = blsize;
759: }
760: PetscFunctionReturn(PETSC_SUCCESS);
761: }
763: /*@
764: NEPCISSSetRefinement - Sets the values of various refinement parameters
765: in the CISS solver.
767: Logically Collective
769: Input Parameters:
770: + nep - the nonlinear eigensolver context
771: . inner - number of iterative refinement iterations (inner loop)
772: - blsize - number of iterative refinement iterations (blocksize loop)
774: Options Database Keys:
775: + -nep_ciss_refine_inner inner - sets number of inner iterations
776: - -nep_ciss_refine_blocksize blsize - sets number of blocksize iterations
778: Note:
779: `PETSC_CURRENT` can be used to preserve the current value of any of the
780: arguments, and `PETSC_DETERMINE` to set them to a default of 0 (no refinement).
782: Level: advanced
784: .seealso: [](ch:nep), `NEPCISS`, `NEPCISSGetRefinement()`
785: @*/
786: PetscErrorCode NEPCISSSetRefinement(NEP nep,PetscInt inner,PetscInt blsize)
787: {
788: PetscFunctionBegin;
792: PetscTryMethod(nep,"NEPCISSSetRefinement_C",(NEP,PetscInt,PetscInt),(nep,inner,blsize));
793: PetscFunctionReturn(PETSC_SUCCESS);
794: }
796: static PetscErrorCode NEPCISSGetRefinement_CISS(NEP nep,PetscInt *inner,PetscInt *blsize)
797: {
798: NEP_CISS *ctx = (NEP_CISS*)nep->data;
800: PetscFunctionBegin;
801: if (inner) *inner = ctx->refine_inner;
802: if (blsize) *blsize = ctx->refine_blocksize;
803: PetscFunctionReturn(PETSC_SUCCESS);
804: }
806: /*@
807: NEPCISSGetRefinement - Gets the values of various refinement parameters
808: in the CISS solver.
810: Not Collective
812: Input Parameter:
813: . nep - the nonlinear eigensolver context
815: Output Parameters:
816: + inner - number of iterative refinement iterations (inner loop)
817: - blsize - number of iterative refinement iterations (blocksize loop)
819: Level: advanced
821: .seealso: [](ch:nep), `NEPCISS`, `NEPCISSSetRefinement()`
822: @*/
823: PetscErrorCode NEPCISSGetRefinement(NEP nep, PetscInt *inner, PetscInt *blsize)
824: {
825: PetscFunctionBegin;
827: PetscUseMethod(nep,"NEPCISSGetRefinement_C",(NEP,PetscInt*,PetscInt*),(nep,inner,blsize));
828: PetscFunctionReturn(PETSC_SUCCESS);
829: }
831: static PetscErrorCode NEPCISSSetExtraction_CISS(NEP nep,NEPCISSExtraction extraction)
832: {
833: NEP_CISS *ctx = (NEP_CISS*)nep->data;
835: PetscFunctionBegin;
836: if (ctx->extraction != extraction) {
837: ctx->extraction = extraction;
838: nep->state = NEP_STATE_INITIAL;
839: }
840: PetscFunctionReturn(PETSC_SUCCESS);
841: }
843: /*@
844: NEPCISSSetExtraction - Sets the extraction technique used in the CISS solver.
846: Logically Collective
848: Input Parameters:
849: + nep - the nonlinear eigensolver context
850: - extraction - the extraction technique, see `NEPCISSExtraction` for possible values
852: Options Database Key:
853: . -nep_ciss_extraction (ritz|hankel|caa) - sets the extraction technique
855: Notes:
856: By default, the Rayleigh-Ritz extraction is used (`NEP_CISS_EXTRACTION_RITZ`).
858: If the `hankel` or the `caa` option is specified (`NEP_CISS_EXTRACTION_HANKEL` or
859: `NEP_CISS_EXTRACTION_CAA`), then the block Hankel method, or the communication-avoiding
860: Arnoldi method, respectively, is used for extracting eigenpairs.
862: Level: advanced
864: .seealso: [](ch:nep), `NEPCISS`, `NEPCISSGetExtraction()`, `NEPCISSExtraction`
865: @*/
866: PetscErrorCode NEPCISSSetExtraction(NEP nep,NEPCISSExtraction extraction)
867: {
868: PetscFunctionBegin;
871: PetscTryMethod(nep,"NEPCISSSetExtraction_C",(NEP,NEPCISSExtraction),(nep,extraction));
872: PetscFunctionReturn(PETSC_SUCCESS);
873: }
875: static PetscErrorCode NEPCISSGetExtraction_CISS(NEP nep,NEPCISSExtraction *extraction)
876: {
877: NEP_CISS *ctx = (NEP_CISS*)nep->data;
879: PetscFunctionBegin;
880: *extraction = ctx->extraction;
881: PetscFunctionReturn(PETSC_SUCCESS);
882: }
884: /*@
885: NEPCISSGetExtraction - Gets the extraction technique used in the CISS solver.
887: Not Collective
889: Input Parameter:
890: . nep - the nonlinear eigensolver context
892: Output Parameter:
893: . extraction - extraction technique
895: Level: advanced
897: .seealso: [](ch:nep), `NEPCISS`, `NEPCISSSetExtraction()`, `NEPCISSExtraction`
898: @*/
899: PetscErrorCode NEPCISSGetExtraction(NEP nep,NEPCISSExtraction *extraction)
900: {
901: PetscFunctionBegin;
903: PetscAssertPointer(extraction,2);
904: PetscUseMethod(nep,"NEPCISSGetExtraction_C",(NEP,NEPCISSExtraction*),(nep,extraction));
905: PetscFunctionReturn(PETSC_SUCCESS);
906: }
908: static PetscErrorCode NEPCISSGetKSPs_CISS(NEP nep,PetscInt *nsolve,KSP **ksp)
909: {
910: NEP_CISS *ctx = (NEP_CISS*)nep->data;
911: SlepcContourData contour;
912: PetscInt i;
913: PC pc;
914: MPI_Comm child;
916: PetscFunctionBegin;
917: if (!ctx->contour) { /* initialize contour data structure first */
918: PetscCall(RGCanUseConjugates(nep->rg,ctx->isreal,&ctx->useconj));
919: PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)nep,&ctx->contour));
920: }
921: contour = ctx->contour;
922: if (!contour->ksp) {
923: PetscCall(PetscMalloc1(contour->npoints,&contour->ksp));
924: PetscCall(PetscSubcommGetChild(contour->subcomm,&child));
925: for (i=0;i<contour->npoints;i++) {
926: PetscCall(KSPCreate(child,&contour->ksp[i]));
927: PetscCall(PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)nep,1));
928: PetscCall(KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)nep)->prefix));
929: PetscCall(KSPAppendOptionsPrefix(contour->ksp[i],"nep_ciss_"));
930: PetscCall(PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)nep)->options));
931: PetscCall(KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE));
932: PetscCall(KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(nep->tol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
933: PetscCall(KSPGetPC(contour->ksp[i],&pc));
934: if ((nep->fui==NEP_USER_INTERFACE_SPLIT && nep->P) || (nep->fui==NEP_USER_INTERFACE_CALLBACK && nep->function_pre!=nep->function)) {
935: PetscCall(KSPSetType(contour->ksp[i],KSPBCGS));
936: PetscCall(PCSetType(pc,PCBJACOBI));
937: } else {
938: PetscCall(KSPSetType(contour->ksp[i],KSPPREONLY));
939: PetscCall(PCSetType(pc,PCLU));
940: }
941: }
942: }
943: if (nsolve) *nsolve = contour->npoints;
944: if (ksp) *ksp = contour->ksp;
945: PetscFunctionReturn(PETSC_SUCCESS);
946: }
948: /*@C
949: NEPCISSGetKSPs - Retrieve the array of linear solver objects associated with
950: the CISS solver.
952: Collective
954: Input Parameter:
955: . nep - the nonlinear eigensolver context
957: Output Parameters:
958: + nsolve - number of solver objects
959: - ksp - array of linear solver object
961: Notes:
962: The number of `KSP` solvers is equal to the number of integration points divided by
963: the number of partitions. This value is halved in the case of real matrices with
964: a region centered at the real axis.
966: Level: advanced
968: .seealso: [](ch:nep), `NEPCISS`, `NEPCISSSetSizes()`
969: @*/
970: PetscErrorCode NEPCISSGetKSPs(NEP nep,PetscInt *nsolve,KSP **ksp)
971: {
972: PetscFunctionBegin;
974: PetscUseMethod(nep,"NEPCISSGetKSPs_C",(NEP,PetscInt*,KSP**),(nep,nsolve,ksp));
975: PetscFunctionReturn(PETSC_SUCCESS);
976: }
978: static PetscErrorCode NEPReset_CISS(NEP nep)
979: {
980: NEP_CISS *ctx = (NEP_CISS*)nep->data;
982: PetscFunctionBegin;
983: PetscCall(BVDestroy(&ctx->S));
984: PetscCall(BVDestroy(&ctx->V));
985: PetscCall(BVDestroy(&ctx->Y));
986: PetscCall(SlepcContourDataReset(ctx->contour));
987: PetscCall(MatDestroy(&ctx->J));
988: PetscCall(BVDestroy(&ctx->pV));
989: PetscFunctionReturn(PETSC_SUCCESS);
990: }
992: static PetscErrorCode NEPSetFromOptions_CISS(NEP nep,PetscOptionItems PetscOptionsObject)
993: {
994: NEP_CISS *ctx = (NEP_CISS*)nep->data;
995: PetscReal r1,r2;
996: PetscInt i,i1,i2,i3,i4,i5,i6,i7;
997: PetscBool b1,flg,flg2,flg3,flg4,flg5,flg6;
998: NEPCISSExtraction extraction;
1000: PetscFunctionBegin;
1001: PetscOptionsHeadBegin(PetscOptionsObject,"NEP CISS Options");
1003: PetscCall(NEPCISSGetSizes(nep,&i1,&i2,&i3,&i4,&i5,&b1));
1004: PetscCall(PetscOptionsInt("-nep_ciss_integration_points","Number of integration points","NEPCISSSetSizes",i1,&i1,&flg));
1005: PetscCall(PetscOptionsInt("-nep_ciss_blocksize","Block size","NEPCISSSetSizes",i2,&i2,&flg2));
1006: PetscCall(PetscOptionsInt("-nep_ciss_moments","Moment size","NEPCISSSetSizes",i3,&i3,&flg3));
1007: PetscCall(PetscOptionsInt("-nep_ciss_partitions","Number of partitions","NEPCISSSetSizes",i4,&i4,&flg4));
1008: PetscCall(PetscOptionsInt("-nep_ciss_maxblocksize","Maximum block size","NEPCISSSetSizes",i5,&i5,&flg5));
1009: PetscCall(PetscOptionsBool("-nep_ciss_realmats","True if T(z) is real for real z","NEPCISSSetSizes",b1,&b1,&flg6));
1010: if (flg || flg2 || flg3 || flg4 || flg5 || flg6) PetscCall(NEPCISSSetSizes(nep,i1,i2,i3,i4,i5,b1));
1012: PetscCall(NEPCISSGetThreshold(nep,&r1,&r2));
1013: PetscCall(PetscOptionsReal("-nep_ciss_delta","Threshold for numerical rank","NEPCISSSetThreshold",r1,&r1,&flg));
1014: PetscCall(PetscOptionsReal("-nep_ciss_spurious_threshold","Threshold for the spurious eigenpairs","NEPCISSSetThreshold",r2,&r2,&flg2));
1015: if (flg || flg2) PetscCall(NEPCISSSetThreshold(nep,r1,r2));
1017: PetscCall(NEPCISSGetRefinement(nep,&i6,&i7));
1018: PetscCall(PetscOptionsInt("-nep_ciss_refine_inner","Number of inner iterative refinement iterations","NEPCISSSetRefinement",i6,&i6,&flg));
1019: PetscCall(PetscOptionsInt("-nep_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","NEPCISSSetRefinement",i7,&i7,&flg2));
1020: if (flg || flg2) PetscCall(NEPCISSSetRefinement(nep,i6,i7));
1022: PetscCall(PetscOptionsEnum("-nep_ciss_extraction","Extraction technique","NEPCISSSetExtraction",NEPCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg));
1023: if (flg) PetscCall(NEPCISSSetExtraction(nep,extraction));
1025: PetscOptionsHeadEnd();
1027: if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg));
1028: PetscCall(RGSetFromOptions(nep->rg)); /* this is necessary here to set useconj */
1029: if (!ctx->contour || !ctx->contour->ksp) PetscCall(NEPCISSGetKSPs(nep,NULL,NULL));
1030: PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)nep),PETSC_ERR_PLIB,"Something went wrong with NEPCISSGetKSPs()");
1031: for (i=0;i<ctx->contour->npoints;i++) PetscCall(KSPSetFromOptions(ctx->contour->ksp[i]));
1032: PetscCall(PetscSubcommSetFromOptions(ctx->contour->subcomm));
1033: PetscFunctionReturn(PETSC_SUCCESS);
1034: }
1036: static PetscErrorCode NEPDestroy_CISS(NEP nep)
1037: {
1038: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1040: PetscFunctionBegin;
1041: PetscCall(SlepcContourDataDestroy(&ctx->contour));
1042: PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
1043: PetscCall(PetscFree(nep->data));
1044: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NULL));
1045: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NULL));
1046: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NULL));
1047: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NULL));
1048: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NULL));
1049: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NULL));
1050: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetExtraction_C",NULL));
1051: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetExtraction_C",NULL));
1052: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetKSPs_C",NULL));
1053: PetscFunctionReturn(PETSC_SUCCESS);
1054: }
1056: static PetscErrorCode NEPView_CISS(NEP nep,PetscViewer viewer)
1057: {
1058: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1059: PetscBool isascii;
1060: PetscViewer sviewer;
1062: PetscFunctionBegin;
1063: PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
1064: if (isascii) {
1065: PetscCall(PetscViewerASCIIPrintf(viewer," sizes { integration points: %" PetscInt_FMT ", block size: %" PetscInt_FMT ", moment size: %" PetscInt_FMT ", partitions: %" PetscInt_FMT ", maximum block size: %" PetscInt_FMT " }\n",ctx->N,ctx->L,ctx->M,ctx->npart,ctx->L_max));
1066: if (ctx->isreal) PetscCall(PetscViewerASCIIPrintf(viewer," exploiting symmetry of integration points\n"));
1067: PetscCall(PetscViewerASCIIPrintf(viewer," threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold));
1068: PetscCall(PetscViewerASCIIPrintf(viewer," iterative refinement { inner: %" PetscInt_FMT ", blocksize: %" PetscInt_FMT " }\n",ctx->refine_inner, ctx->refine_blocksize));
1069: PetscCall(PetscViewerASCIIPrintf(viewer," extraction: %s\n",NEPCISSExtractions[ctx->extraction]));
1070: if (!ctx->contour || !ctx->contour->ksp) PetscCall(NEPCISSGetKSPs(nep,NULL,NULL));
1071: PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)nep),PETSC_ERR_PLIB,"Something went wrong with NEPCISSGetKSPs()");
1072: PetscCall(PetscViewerASCIIPushTab(viewer));
1073: if (ctx->npart>1 && ctx->contour->subcomm) {
1074: PetscCall(PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
1075: if (!ctx->contour->subcomm->color) PetscCall(KSPView(ctx->contour->ksp[0],sviewer));
1076: PetscCall(PetscViewerFlush(sviewer));
1077: PetscCall(PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
1078: /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1079: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1080: } else PetscCall(KSPView(ctx->contour->ksp[0],viewer));
1081: PetscCall(PetscViewerASCIIPopTab(viewer));
1082: }
1083: PetscFunctionReturn(PETSC_SUCCESS);
1084: }
1086: static PetscErrorCode NEPSetDSType_CISS(NEP nep)
1087: {
1088: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1090: PetscFunctionBegin;
1091: switch (ctx->extraction) {
1092: case NEP_CISS_EXTRACTION_RITZ: PetscCall(DSSetType(nep->ds,DSNEP)); break;
1093: case NEP_CISS_EXTRACTION_HANKEL: PetscCall(DSSetType(nep->ds,DSGNHEP)); break;
1094: case NEP_CISS_EXTRACTION_CAA: PetscCall(DSSetType(nep->ds,DSNHEP)); break;
1095: }
1096: PetscFunctionReturn(PETSC_SUCCESS);
1097: }
1099: /*MC
1100: NEPCISS - NEPCISS = "ciss" - A contour integral eigensolver based on the
1101: Sakurai-Sugiura scheme.
1103: Notes:
1104: This solver is based on the numerical contour integration idea
1105: proposed initially for linear problems by {cite:t}`Sak03`. In nonlinear
1106: eigenproblems, a Rayleigh-Ritz projection is done, resulting in
1107: a small dense nonlinear eigenproblem {cite:p}`Asa09,Yok13`.
1109: Contour integral methods are able to compute all eigenvalues
1110: lying inside a region of the complex plane. Use `NEPGetRG()` to
1111: specify the region. However, the computational cost is usually high
1112: because multiple linear systems must be solved. Use `NEPCISSGetKSPs()`
1113: to configure the `KSP` objects for this.
1115: Details of the implementation in SLEPc can be found in {cite:p}`Mae16`.
1117: Level: beginner
1119: .seealso: [](ch:nep), `NEP`, `NEPType`, `NEPSetType()`, `NEPGetRG()`
1120: M*/
1121: SLEPC_EXTERN PetscErrorCode NEPCreate_CISS(NEP nep)
1122: {
1123: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1125: PetscFunctionBegin;
1126: PetscCall(PetscNew(&ctx));
1127: nep->data = ctx;
1128: /* set default values of parameters */
1129: ctx->N = 32;
1130: ctx->L = 16;
1131: ctx->M = ctx->N/4;
1132: ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
1133: ctx->L_max = 64;
1134: ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1135: ctx->isreal = PETSC_FALSE;
1136: ctx->npart = 1;
1138: nep->useds = PETSC_TRUE;
1140: nep->ops->solve = NEPSolve_CISS;
1141: nep->ops->setup = NEPSetUp_CISS;
1142: nep->ops->setfromoptions = NEPSetFromOptions_CISS;
1143: nep->ops->reset = NEPReset_CISS;
1144: nep->ops->destroy = NEPDestroy_CISS;
1145: nep->ops->view = NEPView_CISS;
1146: nep->ops->setdstype = NEPSetDSType_CISS;
1148: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NEPCISSSetSizes_CISS));
1149: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NEPCISSGetSizes_CISS));
1150: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NEPCISSSetThreshold_CISS));
1151: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NEPCISSGetThreshold_CISS));
1152: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NEPCISSSetRefinement_CISS));
1153: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NEPCISSGetRefinement_CISS));
1154: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetExtraction_C",NEPCISSSetExtraction_CISS));
1155: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetExtraction_C",NEPCISSGetExtraction_CISS));
1156: PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetKSPs_C",NEPCISSGetKSPs_CISS));
1157: PetscFunctionReturn(PETSC_SUCCESS);
1158: }