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