Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
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"
12 :
13 : Method: Contour Integral Spectral Slicing
14 :
15 : Algorithm:
16 :
17 : Contour integral based on Sakurai-Sugiura method to construct a
18 : subspace, with various eigenpair extractions (Rayleigh-Ritz,
19 : explicit moment).
20 :
21 : Based on code contributed by Y. Maeda, T. Sakurai.
22 :
23 : References:
24 :
25 : [1] T. Sakurai and H. Sugiura, "A projection method for generalized
26 : eigenvalue problems", J. Comput. Appl. Math. 159:119-128, 2003.
27 :
28 : [2] T. Sakurai and H. Tadano, "CIRR: a Rayleigh-Ritz type method with
29 : contour integral for generalized eigenvalue problems", Hokkaido
30 : Math. J. 36:745-757, 2007.
31 : */
32 :
33 : #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
34 : #include <slepc/private/slepccontour.h>
35 : #include <slepcblaslapack.h>
36 :
37 : typedef struct {
38 : /* user parameters */
39 : PetscInt N; /* number of integration points (32) */
40 : PetscInt L; /* block size (16) */
41 : PetscInt M; /* moment degree (N/4 = 4) */
42 : PetscReal delta; /* threshold of singular value (1e-12) */
43 : PetscInt L_max; /* maximum number of columns of the source matrix V */
44 : PetscReal spurious_threshold; /* discard spurious eigenpairs */
45 : PetscBool isreal; /* A and B are real */
46 : PetscInt npart; /* number of partitions */
47 : PetscInt refine_inner;
48 : PetscInt refine_blocksize;
49 : EPSCISSQuadRule quad;
50 : EPSCISSExtraction extraction;
51 : PetscBool usest;
52 : /* private data */
53 : SlepcContourData contour;
54 : PetscReal *sigma; /* threshold for numerical rank */
55 : PetscScalar *weight;
56 : PetscScalar *omega;
57 : PetscScalar *pp;
58 : BV V;
59 : BV S;
60 : BV pV;
61 : BV Y;
62 : PetscBool useconj;
63 : PetscBool usest_set; /* whether the user set the usest flag or not */
64 : PetscObjectId rgid;
65 : PetscObjectState rgstate;
66 : } EPS_CISS;
67 :
68 : /*
69 : Set up KSP solvers for every integration point, only called if !ctx->usest
70 : */
71 17 : static PetscErrorCode EPSCISSSetUp(EPS eps,Mat A,Mat B,Mat Pa,Mat Pb)
72 : {
73 17 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
74 17 : SlepcContourData contour;
75 17 : PetscInt i,p_id,nsplit;
76 17 : Mat Amat,Pmat;
77 17 : MatStructure str,strp;
78 :
79 17 : PetscFunctionBegin;
80 17 : if (!ctx->contour || !ctx->contour->ksp) PetscCall(EPSCISSGetKSPs(eps,NULL,NULL));
81 17 : PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Something went wrong with EPSCISSGetKSPs()");
82 17 : contour = ctx->contour;
83 17 : PetscCall(STGetMatStructure(eps->st,&str));
84 17 : PetscCall(STGetSplitPreconditionerInfo(eps->st,&nsplit,&strp));
85 389 : for (i=0;i<contour->npoints;i++) {
86 372 : p_id = i*contour->subcomm->n + contour->subcomm->color;
87 372 : PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&Amat));
88 372 : if (B) PetscCall(MatAXPY(Amat,-ctx->omega[p_id],B,str));
89 264 : else PetscCall(MatShift(Amat,-ctx->omega[p_id]));
90 372 : if (nsplit) {
91 64 : PetscCall(MatDuplicate(Pa,MAT_COPY_VALUES,&Pmat));
92 64 : if (Pb) PetscCall(MatAXPY(Pmat,-ctx->omega[p_id],Pb,strp));
93 0 : else PetscCall(MatShift(Pmat,-ctx->omega[p_id]));
94 308 : } else Pmat = Amat;
95 372 : PetscCall(EPS_KSPSetOperators(contour->ksp[i],Amat,Amat));
96 372 : PetscCall(MatDestroy(&Amat));
97 372 : if (nsplit) PetscCall(MatDestroy(&Pmat));
98 : }
99 17 : PetscFunctionReturn(PETSC_SUCCESS);
100 : }
101 :
102 : /*
103 : Y_i = (A-z_i B)^{-1}BV for every integration point, Y=[Y_i] is in the context
104 : */
105 26 : static PetscErrorCode EPSCISSSolve(EPS eps,Mat B,BV V,PetscInt L_start,PetscInt L_end)
106 : {
107 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
108 26 : SlepcContourData contour;
109 26 : PetscInt i,p_id;
110 26 : Mat MV,BMV=NULL,MC;
111 26 : KSP ksp;
112 :
113 26 : PetscFunctionBegin;
114 26 : if (!ctx->contour || !ctx->contour->ksp) PetscCall(EPSCISSGetKSPs(eps,NULL,NULL));
115 26 : contour = ctx->contour;
116 26 : PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Something went wrong with EPSCISSGetKSPs()");
117 26 : PetscCall(BVSetActiveColumns(V,L_start,L_end));
118 26 : PetscCall(BVGetMat(V,&MV));
119 678 : for (i=0;i<contour->npoints;i++) {
120 652 : p_id = i*contour->subcomm->n + contour->subcomm->color;
121 652 : if (ctx->usest) {
122 280 : PetscCall(STSetShift(eps->st,ctx->omega[p_id]));
123 280 : PetscCall(STGetKSP(eps->st,&ksp));
124 372 : } else ksp = contour->ksp[i];
125 652 : PetscCall(BVSetActiveColumns(ctx->Y,i*ctx->L+L_start,i*ctx->L+L_end));
126 652 : PetscCall(BVGetMat(ctx->Y,&MC));
127 652 : if (B) {
128 324 : if (!i) {
129 13 : PetscCall(MatProductCreate(B,MV,NULL,&BMV));
130 13 : PetscCall(MatProductSetType(BMV,MATPRODUCT_AB));
131 13 : PetscCall(MatProductSetFromOptions(BMV));
132 13 : PetscCall(MatProductSymbolic(BMV));
133 : }
134 324 : PetscCall(MatProductNumeric(BMV));
135 324 : PetscCall(KSPMatSolve(ksp,BMV,MC));
136 328 : } else PetscCall(KSPMatSolve(ksp,MV,MC));
137 652 : PetscCall(BVRestoreMat(ctx->Y,&MC));
138 652 : if (ctx->usest && i<contour->npoints-1) PetscCall(KSPReset(ksp));
139 : }
140 26 : PetscCall(MatDestroy(&BMV));
141 26 : PetscCall(BVRestoreMat(V,&MV));
142 26 : PetscFunctionReturn(PETSC_SUCCESS);
143 : }
144 :
145 52 : static PetscErrorCode rescale_eig(EPS eps,PetscInt nv)
146 : {
147 52 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
148 52 : PetscInt i;
149 52 : PetscScalar center;
150 52 : PetscReal radius,a,b,c,d,rgscale;
151 : #if defined(PETSC_USE_COMPLEX)
152 : PetscReal start_ang,end_ang,vscale,theta;
153 : #endif
154 52 : PetscBool isring,isellipse,isinterval;
155 :
156 52 : PetscFunctionBegin;
157 52 : PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse));
158 52 : PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring));
159 52 : PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval));
160 52 : PetscCall(RGGetScale(eps->rg,&rgscale));
161 52 : if (isinterval) {
162 18 : PetscCall(RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d));
163 18 : if (c==d) {
164 242 : for (i=0;i<nv;i++) {
165 : #if defined(PETSC_USE_COMPLEX)
166 : eps->eigr[i] = PetscRealPart(eps->eigr[i]);
167 : #else
168 224 : eps->eigi[i] = 0;
169 : #endif
170 : }
171 : }
172 : }
173 52 : if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
174 10 : if (isellipse) {
175 4 : PetscCall(RGEllipseGetParameters(eps->rg,¢er,&radius,NULL));
176 28 : for (i=0;i<nv;i++) eps->eigr[i] = rgscale*(center + radius*eps->eigr[i]);
177 6 : } else if (isinterval) {
178 6 : PetscCall(RGIntervalGetEndpoints(eps->rg,&a,&b,&c,&d));
179 6 : if (ctx->quad == EPS_CISS_QUADRULE_CHEBYSHEV) {
180 24 : for (i=0;i<nv;i++) {
181 20 : if (c==d) eps->eigr[i] = ((b-a)*(eps->eigr[i]+1.0)/2.0+a)*rgscale;
182 20 : if (a==b) {
183 : #if defined(PETSC_USE_COMPLEX)
184 : eps->eigr[i] = ((d-c)*(eps->eigr[i]+1.0)/2.0+c)*rgscale*PETSC_i;
185 : #else
186 20 : SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Integration points on a vertical line require complex arithmetic");
187 : #endif
188 : }
189 : }
190 : } else {
191 2 : center = (b+a)/2.0+(d+c)/2.0*PETSC_PI;
192 4 : radius = PetscSqrtReal(PetscPowRealInt((b-a)/2.0,2)+PetscPowRealInt((d-c)/2.0,2));
193 14 : for (i=0;i<nv;i++) eps->eigr[i] = center + radius*eps->eigr[i];
194 : }
195 : } else if (isring) { /* only supported in complex scalars */
196 : #if defined(PETSC_USE_COMPLEX)
197 : PetscCall(RGRingGetParameters(eps->rg,¢er,&radius,&vscale,&start_ang,&end_ang,NULL));
198 : if (ctx->quad == EPS_CISS_QUADRULE_CHEBYSHEV) {
199 : for (i=0;i<nv;i++) {
200 : theta = (start_ang*2.0+(end_ang-start_ang)*(PetscRealPart(eps->eigr[i])+1.0))*PETSC_PI;
201 : eps->eigr[i] = rgscale*center + (rgscale*radius+PetscImaginaryPart(eps->eigr[i]))*PetscCMPLX(PetscCosReal(theta),vscale*PetscSinReal(theta));
202 : }
203 : } else {
204 : for (i=0;i<nv;i++) eps->eigr[i] = rgscale*(center + radius*eps->eigr[i]);
205 : }
206 : #endif
207 52 : }
208 : }
209 52 : PetscFunctionReturn(PETSC_SUCCESS);
210 : }
211 :
212 26 : static PetscErrorCode EPSSetUp_CISS(EPS eps)
213 : {
214 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
215 26 : SlepcContourData contour;
216 26 : PetscBool istrivial,isring,isellipse,isinterval,flg;
217 26 : PetscReal c,d;
218 26 : PetscInt nsplit;
219 26 : PetscRandom rand;
220 26 : PetscObjectId id;
221 26 : PetscObjectState state;
222 26 : Mat A[2],Psplit[2];
223 26 : Vec v0;
224 :
225 26 : PetscFunctionBegin;
226 26 : EPSCheckNotStructured(eps);
227 26 : if (eps->ncv==PETSC_DETERMINE) {
228 24 : if (eps->nev==0) eps->nev = 1;
229 24 : eps->ncv = ctx->L_max*ctx->M;
230 24 : if (eps->ncv>eps->n) {
231 15 : eps->ncv = eps->n;
232 15 : ctx->L_max = eps->ncv/ctx->M;
233 15 : PetscCheck(ctx->L_max,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cannot adjust solver parameters, try setting a smaller value of M (moment size)");
234 : }
235 : } else {
236 2 : PetscCall(EPSSetDimensions_Default(eps,&eps->nev,&eps->ncv,&eps->mpd));
237 2 : ctx->L_max = eps->ncv/ctx->M;
238 2 : if (!ctx->L_max) {
239 0 : ctx->L_max = 1;
240 0 : eps->ncv = ctx->L_max*ctx->M;
241 : }
242 : }
243 26 : ctx->L = PetscMin(ctx->L,ctx->L_max);
244 26 : if (eps->max_it==PETSC_DETERMINE) eps->max_it = 5;
245 26 : if (eps->mpd==PETSC_DETERMINE) eps->mpd = eps->ncv;
246 26 : if (!eps->which) eps->which = EPS_ALL;
247 26 : PetscCheck(eps->which==EPS_ALL,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only computing all eigenvalues");
248 26 : EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_EXTRACTION | EPS_FEATURE_STOPPING | EPS_FEATURE_TWOSIDED);
249 :
250 : /* check region */
251 26 : PetscCall(RGIsTrivial(eps->rg,&istrivial));
252 26 : PetscCheck(!istrivial,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
253 26 : PetscCall(RGGetComplement(eps->rg,&flg));
254 26 : PetscCheck(!flg,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
255 26 : PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse));
256 26 : PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring));
257 26 : PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval));
258 26 : PetscCheck(isellipse || isring || isinterval,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Currently only implemented for interval, elliptic or ring regions");
259 :
260 : /* if the region has changed, then reset contour data */
261 26 : PetscCall(PetscObjectGetId((PetscObject)eps->rg,&id));
262 26 : PetscCall(PetscObjectStateGet((PetscObject)eps->rg,&state));
263 26 : if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
264 0 : PetscCall(SlepcContourDataDestroy(&ctx->contour));
265 0 : PetscCall(PetscInfo(eps,"Resetting the contour data structure due to a change of region\n"));
266 0 : ctx->rgid = id; ctx->rgstate = state;
267 : }
268 :
269 : #if !defined(PETSC_USE_COMPLEX)
270 26 : PetscCheck(!isring,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Ring region only supported for complex scalars");
271 : #endif
272 26 : if (isinterval) {
273 9 : PetscCall(RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d));
274 : #if !defined(PETSC_USE_COMPLEX)
275 9 : PetscCheck(c==d && c==0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"In real scalars, endpoints of the imaginary axis must be both zero");
276 : #endif
277 9 : if (!ctx->quad && c==d) ctx->quad = EPS_CISS_QUADRULE_CHEBYSHEV;
278 : }
279 26 : if (!ctx->quad) ctx->quad = EPS_CISS_QUADRULE_TRAPEZOIDAL;
280 :
281 : /* create contour data structure */
282 26 : if (!ctx->contour) {
283 0 : PetscCall(RGCanUseConjugates(eps->rg,ctx->isreal,&ctx->useconj));
284 0 : PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)eps,&ctx->contour));
285 : }
286 :
287 26 : PetscCall(EPSAllocateSolution(eps,0));
288 26 : PetscCall(BVGetRandomContext(eps->V,&rand)); /* make sure the random context is available when duplicating */
289 26 : if (ctx->weight) PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
290 26 : PetscCall(PetscMalloc4(ctx->N,&ctx->weight,ctx->N+1,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma));
291 :
292 : /* allocate basis vectors */
293 26 : PetscCall(BVDestroy(&ctx->S));
294 26 : PetscCall(BVDuplicateResize(eps->V,ctx->L*ctx->M,&ctx->S));
295 26 : PetscCall(BVDestroy(&ctx->V));
296 26 : PetscCall(BVDuplicateResize(eps->V,ctx->L,&ctx->V));
297 :
298 26 : PetscCall(STGetMatrix(eps->st,0,&A[0]));
299 26 : PetscCall(MatIsShell(A[0],&flg));
300 26 : PetscCheck(!flg,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Matrix type shell is not supported in this solver");
301 26 : if (eps->isgeneralized) PetscCall(STGetMatrix(eps->st,1,&A[1]));
302 :
303 26 : if (!ctx->usest_set) ctx->usest = (ctx->npart>1)? PETSC_FALSE: PETSC_TRUE;
304 26 : PetscCheck(!ctx->usest || ctx->npart==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The usest flag is not supported when partitions > 1");
305 :
306 : /* check if a user-defined split preconditioner has been set */
307 26 : PetscCall(STGetSplitPreconditionerInfo(eps->st,&nsplit,NULL));
308 26 : if (nsplit) {
309 4 : PetscCall(STGetSplitPreconditionerTerm(eps->st,0,&Psplit[0]));
310 4 : if (eps->isgeneralized) PetscCall(STGetSplitPreconditionerTerm(eps->st,1,&Psplit[1]));
311 : }
312 :
313 26 : contour = ctx->contour;
314 61 : PetscCall(SlepcContourRedundantMat(contour,eps->isgeneralized?2:1,A,nsplit?Psplit:NULL));
315 26 : if (contour->pA) {
316 8 : PetscCall(BVGetColumn(ctx->V,0,&v0));
317 8 : PetscCall(SlepcContourScatterCreate(contour,v0));
318 8 : PetscCall(BVRestoreColumn(ctx->V,0,&v0));
319 8 : PetscCall(BVDestroy(&ctx->pV));
320 8 : PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV));
321 8 : PetscCall(BVSetSizesFromVec(ctx->pV,contour->xsub,eps->n));
322 8 : PetscCall(BVSetFromOptions(ctx->pV));
323 8 : PetscCall(BVResize(ctx->pV,ctx->L,PETSC_FALSE));
324 : }
325 :
326 26 : EPSCheckDefinite(eps);
327 26 : EPSCheckSinvertCondition(eps,ctx->usest," (with the usest flag set)");
328 :
329 26 : PetscCall(BVDestroy(&ctx->Y));
330 26 : if (contour->pA) {
331 8 : PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y));
332 8 : PetscCall(BVSetSizesFromVec(ctx->Y,contour->xsub,eps->n));
333 8 : PetscCall(BVSetFromOptions(ctx->Y));
334 8 : PetscCall(BVResize(ctx->Y,contour->npoints*ctx->L,PETSC_FALSE));
335 18 : } else PetscCall(BVDuplicateResize(eps->V,contour->npoints*ctx->L,&ctx->Y));
336 :
337 26 : if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) PetscCall(DSSetType(eps->ds,DSGNHEP));
338 21 : else if (eps->isgeneralized) {
339 10 : if (eps->ishermitian && eps->ispositive) PetscCall(DSSetType(eps->ds,DSGHEP));
340 1 : else PetscCall(DSSetType(eps->ds,DSGNHEP));
341 : } else {
342 11 : if (eps->ishermitian) PetscCall(DSSetType(eps->ds,DSHEP));
343 6 : else PetscCall(DSSetType(eps->ds,DSNHEP));
344 : }
345 26 : PetscCall(DSAllocate(eps->ds,eps->ncv));
346 :
347 : #if !defined(PETSC_USE_COMPLEX)
348 26 : PetscCall(EPSSetWorkVecs(eps,3));
349 26 : if (!eps->ishermitian) PetscCall(PetscInfo(eps,"Warning: complex eigenvalues are not calculated exactly without --with-scalar-type=complex in PETSc\n"));
350 : #else
351 : PetscCall(EPSSetWorkVecs(eps,2));
352 : #endif
353 26 : PetscFunctionReturn(PETSC_SUCCESS);
354 : }
355 :
356 26 : static PetscErrorCode EPSSetUpSort_CISS(EPS eps)
357 : {
358 26 : SlepcSC sc;
359 :
360 26 : PetscFunctionBegin;
361 : /* fill sorting criterion context */
362 26 : eps->sc->comparison = SlepcCompareSmallestReal;
363 26 : eps->sc->comparisonctx = NULL;
364 26 : eps->sc->map = NULL;
365 26 : eps->sc->mapobj = NULL;
366 :
367 : /* fill sorting criterion for DS */
368 26 : PetscCall(DSGetSlepcSC(eps->ds,&sc));
369 26 : sc->comparison = SlepcCompareLargestMagnitude;
370 26 : sc->comparisonctx = NULL;
371 26 : sc->map = NULL;
372 26 : sc->mapobj = NULL;
373 26 : PetscFunctionReturn(PETSC_SUCCESS);
374 : }
375 :
376 26 : static PetscErrorCode EPSSolve_CISS(EPS eps)
377 : {
378 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
379 26 : SlepcContourData contour = ctx->contour;
380 26 : Mat A,B,X,M,pA,pB,T,J,Pa=NULL,Pb=NULL;
381 26 : BV V;
382 26 : PetscInt i,j,ld,nmat,L_add=0,nv=0,L_base=ctx->L,inner,nlocal,*inside,nsplit;
383 26 : PetscScalar *Mu,*H0,*H1=NULL,*rr,*temp;
384 26 : PetscReal error,max_error,norm;
385 26 : PetscBool *fl1;
386 26 : Vec si,si1=NULL,w[3];
387 26 : PetscRandom rand;
388 : #if defined(PETSC_USE_COMPLEX)
389 : PetscBool isellipse;
390 : PetscReal est_eig,eta;
391 : #else
392 26 : PetscReal normi;
393 : #endif
394 :
395 26 : PetscFunctionBegin;
396 26 : w[0] = eps->work[0];
397 : #if defined(PETSC_USE_COMPLEX)
398 : w[1] = NULL;
399 : #else
400 26 : w[1] = eps->work[2];
401 : #endif
402 26 : w[2] = eps->work[1];
403 26 : PetscCall(VecGetLocalSize(w[0],&nlocal));
404 26 : PetscCall(DSGetLeadingDimension(eps->ds,&ld));
405 42 : PetscCall(RGComputeQuadrature(eps->rg,ctx->quad==EPS_CISS_QUADRULE_CHEBYSHEV?RG_QUADRULE_CHEBYSHEV:RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight));
406 26 : PetscCall(STGetNumMatrices(eps->st,&nmat));
407 26 : PetscCall(STGetMatrix(eps->st,0,&A));
408 26 : if (nmat>1) PetscCall(STGetMatrix(eps->st,1,&B));
409 13 : else B = NULL;
410 26 : J = (contour->pA && nmat>1)? contour->pA[1]: B;
411 26 : V = contour->pA? ctx->pV: ctx->V;
412 26 : if (!ctx->usest) {
413 17 : T = contour->pA? contour->pA[0]: A;
414 17 : PetscCall(STGetSplitPreconditionerInfo(eps->st,&nsplit,NULL));
415 17 : if (nsplit) {
416 3 : if (contour->pA) {
417 2 : Pa = contour->pP[0];
418 2 : if (nsplit>1) Pb = contour->pP[1];
419 : } else {
420 1 : PetscCall(STGetSplitPreconditionerTerm(eps->st,0,&Pa));
421 1 : if (nsplit>1) PetscCall(STGetSplitPreconditionerTerm(eps->st,1,&Pb));
422 : }
423 : }
424 17 : PetscCall(EPSCISSSetUp(eps,T,J,Pa,Pb));
425 : }
426 26 : PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
427 26 : PetscCall(BVSetRandomSign(ctx->V));
428 26 : PetscCall(BVGetRandomContext(ctx->V,&rand));
429 :
430 26 : if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
431 26 : PetscCall(EPSCISSSolve(eps,J,V,0,ctx->L));
432 : #if defined(PETSC_USE_COMPLEX)
433 : PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse));
434 : if (isellipse) {
435 : PetscCall(BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig));
436 : PetscCall(PetscInfo(eps,"Estimated eigenvalue count: %f\n",(double)est_eig));
437 : eta = PetscPowReal(10.0,-PetscLog10Real(eps->tol)/ctx->N);
438 : L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
439 : if (L_add>ctx->L_max-ctx->L) {
440 : PetscCall(PetscInfo(eps,"Number of eigenvalues inside the contour path may be too large\n"));
441 : L_add = ctx->L_max-ctx->L;
442 : }
443 : }
444 : #endif
445 26 : if (L_add>0) {
446 : PetscCall(PetscInfo(eps,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by Estimate #Eig\n",ctx->L,ctx->L+L_add));
447 : PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
448 : PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
449 : PetscCall(BVSetRandomSign(ctx->V));
450 : if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
451 : ctx->L += L_add;
452 26 : PetscCall(EPSCISSSolve(eps,J,V,ctx->L-L_add,ctx->L));
453 : }
454 26 : PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
455 26 : for (i=0;i<ctx->refine_blocksize;i++) {
456 1 : PetscCall(BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj));
457 1 : PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
458 1 : PetscCall(PetscLogEventBegin(EPS_CISS_SVD,eps,0,0,0));
459 1 : PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
460 1 : PetscCall(PetscLogEventEnd(EPS_CISS_SVD,eps,0,0,0));
461 1 : if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
462 0 : L_add = L_base;
463 0 : if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
464 0 : PetscCall(PetscInfo(eps,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by SVD(H0)\n",ctx->L,ctx->L+L_add));
465 0 : PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
466 0 : PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
467 0 : PetscCall(BVSetRandomSign(ctx->V));
468 0 : if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
469 0 : ctx->L += L_add;
470 0 : PetscCall(EPSCISSSolve(eps,J,V,ctx->L-L_add,ctx->L));
471 0 : if (L_add) {
472 0 : PetscCall(PetscFree2(Mu,H0));
473 0 : PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
474 : }
475 : }
476 26 : if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) PetscCall(PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1));
477 :
478 52 : while (eps->reason == EPS_CONVERGED_ITERATING) {
479 26 : eps->its++;
480 26 : for (inner=0;inner<=ctx->refine_inner;inner++) {
481 26 : if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
482 5 : PetscCall(BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj));
483 5 : PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
484 5 : PetscCall(PetscLogEventBegin(EPS_CISS_SVD,eps,0,0,0));
485 5 : PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
486 5 : PetscCall(PetscLogEventEnd(EPS_CISS_SVD,eps,0,0,0));
487 : break;
488 : } else {
489 21 : PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
490 21 : PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
491 21 : PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
492 21 : PetscCall(BVCopy(ctx->S,ctx->V));
493 21 : PetscCall(BVSVDAndRank(ctx->S,ctx->M,ctx->L,ctx->delta,BV_SVD_METHOD_REFINE,H0,ctx->sigma,&nv));
494 21 : if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
495 0 : if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
496 0 : PetscCall(EPSCISSSolve(eps,J,V,0,ctx->L));
497 : } else break;
498 : }
499 : }
500 26 : eps->nconv = 0;
501 26 : if (nv == 0) eps->reason = EPS_CONVERGED_TOL;
502 : else {
503 26 : PetscCall(DSSetDimensions(eps->ds,nv,0,0));
504 26 : PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
505 :
506 26 : if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
507 5 : PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
508 5 : PetscCall(CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1));
509 5 : PetscCall(DSGetArray(eps->ds,DS_MAT_A,&temp));
510 33 : for (j=0;j<nv;j++) {
511 188 : for (i=0;i<nv;i++) {
512 160 : temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
513 : }
514 : }
515 5 : PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&temp));
516 5 : PetscCall(DSGetArray(eps->ds,DS_MAT_B,&temp));
517 33 : for (j=0;j<nv;j++) {
518 188 : for (i=0;i<nv;i++) {
519 160 : temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
520 : }
521 : }
522 5 : PetscCall(DSRestoreArray(eps->ds,DS_MAT_B,&temp));
523 : } else {
524 21 : PetscCall(BVSetActiveColumns(ctx->S,0,nv));
525 21 : PetscCall(DSGetMat(eps->ds,DS_MAT_A,&pA));
526 21 : PetscCall(MatZeroEntries(pA));
527 21 : PetscCall(BVMatProject(ctx->S,A,ctx->S,pA));
528 21 : PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&pA));
529 21 : if (B) {
530 10 : PetscCall(DSGetMat(eps->ds,DS_MAT_B,&pB));
531 10 : PetscCall(MatZeroEntries(pB));
532 10 : PetscCall(BVMatProject(ctx->S,B,ctx->S,pB));
533 10 : PetscCall(DSRestoreMat(eps->ds,DS_MAT_B,&pB));
534 : }
535 : }
536 :
537 26 : PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi));
538 26 : PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));
539 :
540 26 : PetscCall(PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr));
541 26 : PetscCall(rescale_eig(eps,nv));
542 26 : PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
543 26 : PetscCall(DSGetMat(eps->ds,DS_MAT_X,&X));
544 26 : PetscCall(SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1));
545 26 : PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&X));
546 26 : PetscCall(RGCheckInside(eps->rg,nv,eps->eigr,eps->eigi,inside));
547 549 : for (i=0;i<nv;i++) {
548 523 : if (fl1[i] && inside[i]>=0) {
549 174 : rr[i] = 1.0;
550 174 : eps->nconv++;
551 349 : } else rr[i] = 0.0;
552 : }
553 26 : PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,rr,NULL,&eps->nconv));
554 26 : PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));
555 26 : PetscCall(rescale_eig(eps,nv));
556 26 : PetscCall(PetscFree3(fl1,inside,rr));
557 26 : PetscCall(BVSetActiveColumns(eps->V,0,nv));
558 26 : if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
559 5 : PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
560 5 : PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
561 5 : PetscCall(BVCopy(ctx->S,ctx->V));
562 5 : PetscCall(BVSetActiveColumns(ctx->S,0,nv));
563 : }
564 26 : PetscCall(BVCopy(ctx->S,eps->V));
565 :
566 26 : PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
567 26 : PetscCall(DSGetMat(eps->ds,DS_MAT_X,&X));
568 26 : PetscCall(BVMultInPlace(ctx->S,X,0,eps->nconv));
569 26 : if (eps->ishermitian) PetscCall(BVMultInPlace(eps->V,X,0,eps->nconv));
570 26 : PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&X));
571 : max_error = 0.0;
572 196 : for (i=0;i<eps->nconv;i++) {
573 170 : PetscCall(BVGetColumn(ctx->S,i,&si));
574 : #if !defined(PETSC_USE_COMPLEX)
575 170 : if (eps->eigi[i]!=0.0) PetscCall(BVGetColumn(ctx->S,i+1,&si1));
576 : #endif
577 170 : PetscCall(EPSComputeResidualNorm_Private(eps,PETSC_FALSE,eps->eigr[i],eps->eigi[i],si,si1,w,&error));
578 170 : if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) { /* vector is not normalized */
579 21 : PetscCall(VecNorm(si,NORM_2,&norm));
580 : #if !defined(PETSC_USE_COMPLEX)
581 21 : if (eps->eigi[i]!=0.0) {
582 0 : PetscCall(VecNorm(si1,NORM_2,&normi));
583 0 : norm = SlepcAbsEigenvalue(norm,normi);
584 : }
585 : #endif
586 21 : error /= norm;
587 : }
588 170 : PetscCall((*eps->converged)(eps,eps->eigr[i],eps->eigi[i],error,&error,eps->convergedctx));
589 170 : PetscCall(BVRestoreColumn(ctx->S,i,&si));
590 : #if !defined(PETSC_USE_COMPLEX)
591 170 : if (eps->eigi[i]!=0.0) {
592 4 : PetscCall(BVRestoreColumn(ctx->S,i+1,&si1));
593 : i++;
594 : }
595 : #endif
596 246 : max_error = PetscMax(max_error,error);
597 : }
598 :
599 26 : if (max_error <= eps->tol) eps->reason = EPS_CONVERGED_TOL;
600 0 : else if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
601 : else {
602 0 : if (eps->nconv > ctx->L) nv = eps->nconv;
603 0 : else if (ctx->L > nv) nv = ctx->L;
604 0 : nv = PetscMin(nv,ctx->L*ctx->M);
605 0 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M));
606 0 : PetscCall(MatSetRandom(M,rand));
607 0 : PetscCall(BVSetActiveColumns(ctx->S,0,nv));
608 0 : PetscCall(BVMultInPlace(ctx->S,M,0,ctx->L));
609 0 : PetscCall(MatDestroy(&M));
610 0 : PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
611 0 : PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
612 0 : PetscCall(BVCopy(ctx->S,ctx->V));
613 0 : if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
614 52 : PetscCall(EPSCISSSolve(eps,J,V,0,ctx->L));
615 : }
616 : }
617 : }
618 26 : if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) PetscCall(PetscFree(H1));
619 26 : PetscCall(PetscFree2(Mu,H0));
620 26 : PetscFunctionReturn(PETSC_SUCCESS);
621 : }
622 :
623 26 : static PetscErrorCode EPSComputeVectors_CISS(EPS eps)
624 : {
625 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
626 26 : PetscInt n;
627 26 : Mat Z,B=NULL;
628 :
629 26 : PetscFunctionBegin;
630 26 : if (eps->ishermitian) {
631 16 : if (eps->isgeneralized && !eps->ispositive) PetscCall(EPSComputeVectors_Indefinite(eps));
632 16 : else PetscCall(EPSComputeVectors_Hermitian(eps));
633 16 : if (eps->isgeneralized && eps->ispositive && ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
634 : /* normalize to have unit B-norm */
635 2 : PetscCall(STGetMatrix(eps->st,1,&B));
636 2 : PetscCall(BVSetMatrix(eps->V,B,PETSC_FALSE));
637 2 : PetscCall(BVNormalize(eps->V,NULL));
638 2 : PetscCall(BVSetMatrix(eps->V,NULL,PETSC_FALSE));
639 : }
640 16 : PetscFunctionReturn(PETSC_SUCCESS);
641 : }
642 10 : PetscCall(DSGetDimensions(eps->ds,&n,NULL,NULL,NULL));
643 10 : PetscCall(BVSetActiveColumns(eps->V,0,n));
644 :
645 : /* right eigenvectors */
646 10 : PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
647 :
648 : /* V = V * Z */
649 10 : PetscCall(DSGetMat(eps->ds,DS_MAT_X,&Z));
650 10 : PetscCall(BVMultInPlace(eps->V,Z,0,n));
651 10 : PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&Z));
652 10 : PetscCall(BVSetActiveColumns(eps->V,0,eps->nconv));
653 :
654 : /* normalize */
655 10 : if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) PetscCall(BVNormalize(eps->V,NULL));
656 10 : PetscFunctionReturn(PETSC_SUCCESS);
657 : }
658 :
659 13 : static PetscErrorCode EPSCISSSetSizes_CISS(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
660 : {
661 13 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
662 13 : PetscInt oN,oL,oM,oLmax,onpart;
663 13 : PetscMPIInt size;
664 :
665 13 : PetscFunctionBegin;
666 13 : oN = ctx->N;
667 13 : if (ip == PETSC_DETERMINE) {
668 0 : if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
669 13 : } else if (ip != PETSC_CURRENT) {
670 13 : PetscCheck(ip>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
671 13 : PetscCheck(ip%2==0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
672 13 : if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
673 : }
674 13 : oL = ctx->L;
675 13 : if (bs == PETSC_DETERMINE) {
676 0 : ctx->L = 16;
677 13 : } else if (bs != PETSC_CURRENT) {
678 13 : PetscCheck(bs>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
679 13 : ctx->L = bs;
680 : }
681 13 : oM = ctx->M;
682 13 : if (ms == PETSC_DETERMINE) {
683 0 : ctx->M = ctx->N/4;
684 13 : } else if (ms != PETSC_CURRENT) {
685 13 : PetscCheck(ms>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
686 13 : PetscCheck(ms<=ctx->N,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be less than or equal to the number of integration points");
687 13 : ctx->M = ms;
688 : }
689 13 : onpart = ctx->npart;
690 13 : if (npart == PETSC_DETERMINE) {
691 0 : ctx->npart = 1;
692 13 : } else if (npart != PETSC_CURRENT) {
693 13 : PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size));
694 13 : PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
695 13 : ctx->npart = npart;
696 : }
697 13 : oLmax = ctx->L_max;
698 13 : if (bsmax == PETSC_DETERMINE) {
699 0 : ctx->L_max = 64;
700 13 : } else if (bsmax != PETSC_CURRENT) {
701 13 : PetscCheck(bsmax>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
702 13 : ctx->L_max = PetscMax(bsmax,ctx->L);
703 : }
704 13 : if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
705 12 : PetscCall(SlepcContourDataDestroy(&ctx->contour));
706 12 : PetscCall(PetscInfo(eps,"Resetting the contour data structure due to a change of parameters\n"));
707 12 : eps->state = EPS_STATE_INITIAL;
708 : }
709 13 : ctx->isreal = realmats;
710 13 : if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) eps->state = EPS_STATE_INITIAL;
711 13 : PetscFunctionReturn(PETSC_SUCCESS);
712 : }
713 :
714 : /*@
715 : EPSCISSSetSizes - Sets the values of various size parameters in the CISS solver.
716 :
717 : Logically Collective
718 :
719 : Input Parameters:
720 : + eps - the eigenproblem solver context
721 : . ip - number of integration points
722 : . bs - block size
723 : . ms - moment size
724 : . npart - number of partitions when splitting the communicator
725 : . bsmax - max block size
726 : - realmats - A and B are real
727 :
728 : Options Database Keys:
729 : + -eps_ciss_integration_points - Sets the number of integration points
730 : . -eps_ciss_blocksize - Sets the block size
731 : . -eps_ciss_moments - Sets the moment size
732 : . -eps_ciss_partitions - Sets the number of partitions
733 : . -eps_ciss_maxblocksize - Sets the maximum block size
734 : - -eps_ciss_realmats - A and B are real
735 :
736 : Note:
737 : For all integer arguments, you can use PETSC_CURRENT to keep the current value, and
738 : PETSC_DETERMINE to set them to a default value.
739 :
740 : The default number of partitions is 1. This means the internal KSP object is shared
741 : among all processes of the EPS communicator. Otherwise, the communicator is split
742 : into npart communicators, so that npart KSP solves proceed simultaneously.
743 :
744 : Level: advanced
745 :
746 : .seealso: EPSCISSGetSizes()
747 : @*/
748 13 : PetscErrorCode EPSCISSSetSizes(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
749 : {
750 13 : PetscFunctionBegin;
751 13 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
752 39 : PetscValidLogicalCollectiveInt(eps,ip,2);
753 39 : PetscValidLogicalCollectiveInt(eps,bs,3);
754 39 : PetscValidLogicalCollectiveInt(eps,ms,4);
755 39 : PetscValidLogicalCollectiveInt(eps,npart,5);
756 39 : PetscValidLogicalCollectiveInt(eps,bsmax,6);
757 39 : PetscValidLogicalCollectiveBool(eps,realmats,7);
758 13 : PetscTryMethod(eps,"EPSCISSSetSizes_C",(EPS,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(eps,ip,bs,ms,npart,bsmax,realmats));
759 13 : PetscFunctionReturn(PETSC_SUCCESS);
760 : }
761 :
762 26 : static PetscErrorCode EPSCISSGetSizes_CISS(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
763 : {
764 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
765 :
766 26 : PetscFunctionBegin;
767 26 : if (ip) *ip = ctx->N;
768 26 : if (bs) *bs = ctx->L;
769 26 : if (ms) *ms = ctx->M;
770 26 : if (npart) *npart = ctx->npart;
771 26 : if (bsmax) *bsmax = ctx->L_max;
772 26 : if (realmats) *realmats = ctx->isreal;
773 26 : PetscFunctionReturn(PETSC_SUCCESS);
774 : }
775 :
776 : /*@
777 : EPSCISSGetSizes - Gets the values of various size parameters in the CISS solver.
778 :
779 : Not Collective
780 :
781 : Input Parameter:
782 : . eps - the eigenproblem solver context
783 :
784 : Output Parameters:
785 : + ip - number of integration points
786 : . bs - block size
787 : . ms - moment size
788 : . npart - number of partitions when splitting the communicator
789 : . bsmax - max block size
790 : - realmats - A and B are real
791 :
792 : Level: advanced
793 :
794 : .seealso: EPSCISSSetSizes()
795 : @*/
796 26 : PetscErrorCode EPSCISSGetSizes(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
797 : {
798 26 : PetscFunctionBegin;
799 26 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
800 26 : PetscUseMethod(eps,"EPSCISSGetSizes_C",(EPS,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(eps,ip,bs,ms,npart,bsmax,realmats));
801 26 : PetscFunctionReturn(PETSC_SUCCESS);
802 : }
803 :
804 2 : static PetscErrorCode EPSCISSSetThreshold_CISS(EPS eps,PetscReal delta,PetscReal spur)
805 : {
806 2 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
807 :
808 2 : PetscFunctionBegin;
809 2 : if (delta == (PetscReal)PETSC_DETERMINE) {
810 0 : ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
811 2 : } else if (delta != (PetscReal)PETSC_CURRENT) {
812 2 : PetscCheck(delta>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
813 2 : ctx->delta = delta;
814 : }
815 2 : if (spur == (PetscReal)PETSC_DETERMINE) {
816 0 : ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
817 2 : } else if (spur != (PetscReal)PETSC_CURRENT) {
818 2 : PetscCheck(spur>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
819 2 : ctx->spurious_threshold = spur;
820 : }
821 2 : PetscFunctionReturn(PETSC_SUCCESS);
822 : }
823 :
824 : /*@
825 : EPSCISSSetThreshold - Sets the values of various threshold parameters in
826 : the CISS solver.
827 :
828 : Logically Collective
829 :
830 : Input Parameters:
831 : + eps - the eigenproblem solver context
832 : . delta - threshold for numerical rank
833 : - spur - spurious threshold (to discard spurious eigenpairs)
834 :
835 : Options Database Keys:
836 : + -eps_ciss_delta - Sets the delta
837 : - -eps_ciss_spurious_threshold - Sets the spurious threshold
838 :
839 : Note:
840 : PETSC_CURRENT can be used to preserve the current value of any of the
841 : arguments, and PETSC_DETERMINE to set them to a default value.
842 :
843 : Level: advanced
844 :
845 : .seealso: EPSCISSGetThreshold()
846 : @*/
847 2 : PetscErrorCode EPSCISSSetThreshold(EPS eps,PetscReal delta,PetscReal spur)
848 : {
849 2 : PetscFunctionBegin;
850 2 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
851 6 : PetscValidLogicalCollectiveReal(eps,delta,2);
852 6 : PetscValidLogicalCollectiveReal(eps,spur,3);
853 2 : PetscTryMethod(eps,"EPSCISSSetThreshold_C",(EPS,PetscReal,PetscReal),(eps,delta,spur));
854 2 : PetscFunctionReturn(PETSC_SUCCESS);
855 : }
856 :
857 26 : static PetscErrorCode EPSCISSGetThreshold_CISS(EPS eps,PetscReal *delta,PetscReal *spur)
858 : {
859 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
860 :
861 26 : PetscFunctionBegin;
862 26 : if (delta) *delta = ctx->delta;
863 26 : if (spur) *spur = ctx->spurious_threshold;
864 26 : PetscFunctionReturn(PETSC_SUCCESS);
865 : }
866 :
867 : /*@
868 : EPSCISSGetThreshold - Gets the values of various threshold parameters
869 : in the CISS solver.
870 :
871 : Not Collective
872 :
873 : Input Parameter:
874 : . eps - the eigenproblem solver context
875 :
876 : Output Parameters:
877 : + delta - threshold for numerical rank
878 : - spur - spurious threshold (to discard spurious eigenpairs)
879 :
880 : Level: advanced
881 :
882 : .seealso: EPSCISSSetThreshold()
883 : @*/
884 26 : PetscErrorCode EPSCISSGetThreshold(EPS eps,PetscReal *delta,PetscReal *spur)
885 : {
886 26 : PetscFunctionBegin;
887 26 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
888 26 : PetscUseMethod(eps,"EPSCISSGetThreshold_C",(EPS,PetscReal*,PetscReal*),(eps,delta,spur));
889 26 : PetscFunctionReturn(PETSC_SUCCESS);
890 : }
891 :
892 1 : static PetscErrorCode EPSCISSSetRefinement_CISS(EPS eps,PetscInt inner,PetscInt blsize)
893 : {
894 1 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
895 :
896 1 : PetscFunctionBegin;
897 1 : if (inner == PETSC_DETERMINE) {
898 0 : ctx->refine_inner = 0;
899 1 : } else if (inner != PETSC_CURRENT) {
900 1 : PetscCheck(inner>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
901 1 : ctx->refine_inner = inner;
902 : }
903 1 : if (blsize == PETSC_DETERMINE) {
904 0 : ctx->refine_blocksize = 0;
905 1 : } else if (blsize != PETSC_CURRENT) {
906 1 : PetscCheck(blsize>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
907 1 : ctx->refine_blocksize = blsize;
908 : }
909 1 : PetscFunctionReturn(PETSC_SUCCESS);
910 : }
911 :
912 : /*@
913 : EPSCISSSetRefinement - Sets the values of various refinement parameters
914 : in the CISS solver.
915 :
916 : Logically Collective
917 :
918 : Input Parameters:
919 : + eps - the eigenproblem solver context
920 : . inner - number of iterative refinement iterations (inner loop)
921 : - blsize - number of iterative refinement iterations (blocksize loop)
922 :
923 : Options Database Keys:
924 : + -eps_ciss_refine_inner - Sets number of inner iterations
925 : - -eps_ciss_refine_blocksize - Sets number of blocksize iterations
926 :
927 : Note:
928 : PETSC_CURRENT can be used to preserve the current value of any of the
929 : arguments, and PETSC_DETERMINE to set them to a default of 0.
930 :
931 : Level: advanced
932 :
933 : .seealso: EPSCISSGetRefinement()
934 : @*/
935 1 : PetscErrorCode EPSCISSSetRefinement(EPS eps,PetscInt inner,PetscInt blsize)
936 : {
937 1 : PetscFunctionBegin;
938 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
939 3 : PetscValidLogicalCollectiveInt(eps,inner,2);
940 3 : PetscValidLogicalCollectiveInt(eps,blsize,3);
941 1 : PetscTryMethod(eps,"EPSCISSSetRefinement_C",(EPS,PetscInt,PetscInt),(eps,inner,blsize));
942 1 : PetscFunctionReturn(PETSC_SUCCESS);
943 : }
944 :
945 26 : static PetscErrorCode EPSCISSGetRefinement_CISS(EPS eps,PetscInt *inner,PetscInt *blsize)
946 : {
947 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
948 :
949 26 : PetscFunctionBegin;
950 26 : if (inner) *inner = ctx->refine_inner;
951 26 : if (blsize) *blsize = ctx->refine_blocksize;
952 26 : PetscFunctionReturn(PETSC_SUCCESS);
953 : }
954 :
955 : /*@
956 : EPSCISSGetRefinement - Gets the values of various refinement parameters
957 : in the CISS solver.
958 :
959 : Not Collective
960 :
961 : Input Parameter:
962 : . eps - the eigenproblem solver context
963 :
964 : Output Parameters:
965 : + inner - number of iterative refinement iterations (inner loop)
966 : - blsize - number of iterative refinement iterations (blocksize loop)
967 :
968 : Level: advanced
969 :
970 : .seealso: EPSCISSSetRefinement()
971 : @*/
972 26 : PetscErrorCode EPSCISSGetRefinement(EPS eps, PetscInt *inner, PetscInt *blsize)
973 : {
974 26 : PetscFunctionBegin;
975 26 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
976 26 : PetscUseMethod(eps,"EPSCISSGetRefinement_C",(EPS,PetscInt*,PetscInt*),(eps,inner,blsize));
977 26 : PetscFunctionReturn(PETSC_SUCCESS);
978 : }
979 :
980 10 : static PetscErrorCode EPSCISSSetUseST_CISS(EPS eps,PetscBool usest)
981 : {
982 10 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
983 :
984 10 : PetscFunctionBegin;
985 10 : ctx->usest = usest;
986 10 : ctx->usest_set = PETSC_TRUE;
987 10 : eps->state = EPS_STATE_INITIAL;
988 10 : PetscFunctionReturn(PETSC_SUCCESS);
989 : }
990 :
991 : /*@
992 : EPSCISSSetUseST - Sets a flag indicating that the CISS solver will
993 : use the ST object for the linear solves.
994 :
995 : Logically Collective
996 :
997 : Input Parameters:
998 : + eps - the eigenproblem solver context
999 : - usest - boolean flag to use the ST object or not
1000 :
1001 : Options Database Keys:
1002 : . -eps_ciss_usest <bool> - whether the ST object will be used or not
1003 :
1004 : Level: advanced
1005 :
1006 : .seealso: EPSCISSGetUseST()
1007 : @*/
1008 10 : PetscErrorCode EPSCISSSetUseST(EPS eps,PetscBool usest)
1009 : {
1010 10 : PetscFunctionBegin;
1011 10 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1012 30 : PetscValidLogicalCollectiveBool(eps,usest,2);
1013 10 : PetscTryMethod(eps,"EPSCISSSetUseST_C",(EPS,PetscBool),(eps,usest));
1014 10 : PetscFunctionReturn(PETSC_SUCCESS);
1015 : }
1016 :
1017 26 : static PetscErrorCode EPSCISSGetUseST_CISS(EPS eps,PetscBool *usest)
1018 : {
1019 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1020 :
1021 26 : PetscFunctionBegin;
1022 26 : *usest = ctx->usest;
1023 26 : PetscFunctionReturn(PETSC_SUCCESS);
1024 : }
1025 :
1026 : /*@
1027 : EPSCISSGetUseST - Gets the flag for using the ST object
1028 : in the CISS solver.
1029 :
1030 : Not Collective
1031 :
1032 : Input Parameter:
1033 : . eps - the eigenproblem solver context
1034 :
1035 : Output Parameters:
1036 : . usest - boolean flag indicating if the ST object is being used
1037 :
1038 : Level: advanced
1039 :
1040 : .seealso: EPSCISSSetUseST()
1041 : @*/
1042 26 : PetscErrorCode EPSCISSGetUseST(EPS eps,PetscBool *usest)
1043 : {
1044 26 : PetscFunctionBegin;
1045 26 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1046 26 : PetscAssertPointer(usest,2);
1047 26 : PetscUseMethod(eps,"EPSCISSGetUseST_C",(EPS,PetscBool*),(eps,usest));
1048 26 : PetscFunctionReturn(PETSC_SUCCESS);
1049 : }
1050 :
1051 5 : static PetscErrorCode EPSCISSSetQuadRule_CISS(EPS eps,EPSCISSQuadRule quad)
1052 : {
1053 5 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1054 :
1055 5 : PetscFunctionBegin;
1056 5 : if (ctx->quad != quad) {
1057 5 : ctx->quad = quad;
1058 5 : eps->state = EPS_STATE_INITIAL;
1059 : }
1060 5 : PetscFunctionReturn(PETSC_SUCCESS);
1061 : }
1062 :
1063 : /*@
1064 : EPSCISSSetQuadRule - Sets the quadrature rule used in the CISS solver.
1065 :
1066 : Logically Collective
1067 :
1068 : Input Parameters:
1069 : + eps - the eigenproblem solver context
1070 : - quad - the quadrature rule
1071 :
1072 : Options Database Key:
1073 : . -eps_ciss_quadrule - Sets the quadrature rule (either 'trapezoidal' or
1074 : 'chebyshev')
1075 :
1076 : Notes:
1077 : By default, the trapezoidal rule is used (EPS_CISS_QUADRULE_TRAPEZOIDAL).
1078 :
1079 : If the 'chebyshev' option is specified (EPS_CISS_QUADRULE_CHEBYSHEV), then
1080 : Chebyshev points are used as quadrature points.
1081 :
1082 : Level: advanced
1083 :
1084 : .seealso: EPSCISSGetQuadRule(), EPSCISSQuadRule
1085 : @*/
1086 5 : PetscErrorCode EPSCISSSetQuadRule(EPS eps,EPSCISSQuadRule quad)
1087 : {
1088 5 : PetscFunctionBegin;
1089 5 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1090 15 : PetscValidLogicalCollectiveEnum(eps,quad,2);
1091 5 : PetscTryMethod(eps,"EPSCISSSetQuadRule_C",(EPS,EPSCISSQuadRule),(eps,quad));
1092 5 : PetscFunctionReturn(PETSC_SUCCESS);
1093 : }
1094 :
1095 1 : static PetscErrorCode EPSCISSGetQuadRule_CISS(EPS eps,EPSCISSQuadRule *quad)
1096 : {
1097 1 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1098 :
1099 1 : PetscFunctionBegin;
1100 1 : *quad = ctx->quad;
1101 1 : PetscFunctionReturn(PETSC_SUCCESS);
1102 : }
1103 :
1104 : /*@
1105 : EPSCISSGetQuadRule - Gets the quadrature rule used in the CISS solver.
1106 :
1107 : Not Collective
1108 :
1109 : Input Parameter:
1110 : . eps - the eigenproblem solver context
1111 :
1112 : Output Parameters:
1113 : . quad - quadrature rule
1114 :
1115 : Level: advanced
1116 :
1117 : .seealso: EPSCISSSetQuadRule() EPSCISSQuadRule
1118 : @*/
1119 1 : PetscErrorCode EPSCISSGetQuadRule(EPS eps,EPSCISSQuadRule *quad)
1120 : {
1121 1 : PetscFunctionBegin;
1122 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1123 1 : PetscAssertPointer(quad,2);
1124 1 : PetscUseMethod(eps,"EPSCISSGetQuadRule_C",(EPS,EPSCISSQuadRule*),(eps,quad));
1125 1 : PetscFunctionReturn(PETSC_SUCCESS);
1126 : }
1127 :
1128 6 : static PetscErrorCode EPSCISSSetExtraction_CISS(EPS eps,EPSCISSExtraction extraction)
1129 : {
1130 6 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1131 :
1132 6 : PetscFunctionBegin;
1133 6 : if (ctx->extraction != extraction) {
1134 5 : ctx->extraction = extraction;
1135 5 : eps->state = EPS_STATE_INITIAL;
1136 : }
1137 6 : PetscFunctionReturn(PETSC_SUCCESS);
1138 : }
1139 :
1140 : /*@
1141 : EPSCISSSetExtraction - Sets the extraction technique used in the CISS solver.
1142 :
1143 : Logically Collective
1144 :
1145 : Input Parameters:
1146 : + eps - the eigenproblem solver context
1147 : - extraction - the extraction technique
1148 :
1149 : Options Database Key:
1150 : . -eps_ciss_extraction - Sets the extraction technique (either 'ritz' or
1151 : 'hankel')
1152 :
1153 : Notes:
1154 : By default, the Rayleigh-Ritz extraction is used (EPS_CISS_EXTRACTION_RITZ).
1155 :
1156 : If the 'hankel' option is specified (EPS_CISS_EXTRACTION_HANKEL), then
1157 : the Block Hankel method is used for extracting eigenpairs.
1158 :
1159 : Level: advanced
1160 :
1161 : .seealso: EPSCISSGetExtraction(), EPSCISSExtraction
1162 : @*/
1163 6 : PetscErrorCode EPSCISSSetExtraction(EPS eps,EPSCISSExtraction extraction)
1164 : {
1165 6 : PetscFunctionBegin;
1166 6 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1167 18 : PetscValidLogicalCollectiveEnum(eps,extraction,2);
1168 6 : PetscTryMethod(eps,"EPSCISSSetExtraction_C",(EPS,EPSCISSExtraction),(eps,extraction));
1169 6 : PetscFunctionReturn(PETSC_SUCCESS);
1170 : }
1171 :
1172 1 : static PetscErrorCode EPSCISSGetExtraction_CISS(EPS eps,EPSCISSExtraction *extraction)
1173 : {
1174 1 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1175 :
1176 1 : PetscFunctionBegin;
1177 1 : *extraction = ctx->extraction;
1178 1 : PetscFunctionReturn(PETSC_SUCCESS);
1179 : }
1180 :
1181 : /*@
1182 : EPSCISSGetExtraction - Gets the extraction technique used in the CISS solver.
1183 :
1184 : Not Collective
1185 :
1186 : Input Parameter:
1187 : . eps - the eigenproblem solver context
1188 :
1189 : Output Parameters:
1190 : . extraction - extraction technique
1191 :
1192 : Level: advanced
1193 :
1194 : .seealso: EPSCISSSetExtraction() EPSCISSExtraction
1195 : @*/
1196 1 : PetscErrorCode EPSCISSGetExtraction(EPS eps,EPSCISSExtraction *extraction)
1197 : {
1198 1 : PetscFunctionBegin;
1199 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1200 1 : PetscAssertPointer(extraction,2);
1201 1 : PetscUseMethod(eps,"EPSCISSGetExtraction_C",(EPS,EPSCISSExtraction*),(eps,extraction));
1202 1 : PetscFunctionReturn(PETSC_SUCCESS);
1203 : }
1204 :
1205 26 : static PetscErrorCode EPSCISSGetKSPs_CISS(EPS eps,PetscInt *nsolve,KSP **ksp)
1206 : {
1207 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1208 26 : SlepcContourData contour;
1209 26 : PetscInt i,nsplit;
1210 26 : PC pc;
1211 26 : MPI_Comm child;
1212 :
1213 26 : PetscFunctionBegin;
1214 26 : if (!ctx->contour) { /* initialize contour data structure first */
1215 26 : PetscCall(RGCanUseConjugates(eps->rg,ctx->isreal,&ctx->useconj));
1216 26 : PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)eps,&ctx->contour));
1217 : }
1218 26 : contour = ctx->contour;
1219 26 : if (!contour->ksp) {
1220 26 : PetscCall(PetscMalloc1(contour->npoints,&contour->ksp));
1221 26 : PetscCall(EPSGetST(eps,&eps->st));
1222 26 : PetscCall(STGetSplitPreconditionerInfo(eps->st,&nsplit,NULL));
1223 26 : PetscCall(PetscSubcommGetChild(contour->subcomm,&child));
1224 678 : for (i=0;i<contour->npoints;i++) {
1225 652 : PetscCall(KSPCreate(child,&contour->ksp[i]));
1226 652 : PetscCall(PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)eps,1));
1227 652 : PetscCall(KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)eps)->prefix));
1228 652 : PetscCall(KSPAppendOptionsPrefix(contour->ksp[i],"eps_ciss_"));
1229 652 : PetscCall(PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)eps)->options));
1230 652 : PetscCall(KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE));
1231 908 : PetscCall(KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(eps->tol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
1232 652 : PetscCall(KSPGetPC(contour->ksp[i],&pc));
1233 652 : if (nsplit) {
1234 96 : PetscCall(KSPSetType(contour->ksp[i],KSPBCGS));
1235 96 : PetscCall(PCSetType(pc,PCBJACOBI));
1236 : } else {
1237 556 : PetscCall(KSPSetType(contour->ksp[i],KSPPREONLY));
1238 652 : PetscCall(PCSetType(pc,PCLU));
1239 : }
1240 : }
1241 : }
1242 26 : if (nsolve) *nsolve = contour->npoints;
1243 26 : if (ksp) *ksp = contour->ksp;
1244 26 : PetscFunctionReturn(PETSC_SUCCESS);
1245 : }
1246 :
1247 : /*@C
1248 : EPSCISSGetKSPs - Retrieve the array of linear solver objects associated with
1249 : the CISS solver.
1250 :
1251 : Not Collective
1252 :
1253 : Input Parameter:
1254 : . eps - the eigenproblem solver solver
1255 :
1256 : Output Parameters:
1257 : + nsolve - number of solver objects
1258 : - ksp - array of linear solver object
1259 :
1260 : Notes:
1261 : The number of KSP solvers is equal to the number of integration points divided by
1262 : the number of partitions. This value is halved in the case of real matrices with
1263 : a region centered at the real axis.
1264 :
1265 : Level: advanced
1266 :
1267 : .seealso: EPSCISSSetSizes()
1268 : @*/
1269 26 : PetscErrorCode EPSCISSGetKSPs(EPS eps,PetscInt *nsolve,KSP **ksp)
1270 : {
1271 26 : PetscFunctionBegin;
1272 26 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1273 26 : PetscUseMethod(eps,"EPSCISSGetKSPs_C",(EPS,PetscInt*,KSP**),(eps,nsolve,ksp));
1274 26 : PetscFunctionReturn(PETSC_SUCCESS);
1275 : }
1276 :
1277 26 : static PetscErrorCode EPSReset_CISS(EPS eps)
1278 : {
1279 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1280 :
1281 26 : PetscFunctionBegin;
1282 26 : PetscCall(BVDestroy(&ctx->S));
1283 26 : PetscCall(BVDestroy(&ctx->V));
1284 26 : PetscCall(BVDestroy(&ctx->Y));
1285 26 : if (!ctx->usest) PetscCall(SlepcContourDataReset(ctx->contour));
1286 26 : PetscCall(BVDestroy(&ctx->pV));
1287 26 : PetscFunctionReturn(PETSC_SUCCESS);
1288 : }
1289 :
1290 26 : static PetscErrorCode EPSSetFromOptions_CISS(EPS eps,PetscOptionItems *PetscOptionsObject)
1291 : {
1292 26 : PetscReal r3,r4;
1293 26 : PetscInt i,i1,i2,i3,i4,i5,i6,i7;
1294 26 : PetscBool b1,b2,flg,flg2,flg3,flg4,flg5,flg6;
1295 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1296 26 : EPSCISSQuadRule quad;
1297 26 : EPSCISSExtraction extraction;
1298 :
1299 26 : PetscFunctionBegin;
1300 26 : PetscOptionsHeadBegin(PetscOptionsObject,"EPS CISS Options");
1301 :
1302 26 : PetscCall(EPSCISSGetSizes(eps,&i1,&i2,&i3,&i4,&i5,&b1));
1303 26 : PetscCall(PetscOptionsInt("-eps_ciss_integration_points","Number of integration points","EPSCISSSetSizes",i1,&i1,&flg));
1304 26 : PetscCall(PetscOptionsInt("-eps_ciss_blocksize","Block size","EPSCISSSetSizes",i2,&i2,&flg2));
1305 26 : PetscCall(PetscOptionsInt("-eps_ciss_moments","Moment size","EPSCISSSetSizes",i3,&i3,&flg3));
1306 26 : PetscCall(PetscOptionsInt("-eps_ciss_partitions","Number of partitions","EPSCISSSetSizes",i4,&i4,&flg4));
1307 26 : PetscCall(PetscOptionsInt("-eps_ciss_maxblocksize","Maximum block size","EPSCISSSetSizes",i5,&i5,&flg5));
1308 26 : PetscCall(PetscOptionsBool("-eps_ciss_realmats","True if A and B are real","EPSCISSSetSizes",b1,&b1,&flg6));
1309 26 : if (flg || flg2 || flg3 || flg4 || flg5 || flg6) PetscCall(EPSCISSSetSizes(eps,i1,i2,i3,i4,i5,b1));
1310 :
1311 26 : PetscCall(EPSCISSGetThreshold(eps,&r3,&r4));
1312 26 : PetscCall(PetscOptionsReal("-eps_ciss_delta","Threshold for numerical rank","EPSCISSSetThreshold",r3,&r3,&flg));
1313 26 : PetscCall(PetscOptionsReal("-eps_ciss_spurious_threshold","Threshold for the spurious eigenpairs","EPSCISSSetThreshold",r4,&r4,&flg2));
1314 26 : if (flg || flg2) PetscCall(EPSCISSSetThreshold(eps,r3,r4));
1315 :
1316 26 : PetscCall(EPSCISSGetRefinement(eps,&i6,&i7));
1317 26 : PetscCall(PetscOptionsInt("-eps_ciss_refine_inner","Number of inner iterative refinement iterations","EPSCISSSetRefinement",i6,&i6,&flg));
1318 26 : PetscCall(PetscOptionsInt("-eps_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","EPSCISSSetRefinement",i7,&i7,&flg2));
1319 26 : if (flg || flg2) PetscCall(EPSCISSSetRefinement(eps,i6,i7));
1320 :
1321 26 : PetscCall(EPSCISSGetUseST(eps,&b2));
1322 26 : PetscCall(PetscOptionsBool("-eps_ciss_usest","Use ST for linear solves","EPSCISSSetUseST",b2,&b2,&flg));
1323 26 : if (flg) PetscCall(EPSCISSSetUseST(eps,b2));
1324 :
1325 26 : PetscCall(PetscOptionsEnum("-eps_ciss_quadrule","Quadrature rule","EPSCISSSetQuadRule",EPSCISSQuadRules,(PetscEnum)ctx->quad,(PetscEnum*)&quad,&flg));
1326 26 : if (flg) PetscCall(EPSCISSSetQuadRule(eps,quad));
1327 :
1328 26 : PetscCall(PetscOptionsEnum("-eps_ciss_extraction","Extraction technique","EPSCISSSetExtraction",EPSCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg));
1329 26 : if (flg) PetscCall(EPSCISSSetExtraction(eps,extraction));
1330 :
1331 26 : PetscOptionsHeadEnd();
1332 :
1333 26 : if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
1334 26 : PetscCall(RGSetFromOptions(eps->rg)); /* this is necessary here to set useconj */
1335 26 : if (!ctx->contour || !ctx->contour->ksp) PetscCall(EPSCISSGetKSPs(eps,NULL,NULL));
1336 26 : PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Something went wrong with EPSCISSGetKSPs()");
1337 678 : for (i=0;i<ctx->contour->npoints;i++) PetscCall(KSPSetFromOptions(ctx->contour->ksp[i]));
1338 26 : PetscCall(PetscSubcommSetFromOptions(ctx->contour->subcomm));
1339 26 : PetscFunctionReturn(PETSC_SUCCESS);
1340 : }
1341 :
1342 26 : static PetscErrorCode EPSDestroy_CISS(EPS eps)
1343 : {
1344 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1345 :
1346 26 : PetscFunctionBegin;
1347 26 : PetscCall(SlepcContourDataDestroy(&ctx->contour));
1348 26 : PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
1349 26 : PetscCall(PetscFree(eps->data));
1350 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",NULL));
1351 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",NULL));
1352 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",NULL));
1353 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",NULL));
1354 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",NULL));
1355 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",NULL));
1356 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetUseST_C",NULL));
1357 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetUseST_C",NULL));
1358 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetQuadRule_C",NULL));
1359 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetQuadRule_C",NULL));
1360 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetExtraction_C",NULL));
1361 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetExtraction_C",NULL));
1362 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetKSPs_C",NULL));
1363 26 : PetscFunctionReturn(PETSC_SUCCESS);
1364 : }
1365 :
1366 0 : static PetscErrorCode EPSView_CISS(EPS eps,PetscViewer viewer)
1367 : {
1368 0 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1369 0 : PetscBool isascii;
1370 0 : PetscViewer sviewer;
1371 :
1372 0 : PetscFunctionBegin;
1373 0 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
1374 0 : if (isascii) {
1375 0 : 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));
1376 0 : if (ctx->isreal) PetscCall(PetscViewerASCIIPrintf(viewer," exploiting symmetry of integration points\n"));
1377 0 : PetscCall(PetscViewerASCIIPrintf(viewer," threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold));
1378 0 : PetscCall(PetscViewerASCIIPrintf(viewer," iterative refinement { inner: %" PetscInt_FMT ", blocksize: %" PetscInt_FMT " }\n",ctx->refine_inner, ctx->refine_blocksize));
1379 0 : PetscCall(PetscViewerASCIIPrintf(viewer," extraction: %s\n",EPSCISSExtractions[ctx->extraction]));
1380 0 : PetscCall(PetscViewerASCIIPrintf(viewer," quadrature rule: %s\n",EPSCISSQuadRules[ctx->quad]));
1381 0 : if (ctx->usest) PetscCall(PetscViewerASCIIPrintf(viewer," using ST for linear solves\n"));
1382 : else {
1383 0 : if (!ctx->contour || !ctx->contour->ksp) PetscCall(EPSCISSGetKSPs(eps,NULL,NULL));
1384 0 : PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Something went wrong with EPSCISSGetKSPs()");
1385 0 : PetscCall(PetscViewerASCIIPushTab(viewer));
1386 0 : if (ctx->npart>1 && ctx->contour->subcomm) {
1387 0 : PetscCall(PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
1388 0 : if (!ctx->contour->subcomm->color) PetscCall(KSPView(ctx->contour->ksp[0],sviewer));
1389 0 : PetscCall(PetscViewerFlush(sviewer));
1390 0 : PetscCall(PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
1391 : /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1392 0 : PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1393 0 : } else PetscCall(KSPView(ctx->contour->ksp[0],viewer));
1394 0 : PetscCall(PetscViewerASCIIPopTab(viewer));
1395 : }
1396 : }
1397 0 : PetscFunctionReturn(PETSC_SUCCESS);
1398 : }
1399 :
1400 52 : static PetscErrorCode EPSSetDefaultST_CISS(EPS eps)
1401 : {
1402 52 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1403 52 : PetscBool usest = ctx->usest;
1404 52 : KSP ksp;
1405 52 : PC pc;
1406 :
1407 52 : PetscFunctionBegin;
1408 52 : if (!((PetscObject)eps->st)->type_name) {
1409 21 : if (!ctx->usest_set) usest = (ctx->npart>1)? PETSC_FALSE: PETSC_TRUE;
1410 21 : if (usest) PetscCall(STSetType(eps->st,STSINVERT));
1411 : else {
1412 : /* we are not going to use ST, so avoid factorizing the matrix */
1413 14 : PetscCall(STSetType(eps->st,STSHIFT));
1414 14 : if (eps->isgeneralized) {
1415 3 : PetscCall(STGetKSP(eps->st,&ksp));
1416 3 : PetscCall(KSPGetPC(ksp,&pc));
1417 3 : PetscCall(PCSetType(pc,PCNONE));
1418 : }
1419 : }
1420 : }
1421 52 : PetscFunctionReturn(PETSC_SUCCESS);
1422 : }
1423 :
1424 26 : SLEPC_EXTERN PetscErrorCode EPSCreate_CISS(EPS eps)
1425 : {
1426 26 : EPS_CISS *ctx = (EPS_CISS*)eps->data;
1427 :
1428 26 : PetscFunctionBegin;
1429 26 : PetscCall(PetscNew(&ctx));
1430 26 : eps->data = ctx;
1431 :
1432 26 : eps->useds = PETSC_TRUE;
1433 26 : eps->categ = EPS_CATEGORY_CONTOUR;
1434 :
1435 26 : eps->ops->solve = EPSSolve_CISS;
1436 26 : eps->ops->setup = EPSSetUp_CISS;
1437 26 : eps->ops->setupsort = EPSSetUpSort_CISS;
1438 26 : eps->ops->setfromoptions = EPSSetFromOptions_CISS;
1439 26 : eps->ops->destroy = EPSDestroy_CISS;
1440 26 : eps->ops->reset = EPSReset_CISS;
1441 26 : eps->ops->view = EPSView_CISS;
1442 26 : eps->ops->computevectors = EPSComputeVectors_CISS;
1443 26 : eps->ops->setdefaultst = EPSSetDefaultST_CISS;
1444 :
1445 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",EPSCISSSetSizes_CISS));
1446 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",EPSCISSGetSizes_CISS));
1447 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",EPSCISSSetThreshold_CISS));
1448 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",EPSCISSGetThreshold_CISS));
1449 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",EPSCISSSetRefinement_CISS));
1450 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",EPSCISSGetRefinement_CISS));
1451 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetUseST_C",EPSCISSSetUseST_CISS));
1452 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetUseST_C",EPSCISSGetUseST_CISS));
1453 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetQuadRule_C",EPSCISSSetQuadRule_CISS));
1454 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetQuadRule_C",EPSCISSGetQuadRule_CISS));
1455 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetExtraction_C",EPSCISSSetExtraction_CISS));
1456 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetExtraction_C",EPSCISSGetExtraction_CISS));
1457 26 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetKSPs_C",EPSCISSGetKSPs_CISS));
1458 :
1459 : /* set default values of parameters */
1460 26 : ctx->N = 32;
1461 26 : ctx->L = 16;
1462 26 : ctx->M = ctx->N/4;
1463 26 : ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
1464 26 : ctx->L_max = 64;
1465 26 : ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1466 26 : ctx->usest = PETSC_TRUE;
1467 26 : ctx->usest_set = PETSC_FALSE;
1468 26 : ctx->isreal = PETSC_FALSE;
1469 26 : ctx->refine_inner = 0;
1470 26 : ctx->refine_blocksize = 0;
1471 26 : ctx->npart = 1;
1472 26 : ctx->quad = (EPSCISSQuadRule)0;
1473 26 : ctx->extraction = EPS_CISS_EXTRACTION_RITZ;
1474 26 : PetscFunctionReturn(PETSC_SUCCESS);
1475 : }
|