Actual source code: ciss.c

slepc-main 2024-12-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc eigensolver: "ciss"

 13:    Method: Contour Integral Spectral Slicing

 15:    Algorithm:

 17:        Contour integral based on Sakurai-Sugiura method to construct a
 18:        subspace, with various eigenpair extractions (Rayleigh-Ritz,
 19:        explicit moment).

 21:    Based on code contributed by Y. Maeda, T. Sakurai.

 23:    References:

 25:        [1] T. Sakurai and H. Sugiura, "A projection method for generalized
 26:            eigenvalue problems", J. Comput. Appl. Math. 159:119-128, 2003.

 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: */

 33: #include <slepc/private/epsimpl.h>
 34: #include <slepc/private/slepccontour.h>
 35: #include <slepcblaslapack.h>

 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;

 68: /*
 69:   Set up KSP solvers for every integration point, only called if !ctx->usest
 70: */
 71: static PetscErrorCode EPSCISSSetUp(EPS eps,Mat A,Mat B,Mat Pa,Mat Pb)
 72: {
 73:   EPS_CISS         *ctx = (EPS_CISS*)eps->data;
 74:   SlepcContourData contour;
 75:   PetscInt         i,p_id,nsplit;
 76:   Mat              Amat,Pmat;
 77:   MatStructure     str,strp;

 79:   PetscFunctionBegin;
 80:   if (!ctx->contour || !ctx->contour->ksp) PetscCall(EPSCISSGetKSPs(eps,NULL,NULL));
 81:   PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Something went wrong with EPSCISSGetKSPs()");
 82:   contour = ctx->contour;
 83:   PetscCall(STGetMatStructure(eps->st,&str));
 84:   PetscCall(STGetSplitPreconditionerInfo(eps->st,&nsplit,&strp));
 85:   for (i=0;i<contour->npoints;i++) {
 86:     p_id = i*contour->subcomm->n + contour->subcomm->color;
 87:     PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&Amat));
 88:     if (B) PetscCall(MatAXPY(Amat,-ctx->omega[p_id],B,str));
 89:     else PetscCall(MatShift(Amat,-ctx->omega[p_id]));
 90:     if (nsplit) {
 91:       PetscCall(MatDuplicate(Pa,MAT_COPY_VALUES,&Pmat));
 92:       if (Pb) PetscCall(MatAXPY(Pmat,-ctx->omega[p_id],Pb,strp));
 93:       else PetscCall(MatShift(Pmat,-ctx->omega[p_id]));
 94:     } else Pmat = Amat;
 95:     PetscCall(EPS_KSPSetOperators(contour->ksp[i],Amat,Amat));
 96:     PetscCall(MatDestroy(&Amat));
 97:     if (nsplit) PetscCall(MatDestroy(&Pmat));
 98:   }
 99:   PetscFunctionReturn(PETSC_SUCCESS);
100: }

102: /*
103:   Y_i = (A-z_i B)^{-1}BV for every integration point, Y=[Y_i] is in the context
104: */
105: static PetscErrorCode EPSCISSSolve(EPS eps,Mat B,BV V,PetscInt L_start,PetscInt L_end)
106: {
107:   EPS_CISS         *ctx = (EPS_CISS*)eps->data;
108:   SlepcContourData contour;
109:   PetscInt         i,p_id;
110:   Mat              MV,BMV=NULL,MC;
111:   KSP              ksp;

113:   PetscFunctionBegin;
114:   if (!ctx->contour || !ctx->contour->ksp) PetscCall(EPSCISSGetKSPs(eps,NULL,NULL));
115:   contour = ctx->contour;
116:   PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Something went wrong with EPSCISSGetKSPs()");
117:   PetscCall(BVSetActiveColumns(V,L_start,L_end));
118:   PetscCall(BVGetMat(V,&MV));
119:   for (i=0;i<contour->npoints;i++) {
120:     p_id = i*contour->subcomm->n + contour->subcomm->color;
121:     if (ctx->usest)  {
122:       PetscCall(STSetShift(eps->st,ctx->omega[p_id]));
123:       PetscCall(STGetKSP(eps->st,&ksp));
124:     } else ksp = contour->ksp[i];
125:     PetscCall(BVSetActiveColumns(ctx->Y,i*ctx->L+L_start,i*ctx->L+L_end));
126:     PetscCall(BVGetMat(ctx->Y,&MC));
127:     if (B) {
128:       if (!i) {
129:         PetscCall(MatProductCreate(B,MV,NULL,&BMV));
130:         PetscCall(MatProductSetType(BMV,MATPRODUCT_AB));
131:         PetscCall(MatProductSetFromOptions(BMV));
132:         PetscCall(MatProductSymbolic(BMV));
133:       }
134:       PetscCall(MatProductNumeric(BMV));
135:       PetscCall(KSPMatSolve(ksp,BMV,MC));
136:     } else PetscCall(KSPMatSolve(ksp,MV,MC));
137:     PetscCall(BVRestoreMat(ctx->Y,&MC));
138:     if (ctx->usest && i<contour->npoints-1) PetscCall(KSPReset(ksp));
139:   }
140:   PetscCall(MatDestroy(&BMV));
141:   PetscCall(BVRestoreMat(V,&MV));
142:   PetscFunctionReturn(PETSC_SUCCESS);
143: }

145: static PetscErrorCode rescale_eig(EPS eps,PetscInt nv)
146: {
147:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
148:   PetscInt       i;
149:   PetscScalar    center;
150:   PetscReal      radius,a,b,c,d,rgscale;
151: #if defined(PETSC_USE_COMPLEX)
152:   PetscReal      start_ang,end_ang,vscale,theta;
153: #endif
154:   PetscBool      isring,isellipse,isinterval;

156:   PetscFunctionBegin;
157:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse));
158:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring));
159:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval));
160:   PetscCall(RGGetScale(eps->rg,&rgscale));
161:   if (isinterval) {
162:     PetscCall(RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d));
163:     if (c==d) {
164:       for (i=0;i<nv;i++) {
165: #if defined(PETSC_USE_COMPLEX)
166:         eps->eigr[i] = PetscRealPart(eps->eigr[i]);
167: #else
168:         eps->eigi[i] = 0;
169: #endif
170:       }
171:     }
172:   }
173:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
174:     if (isellipse) {
175:       PetscCall(RGEllipseGetParameters(eps->rg,&center,&radius,NULL));
176:       for (i=0;i<nv;i++) eps->eigr[i] = rgscale*(center + radius*eps->eigr[i]);
177:     } else if (isinterval) {
178:       PetscCall(RGIntervalGetEndpoints(eps->rg,&a,&b,&c,&d));
179:       if (ctx->quad == EPS_CISS_QUADRULE_CHEBYSHEV) {
180:         for (i=0;i<nv;i++) {
181:           if (c==d) eps->eigr[i] = ((b-a)*(eps->eigr[i]+1.0)/2.0+a)*rgscale;
182:           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:             SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Integration points on a vertical line require complex arithmetic");
187: #endif
188:           }
189:         }
190:       } else {
191:         center = (b+a)/2.0+(d+c)/2.0*PETSC_PI;
192:         radius = PetscSqrtReal(PetscPowRealInt((b-a)/2.0,2)+PetscPowRealInt((d-c)/2.0,2));
193:         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,&center,&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:     }
208:   }
209:   PetscFunctionReturn(PETSC_SUCCESS);
210: }

212: static PetscErrorCode EPSSetUp_CISS(EPS eps)
213: {
214:   EPS_CISS         *ctx = (EPS_CISS*)eps->data;
215:   SlepcContourData contour;
216:   PetscBool        istrivial,isring,isellipse,isinterval,flg;
217:   PetscReal        c,d;
218:   PetscInt         nsplit;
219:   PetscRandom      rand;
220:   PetscObjectId    id;
221:   PetscObjectState state;
222:   Mat              A[2],Psplit[2];
223:   Vec              v0;

225:   PetscFunctionBegin;
226:   EPSCheckNotStructured(eps);
227:   if (eps->ncv==PETSC_DETERMINE) {
228:     if (eps->nev==0) eps->nev = 1;
229:     eps->ncv = ctx->L_max*ctx->M;
230:     if (eps->ncv>eps->n) {
231:       eps->ncv = eps->n;
232:       ctx->L_max = eps->ncv/ctx->M;
233:       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:     PetscCall(EPSSetDimensions_Default(eps,&eps->nev,&eps->ncv,&eps->mpd));
237:     ctx->L_max = eps->ncv/ctx->M;
238:     if (!ctx->L_max) {
239:       ctx->L_max = 1;
240:       eps->ncv = ctx->L_max*ctx->M;
241:     }
242:   }
243:   ctx->L = PetscMin(ctx->L,ctx->L_max);
244:   if (eps->max_it==PETSC_DETERMINE) eps->max_it = 5;
245:   if (eps->mpd==PETSC_DETERMINE) eps->mpd = eps->ncv;
246:   if (!eps->which) eps->which = EPS_ALL;
247:   PetscCheck(eps->which==EPS_ALL,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only computing all eigenvalues");
248:   EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_EXTRACTION | EPS_FEATURE_STOPPING | EPS_FEATURE_TWOSIDED);

250:   /* check region */
251:   PetscCall(RGIsTrivial(eps->rg,&istrivial));
252:   PetscCheck(!istrivial,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
253:   PetscCall(RGGetComplement(eps->rg,&flg));
254:   PetscCheck(!flg,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
255:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse));
256:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring));
257:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval));
258:   PetscCheck(isellipse || isring || isinterval,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Currently only implemented for interval, elliptic or ring regions");

260:   /* if the region has changed, then reset contour data */
261:   PetscCall(PetscObjectGetId((PetscObject)eps->rg,&id));
262:   PetscCall(PetscObjectStateGet((PetscObject)eps->rg,&state));
263:   if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
264:     PetscCall(SlepcContourDataDestroy(&ctx->contour));
265:     PetscCall(PetscInfo(eps,"Resetting the contour data structure due to a change of region\n"));
266:     ctx->rgid = id; ctx->rgstate = state;
267:   }

269: #if !defined(PETSC_USE_COMPLEX)
270:   PetscCheck(!isring,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Ring region only supported for complex scalars");
271: #endif
272:   if (isinterval) {
273:     PetscCall(RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d));
274: #if !defined(PETSC_USE_COMPLEX)
275:     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:     if (!ctx->quad && c==d) ctx->quad = EPS_CISS_QUADRULE_CHEBYSHEV;
278:   }
279:   if (!ctx->quad) ctx->quad = EPS_CISS_QUADRULE_TRAPEZOIDAL;

281:   /* create contour data structure */
282:   if (!ctx->contour) {
283:     PetscCall(RGCanUseConjugates(eps->rg,ctx->isreal,&ctx->useconj));
284:     PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)eps,&ctx->contour));
285:   }

287:   PetscCall(EPSAllocateSolution(eps,0));
288:   PetscCall(BVGetRandomContext(eps->V,&rand));  /* make sure the random context is available when duplicating */
289:   if (ctx->weight) PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
290:   PetscCall(PetscMalloc4(ctx->N,&ctx->weight,ctx->N+1,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma));

292:   /* allocate basis vectors */
293:   PetscCall(BVDestroy(&ctx->S));
294:   PetscCall(BVDuplicateResize(eps->V,ctx->L*ctx->M,&ctx->S));
295:   PetscCall(BVDestroy(&ctx->V));
296:   PetscCall(BVDuplicateResize(eps->V,ctx->L,&ctx->V));

298:   PetscCall(STGetMatrix(eps->st,0,&A[0]));
299:   PetscCall(MatIsShell(A[0],&flg));
300:   PetscCheck(!flg,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Matrix type shell is not supported in this solver");
301:   if (eps->isgeneralized) PetscCall(STGetMatrix(eps->st,1,&A[1]));

303:   if (!ctx->usest_set) ctx->usest = (ctx->npart>1)? PETSC_FALSE: PETSC_TRUE;
304:   PetscCheck(!ctx->usest || ctx->npart==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The usest flag is not supported when partitions > 1");

306:   /* check if a user-defined split preconditioner has been set */
307:   PetscCall(STGetSplitPreconditionerInfo(eps->st,&nsplit,NULL));
308:   if (nsplit) {
309:     PetscCall(STGetSplitPreconditionerTerm(eps->st,0,&Psplit[0]));
310:     if (eps->isgeneralized) PetscCall(STGetSplitPreconditionerTerm(eps->st,1,&Psplit[1]));
311:   }

313:   contour = ctx->contour;
314:   PetscCall(SlepcContourRedundantMat(contour,eps->isgeneralized?2:1,A,nsplit?Psplit:NULL));
315:   if (contour->pA) {
316:     PetscCall(BVGetColumn(ctx->V,0,&v0));
317:     PetscCall(SlepcContourScatterCreate(contour,v0));
318:     PetscCall(BVRestoreColumn(ctx->V,0,&v0));
319:     PetscCall(BVDestroy(&ctx->pV));
320:     PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV));
321:     PetscCall(BVSetSizesFromVec(ctx->pV,contour->xsub,eps->n));
322:     PetscCall(BVSetFromOptions(ctx->pV));
323:     PetscCall(BVResize(ctx->pV,ctx->L,PETSC_FALSE));
324:   }

326:   EPSCheckDefinite(eps);
327:   EPSCheckSinvertCondition(eps,ctx->usest," (with the usest flag set)");

329:   PetscCall(BVDestroy(&ctx->Y));
330:   if (contour->pA) {
331:     PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y));
332:     PetscCall(BVSetSizesFromVec(ctx->Y,contour->xsub,eps->n));
333:     PetscCall(BVSetFromOptions(ctx->Y));
334:     PetscCall(BVResize(ctx->Y,contour->npoints*ctx->L,PETSC_FALSE));
335:   } else PetscCall(BVDuplicateResize(eps->V,contour->npoints*ctx->L,&ctx->Y));

337:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) PetscCall(DSSetType(eps->ds,DSGNHEP));
338:   else if (eps->isgeneralized) {
339:     if (eps->ishermitian && eps->ispositive) PetscCall(DSSetType(eps->ds,DSGHEP));
340:     else PetscCall(DSSetType(eps->ds,DSGNHEP));
341:   } else {
342:     if (eps->ishermitian) PetscCall(DSSetType(eps->ds,DSHEP));
343:     else PetscCall(DSSetType(eps->ds,DSNHEP));
344:   }
345:   PetscCall(DSAllocate(eps->ds,eps->ncv));

347: #if !defined(PETSC_USE_COMPLEX)
348:   PetscCall(EPSSetWorkVecs(eps,3));
349:   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:   PetscFunctionReturn(PETSC_SUCCESS);
354: }

356: static PetscErrorCode EPSSetUpSort_CISS(EPS eps)
357: {
358:   SlepcSC        sc;

360:   PetscFunctionBegin;
361:   /* fill sorting criterion context */
362:   eps->sc->comparison    = SlepcCompareSmallestReal;
363:   eps->sc->comparisonctx = NULL;
364:   eps->sc->map           = NULL;
365:   eps->sc->mapobj        = NULL;

367:   /* fill sorting criterion for DS */
368:   PetscCall(DSGetSlepcSC(eps->ds,&sc));
369:   sc->comparison    = SlepcCompareLargestMagnitude;
370:   sc->comparisonctx = NULL;
371:   sc->map           = NULL;
372:   sc->mapobj        = NULL;
373:   PetscFunctionReturn(PETSC_SUCCESS);
374: }

376: static PetscErrorCode EPSSolve_CISS(EPS eps)
377: {
378:   EPS_CISS         *ctx = (EPS_CISS*)eps->data;
379:   SlepcContourData contour = ctx->contour;
380:   Mat              A,B,X,M,pA,pB,T,J,Pa=NULL,Pb=NULL;
381:   BV               V;
382:   PetscInt         i,j,ld,nmat,L_add=0,nv=0,L_base=ctx->L,inner,nlocal,*inside,nsplit;
383:   PetscScalar      *Mu,*H0,*H1=NULL,*rr,*temp;
384:   PetscReal        error,max_error,norm;
385:   PetscBool        *fl1;
386:   Vec              si,si1=NULL,w[3];
387:   PetscRandom      rand;
388: #if defined(PETSC_USE_COMPLEX)
389:   PetscBool        isellipse;
390:   PetscReal        est_eig,eta;
391: #else
392:   PetscReal        normi;
393: #endif

395:   PetscFunctionBegin;
396:   w[0] = eps->work[0];
397: #if defined(PETSC_USE_COMPLEX)
398:   w[1] = NULL;
399: #else
400:   w[1] = eps->work[2];
401: #endif
402:   w[2] = eps->work[1];
403:   PetscCall(VecGetLocalSize(w[0],&nlocal));
404:   PetscCall(DSGetLeadingDimension(eps->ds,&ld));
405:   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:   PetscCall(STGetNumMatrices(eps->st,&nmat));
407:   PetscCall(STGetMatrix(eps->st,0,&A));
408:   if (nmat>1) PetscCall(STGetMatrix(eps->st,1,&B));
409:   else B = NULL;
410:   J = (contour->pA && nmat>1)? contour->pA[1]: B;
411:   V = contour->pA? ctx->pV: ctx->V;
412:   if (!ctx->usest) {
413:     T = contour->pA? contour->pA[0]: A;
414:     PetscCall(STGetSplitPreconditionerInfo(eps->st,&nsplit,NULL));
415:     if (nsplit) {
416:       if (contour->pA) {
417:         Pa = contour->pP[0];
418:         if (nsplit>1) Pb = contour->pP[1];
419:       } else {
420:         PetscCall(STGetSplitPreconditionerTerm(eps->st,0,&Pa));
421:         if (nsplit>1) PetscCall(STGetSplitPreconditionerTerm(eps->st,1,&Pb));
422:       }
423:     }
424:     PetscCall(EPSCISSSetUp(eps,T,J,Pa,Pb));
425:   }
426:   PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
427:   PetscCall(BVSetRandomSign(ctx->V));
428:   PetscCall(BVGetRandomContext(ctx->V,&rand));

430:   if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
431:   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:   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:     PetscCall(EPSCISSSolve(eps,J,V,ctx->L-L_add,ctx->L));
453:   }
454:   PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
455:   for (i=0;i<ctx->refine_blocksize;i++) {
456:     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:     PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
458:     PetscCall(PetscLogEventBegin(EPS_CISS_SVD,eps,0,0,0));
459:     PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
460:     PetscCall(PetscLogEventEnd(EPS_CISS_SVD,eps,0,0,0));
461:     if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
462:     L_add = L_base;
463:     if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
464:     PetscCall(PetscInfo(eps,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by SVD(H0)\n",ctx->L,ctx->L+L_add));
465:     PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
466:     PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
467:     PetscCall(BVSetRandomSign(ctx->V));
468:     if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
469:     ctx->L += L_add;
470:     PetscCall(EPSCISSSolve(eps,J,V,ctx->L-L_add,ctx->L));
471:     if (L_add) {
472:       PetscCall(PetscFree2(Mu,H0));
473:       PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
474:     }
475:   }
476:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) PetscCall(PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1));

478:   while (eps->reason == EPS_CONVERGED_ITERATING) {
479:     eps->its++;
480:     for (inner=0;inner<=ctx->refine_inner;inner++) {
481:       if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
482:         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:         PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
484:         PetscCall(PetscLogEventBegin(EPS_CISS_SVD,eps,0,0,0));
485:         PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
486:         PetscCall(PetscLogEventEnd(EPS_CISS_SVD,eps,0,0,0));
487:         break;
488:       } else {
489:         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:         PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
491:         PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
492:         PetscCall(BVCopy(ctx->S,ctx->V));
493:         PetscCall(BVSVDAndRank(ctx->S,ctx->M,ctx->L,ctx->delta,BV_SVD_METHOD_REFINE,H0,ctx->sigma,&nv));
494:         if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
495:           if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
496:           PetscCall(EPSCISSSolve(eps,J,V,0,ctx->L));
497:         } else break;
498:       }
499:     }
500:     eps->nconv = 0;
501:     if (nv == 0) eps->reason = EPS_CONVERGED_TOL;
502:     else {
503:       PetscCall(DSSetDimensions(eps->ds,nv,0,0));
504:       PetscCall(DSSetState(eps->ds,DS_STATE_RAW));

506:       if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
507:         PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
508:         PetscCall(CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1));
509:         PetscCall(DSGetArray(eps->ds,DS_MAT_A,&temp));
510:         for (j=0;j<nv;j++) {
511:           for (i=0;i<nv;i++) {
512:             temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
513:           }
514:         }
515:         PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&temp));
516:         PetscCall(DSGetArray(eps->ds,DS_MAT_B,&temp));
517:         for (j=0;j<nv;j++) {
518:           for (i=0;i<nv;i++) {
519:             temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
520:           }
521:         }
522:         PetscCall(DSRestoreArray(eps->ds,DS_MAT_B,&temp));
523:       } else {
524:         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
525:         PetscCall(DSGetMat(eps->ds,DS_MAT_A,&pA));
526:         PetscCall(MatZeroEntries(pA));
527:         PetscCall(BVMatProject(ctx->S,A,ctx->S,pA));
528:         PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&pA));
529:         if (B) {
530:           PetscCall(DSGetMat(eps->ds,DS_MAT_B,&pB));
531:           PetscCall(MatZeroEntries(pB));
532:           PetscCall(BVMatProject(ctx->S,B,ctx->S,pB));
533:           PetscCall(DSRestoreMat(eps->ds,DS_MAT_B,&pB));
534:         }
535:       }

537:       PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi));
538:       PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));

540:       PetscCall(PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr));
541:       PetscCall(rescale_eig(eps,nv));
542:       PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
543:       PetscCall(DSGetMat(eps->ds,DS_MAT_X,&X));
544:       PetscCall(SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1));
545:       PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&X));
546:       PetscCall(RGCheckInside(eps->rg,nv,eps->eigr,eps->eigi,inside));
547:       for (i=0;i<nv;i++) {
548:         if (fl1[i] && inside[i]>=0) {
549:           rr[i] = 1.0;
550:           eps->nconv++;
551:         } else rr[i] = 0.0;
552:       }
553:       PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,rr,NULL,&eps->nconv));
554:       PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));
555:       PetscCall(rescale_eig(eps,nv));
556:       PetscCall(PetscFree3(fl1,inside,rr));
557:       PetscCall(BVSetActiveColumns(eps->V,0,nv));
558:       if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
559:         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:         PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
561:         PetscCall(BVCopy(ctx->S,ctx->V));
562:         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
563:       }
564:       PetscCall(BVCopy(ctx->S,eps->V));

566:       PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
567:       PetscCall(DSGetMat(eps->ds,DS_MAT_X,&X));
568:       PetscCall(BVMultInPlace(ctx->S,X,0,eps->nconv));
569:       if (eps->ishermitian) PetscCall(BVMultInPlace(eps->V,X,0,eps->nconv));
570:       PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&X));
571:       max_error = 0.0;
572:       for (i=0;i<eps->nconv;i++) {
573:         PetscCall(BVGetColumn(ctx->S,i,&si));
574: #if !defined(PETSC_USE_COMPLEX)
575:         if (eps->eigi[i]!=0.0) PetscCall(BVGetColumn(ctx->S,i+1,&si1));
576: #endif
577:         PetscCall(EPSComputeResidualNorm_Private(eps,PETSC_FALSE,eps->eigr[i],eps->eigi[i],si,si1,w,&error));
578:         if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {  /* vector is not normalized */
579:           PetscCall(VecNorm(si,NORM_2,&norm));
580: #if !defined(PETSC_USE_COMPLEX)
581:           if (eps->eigi[i]!=0.0) {
582:             PetscCall(VecNorm(si1,NORM_2,&normi));
583:             norm = SlepcAbsEigenvalue(norm,normi);
584:           }
585: #endif
586:           error /= norm;
587:         }
588:         PetscCall((*eps->converged)(eps,eps->eigr[i],eps->eigi[i],error,&error,eps->convergedctx));
589:         PetscCall(BVRestoreColumn(ctx->S,i,&si));
590: #if !defined(PETSC_USE_COMPLEX)
591:         if (eps->eigi[i]!=0.0) {
592:           PetscCall(BVRestoreColumn(ctx->S,i+1,&si1));
593:           i++;
594:         }
595: #endif
596:         max_error = PetscMax(max_error,error);
597:       }

599:       if (max_error <= eps->tol) eps->reason = EPS_CONVERGED_TOL;
600:       else if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
601:       else {
602:         if (eps->nconv > ctx->L) nv = eps->nconv;
603:         else if (ctx->L > nv) nv = ctx->L;
604:         nv = PetscMin(nv,ctx->L*ctx->M);
605:         PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M));
606:         PetscCall(MatSetRandom(M,rand));
607:         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
608:         PetscCall(BVMultInPlace(ctx->S,M,0,ctx->L));
609:         PetscCall(MatDestroy(&M));
610:         PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
611:         PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
612:         PetscCall(BVCopy(ctx->S,ctx->V));
613:         if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
614:         PetscCall(EPSCISSSolve(eps,J,V,0,ctx->L));
615:       }
616:     }
617:   }
618:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) PetscCall(PetscFree(H1));
619:   PetscCall(PetscFree2(Mu,H0));
620:   PetscFunctionReturn(PETSC_SUCCESS);
621: }

623: static PetscErrorCode EPSComputeVectors_CISS(EPS eps)
624: {
625:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
626:   PetscInt       n;
627:   Mat            Z,B=NULL;

629:   PetscFunctionBegin;
630:   if (eps->ishermitian) {
631:     if (eps->isgeneralized && !eps->ispositive) PetscCall(EPSComputeVectors_Indefinite(eps));
632:     else PetscCall(EPSComputeVectors_Hermitian(eps));
633:     if (eps->isgeneralized && eps->ispositive && ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
634:       /* normalize to have unit B-norm */
635:       PetscCall(STGetMatrix(eps->st,1,&B));
636:       PetscCall(BVSetMatrix(eps->V,B,PETSC_FALSE));
637:       PetscCall(BVNormalize(eps->V,NULL));
638:       PetscCall(BVSetMatrix(eps->V,NULL,PETSC_FALSE));
639:     }
640:     PetscFunctionReturn(PETSC_SUCCESS);
641:   }
642:   PetscCall(DSGetDimensions(eps->ds,&n,NULL,NULL,NULL));
643:   PetscCall(BVSetActiveColumns(eps->V,0,n));

645:   /* right eigenvectors */
646:   PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));

648:   /* V = V * Z */
649:   PetscCall(DSGetMat(eps->ds,DS_MAT_X,&Z));
650:   PetscCall(BVMultInPlace(eps->V,Z,0,n));
651:   PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&Z));
652:   PetscCall(BVSetActiveColumns(eps->V,0,eps->nconv));

654:   /* normalize */
655:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) PetscCall(BVNormalize(eps->V,NULL));
656:   PetscFunctionReturn(PETSC_SUCCESS);
657: }

659: static PetscErrorCode EPSCISSSetSizes_CISS(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
660: {
661:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
662:   PetscInt       oN,oL,oM,oLmax,onpart;
663:   PetscMPIInt    size;

665:   PetscFunctionBegin;
666:   oN = ctx->N;
667:   if (ip == PETSC_DETERMINE) {
668:     if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
669:   } else if (ip != PETSC_CURRENT) {
670:     PetscCheck(ip>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
671:     PetscCheck(ip%2==0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
672:     if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
673:   }
674:   oL = ctx->L;
675:   if (bs == PETSC_DETERMINE) {
676:     ctx->L = 16;
677:   } else if (bs != PETSC_CURRENT) {
678:     PetscCheck(bs>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
679:     ctx->L = bs;
680:   }
681:   oM = ctx->M;
682:   if (ms == PETSC_DETERMINE) {
683:     ctx->M = ctx->N/4;
684:   } else if (ms != PETSC_CURRENT) {
685:     PetscCheck(ms>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
686:     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:     ctx->M = ms;
688:   }
689:   onpart = ctx->npart;
690:   if (npart == PETSC_DETERMINE) {
691:     ctx->npart = 1;
692:   } else if (npart != PETSC_CURRENT) {
693:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size));
694:     PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
695:     ctx->npart = npart;
696:   }
697:   oLmax = ctx->L_max;
698:   if (bsmax == PETSC_DETERMINE) {
699:     ctx->L_max = 64;
700:   } else if (bsmax != PETSC_CURRENT) {
701:     PetscCheck(bsmax>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
702:     ctx->L_max = PetscMax(bsmax,ctx->L);
703:   }
704:   if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
705:     PetscCall(SlepcContourDataDestroy(&ctx->contour));
706:     PetscCall(PetscInfo(eps,"Resetting the contour data structure due to a change of parameters\n"));
707:     eps->state = EPS_STATE_INITIAL;
708:   }
709:   ctx->isreal = realmats;
710:   if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) eps->state = EPS_STATE_INITIAL;
711:   PetscFunctionReturn(PETSC_SUCCESS);
712: }

714: /*@
715:    EPSCISSSetSizes - Sets the values of various size parameters in the CISS solver.

717:    Logically Collective

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

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

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.

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.

744:    Level: advanced

746: .seealso: EPSCISSGetSizes()
747: @*/
748: PetscErrorCode EPSCISSSetSizes(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
749: {
750:   PetscFunctionBegin;
758:   PetscTryMethod(eps,"EPSCISSSetSizes_C",(EPS,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(eps,ip,bs,ms,npart,bsmax,realmats));
759:   PetscFunctionReturn(PETSC_SUCCESS);
760: }

762: static PetscErrorCode EPSCISSGetSizes_CISS(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
763: {
764:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

766:   PetscFunctionBegin;
767:   if (ip) *ip = ctx->N;
768:   if (bs) *bs = ctx->L;
769:   if (ms) *ms = ctx->M;
770:   if (npart) *npart = ctx->npart;
771:   if (bsmax) *bsmax = ctx->L_max;
772:   if (realmats) *realmats = ctx->isreal;
773:   PetscFunctionReturn(PETSC_SUCCESS);
774: }

776: /*@
777:    EPSCISSGetSizes - Gets the values of various size parameters in the CISS solver.

779:    Not Collective

781:    Input Parameter:
782: .  eps - the eigenproblem solver context

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

792:    Level: advanced

794: .seealso: EPSCISSSetSizes()
795: @*/
796: PetscErrorCode EPSCISSGetSizes(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
797: {
798:   PetscFunctionBegin;
800:   PetscUseMethod(eps,"EPSCISSGetSizes_C",(EPS,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(eps,ip,bs,ms,npart,bsmax,realmats));
801:   PetscFunctionReturn(PETSC_SUCCESS);
802: }

804: static PetscErrorCode EPSCISSSetThreshold_CISS(EPS eps,PetscReal delta,PetscReal spur)
805: {
806:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

808:   PetscFunctionBegin;
809:   if (delta == (PetscReal)PETSC_DETERMINE) {
810:     ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
811:   } else if (delta != (PetscReal)PETSC_CURRENT) {
812:     PetscCheck(delta>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
813:     ctx->delta = delta;
814:   }
815:   if (spur == (PetscReal)PETSC_DETERMINE) {
816:     ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
817:   } else if (spur != (PetscReal)PETSC_CURRENT) {
818:     PetscCheck(spur>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
819:     ctx->spurious_threshold = spur;
820:   }
821:   PetscFunctionReturn(PETSC_SUCCESS);
822: }

824: /*@
825:    EPSCISSSetThreshold - Sets the values of various threshold parameters in
826:    the CISS solver.

828:    Logically Collective

830:    Input Parameters:
831: +  eps   - the eigenproblem solver context
832: .  delta - threshold for numerical rank
833: -  spur  - spurious threshold (to discard spurious eigenpairs)

835:    Options Database Keys:
836: +  -eps_ciss_delta - Sets the delta
837: -  -eps_ciss_spurious_threshold - Sets the spurious threshold

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.

843:    Level: advanced

845: .seealso: EPSCISSGetThreshold()
846: @*/
847: PetscErrorCode EPSCISSSetThreshold(EPS eps,PetscReal delta,PetscReal spur)
848: {
849:   PetscFunctionBegin;
853:   PetscTryMethod(eps,"EPSCISSSetThreshold_C",(EPS,PetscReal,PetscReal),(eps,delta,spur));
854:   PetscFunctionReturn(PETSC_SUCCESS);
855: }

857: static PetscErrorCode EPSCISSGetThreshold_CISS(EPS eps,PetscReal *delta,PetscReal *spur)
858: {
859:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

861:   PetscFunctionBegin;
862:   if (delta) *delta = ctx->delta;
863:   if (spur)  *spur = ctx->spurious_threshold;
864:   PetscFunctionReturn(PETSC_SUCCESS);
865: }

867: /*@
868:    EPSCISSGetThreshold - Gets the values of various threshold parameters
869:    in the CISS solver.

871:    Not Collective

873:    Input Parameter:
874: .  eps - the eigenproblem solver context

876:    Output Parameters:
877: +  delta - threshold for numerical rank
878: -  spur  - spurious threshold (to discard spurious eigenpairs)

880:    Level: advanced

882: .seealso: EPSCISSSetThreshold()
883: @*/
884: PetscErrorCode EPSCISSGetThreshold(EPS eps,PetscReal *delta,PetscReal *spur)
885: {
886:   PetscFunctionBegin;
888:   PetscUseMethod(eps,"EPSCISSGetThreshold_C",(EPS,PetscReal*,PetscReal*),(eps,delta,spur));
889:   PetscFunctionReturn(PETSC_SUCCESS);
890: }

892: static PetscErrorCode EPSCISSSetRefinement_CISS(EPS eps,PetscInt inner,PetscInt blsize)
893: {
894:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

896:   PetscFunctionBegin;
897:   if (inner == PETSC_DETERMINE) {
898:     ctx->refine_inner = 0;
899:   } else if (inner != PETSC_CURRENT) {
900:     PetscCheck(inner>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
901:     ctx->refine_inner = inner;
902:   }
903:   if (blsize == PETSC_DETERMINE) {
904:     ctx->refine_blocksize = 0;
905:   } else if (blsize != PETSC_CURRENT) {
906:     PetscCheck(blsize>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
907:     ctx->refine_blocksize = blsize;
908:   }
909:   PetscFunctionReturn(PETSC_SUCCESS);
910: }

912: /*@
913:    EPSCISSSetRefinement - Sets the values of various refinement parameters
914:    in the CISS solver.

916:    Logically Collective

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)

923:    Options Database Keys:
924: +  -eps_ciss_refine_inner - Sets number of inner iterations
925: -  -eps_ciss_refine_blocksize - Sets number of blocksize iterations

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.

931:    Level: advanced

933: .seealso: EPSCISSGetRefinement()
934: @*/
935: PetscErrorCode EPSCISSSetRefinement(EPS eps,PetscInt inner,PetscInt blsize)
936: {
937:   PetscFunctionBegin;
941:   PetscTryMethod(eps,"EPSCISSSetRefinement_C",(EPS,PetscInt,PetscInt),(eps,inner,blsize));
942:   PetscFunctionReturn(PETSC_SUCCESS);
943: }

945: static PetscErrorCode EPSCISSGetRefinement_CISS(EPS eps,PetscInt *inner,PetscInt *blsize)
946: {
947:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

949:   PetscFunctionBegin;
950:   if (inner)  *inner = ctx->refine_inner;
951:   if (blsize) *blsize = ctx->refine_blocksize;
952:   PetscFunctionReturn(PETSC_SUCCESS);
953: }

955: /*@
956:    EPSCISSGetRefinement - Gets the values of various refinement parameters
957:    in the CISS solver.

959:    Not Collective

961:    Input Parameter:
962: .  eps - the eigenproblem solver context

964:    Output Parameters:
965: +  inner  - number of iterative refinement iterations (inner loop)
966: -  blsize - number of iterative refinement iterations (blocksize loop)

968:    Level: advanced

970: .seealso: EPSCISSSetRefinement()
971: @*/
972: PetscErrorCode EPSCISSGetRefinement(EPS eps, PetscInt *inner, PetscInt *blsize)
973: {
974:   PetscFunctionBegin;
976:   PetscUseMethod(eps,"EPSCISSGetRefinement_C",(EPS,PetscInt*,PetscInt*),(eps,inner,blsize));
977:   PetscFunctionReturn(PETSC_SUCCESS);
978: }

980: static PetscErrorCode EPSCISSSetUseST_CISS(EPS eps,PetscBool usest)
981: {
982:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

984:   PetscFunctionBegin;
985:   ctx->usest     = usest;
986:   ctx->usest_set = PETSC_TRUE;
987:   eps->state     = EPS_STATE_INITIAL;
988:   PetscFunctionReturn(PETSC_SUCCESS);
989: }

991: /*@
992:    EPSCISSSetUseST - Sets a flag indicating that the CISS solver will
993:    use the ST object for the linear solves.

995:    Logically Collective

997:    Input Parameters:
998: +  eps    - the eigenproblem solver context
999: -  usest  - boolean flag to use the ST object or not

1001:    Options Database Keys:
1002: .  -eps_ciss_usest <bool> - whether the ST object will be used or not

1004:    Level: advanced

1006: .seealso: EPSCISSGetUseST()
1007: @*/
1008: PetscErrorCode EPSCISSSetUseST(EPS eps,PetscBool usest)
1009: {
1010:   PetscFunctionBegin;
1013:   PetscTryMethod(eps,"EPSCISSSetUseST_C",(EPS,PetscBool),(eps,usest));
1014:   PetscFunctionReturn(PETSC_SUCCESS);
1015: }

1017: static PetscErrorCode EPSCISSGetUseST_CISS(EPS eps,PetscBool *usest)
1018: {
1019:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1021:   PetscFunctionBegin;
1022:   *usest = ctx->usest;
1023:   PetscFunctionReturn(PETSC_SUCCESS);
1024: }

1026: /*@
1027:    EPSCISSGetUseST - Gets the flag for using the ST object
1028:    in the CISS solver.

1030:    Not Collective

1032:    Input Parameter:
1033: .  eps - the eigenproblem solver context

1035:    Output Parameters:
1036: .  usest - boolean flag indicating if the ST object is being used

1038:    Level: advanced

1040: .seealso: EPSCISSSetUseST()
1041: @*/
1042: PetscErrorCode EPSCISSGetUseST(EPS eps,PetscBool *usest)
1043: {
1044:   PetscFunctionBegin;
1046:   PetscAssertPointer(usest,2);
1047:   PetscUseMethod(eps,"EPSCISSGetUseST_C",(EPS,PetscBool*),(eps,usest));
1048:   PetscFunctionReturn(PETSC_SUCCESS);
1049: }

1051: static PetscErrorCode EPSCISSSetQuadRule_CISS(EPS eps,EPSCISSQuadRule quad)
1052: {
1053:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1055:   PetscFunctionBegin;
1056:   if (ctx->quad != quad) {
1057:     ctx->quad  = quad;
1058:     eps->state = EPS_STATE_INITIAL;
1059:   }
1060:   PetscFunctionReturn(PETSC_SUCCESS);
1061: }

1063: /*@
1064:    EPSCISSSetQuadRule - Sets the quadrature rule used in the CISS solver.

1066:    Logically Collective

1068:    Input Parameters:
1069: +  eps  - the eigenproblem solver context
1070: -  quad - the quadrature rule

1072:    Options Database Key:
1073: .  -eps_ciss_quadrule - Sets the quadrature rule (either 'trapezoidal' or
1074:                            'chebyshev')

1076:    Notes:
1077:    By default, the trapezoidal rule is used (EPS_CISS_QUADRULE_TRAPEZOIDAL).

1079:    If the 'chebyshev' option is specified (EPS_CISS_QUADRULE_CHEBYSHEV), then
1080:    Chebyshev points are used as quadrature points.

1082:    Level: advanced

1084: .seealso: EPSCISSGetQuadRule(), EPSCISSQuadRule
1085: @*/
1086: PetscErrorCode EPSCISSSetQuadRule(EPS eps,EPSCISSQuadRule quad)
1087: {
1088:   PetscFunctionBegin;
1091:   PetscTryMethod(eps,"EPSCISSSetQuadRule_C",(EPS,EPSCISSQuadRule),(eps,quad));
1092:   PetscFunctionReturn(PETSC_SUCCESS);
1093: }

1095: static PetscErrorCode EPSCISSGetQuadRule_CISS(EPS eps,EPSCISSQuadRule *quad)
1096: {
1097:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1099:   PetscFunctionBegin;
1100:   *quad = ctx->quad;
1101:   PetscFunctionReturn(PETSC_SUCCESS);
1102: }

1104: /*@
1105:    EPSCISSGetQuadRule - Gets the quadrature rule used in the CISS solver.

1107:    Not Collective

1109:    Input Parameter:
1110: .  eps - the eigenproblem solver context

1112:    Output Parameters:
1113: .  quad - quadrature rule

1115:    Level: advanced

1117: .seealso: EPSCISSSetQuadRule() EPSCISSQuadRule
1118: @*/
1119: PetscErrorCode EPSCISSGetQuadRule(EPS eps,EPSCISSQuadRule *quad)
1120: {
1121:   PetscFunctionBegin;
1123:   PetscAssertPointer(quad,2);
1124:   PetscUseMethod(eps,"EPSCISSGetQuadRule_C",(EPS,EPSCISSQuadRule*),(eps,quad));
1125:   PetscFunctionReturn(PETSC_SUCCESS);
1126: }

1128: static PetscErrorCode EPSCISSSetExtraction_CISS(EPS eps,EPSCISSExtraction extraction)
1129: {
1130:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1132:   PetscFunctionBegin;
1133:   if (ctx->extraction != extraction) {
1134:     ctx->extraction = extraction;
1135:     eps->state      = EPS_STATE_INITIAL;
1136:   }
1137:   PetscFunctionReturn(PETSC_SUCCESS);
1138: }

1140: /*@
1141:    EPSCISSSetExtraction - Sets the extraction technique used in the CISS solver.

1143:    Logically Collective

1145:    Input Parameters:
1146: +  eps        - the eigenproblem solver context
1147: -  extraction - the extraction technique

1149:    Options Database Key:
1150: .  -eps_ciss_extraction - Sets the extraction technique (either 'ritz' or
1151:                            'hankel')

1153:    Notes:
1154:    By default, the Rayleigh-Ritz extraction is used (EPS_CISS_EXTRACTION_RITZ).

1156:    If the 'hankel' option is specified (EPS_CISS_EXTRACTION_HANKEL), then
1157:    the Block Hankel method is used for extracting eigenpairs.

1159:    Level: advanced

1161: .seealso: EPSCISSGetExtraction(), EPSCISSExtraction
1162: @*/
1163: PetscErrorCode EPSCISSSetExtraction(EPS eps,EPSCISSExtraction extraction)
1164: {
1165:   PetscFunctionBegin;
1168:   PetscTryMethod(eps,"EPSCISSSetExtraction_C",(EPS,EPSCISSExtraction),(eps,extraction));
1169:   PetscFunctionReturn(PETSC_SUCCESS);
1170: }

1172: static PetscErrorCode EPSCISSGetExtraction_CISS(EPS eps,EPSCISSExtraction *extraction)
1173: {
1174:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1176:   PetscFunctionBegin;
1177:   *extraction = ctx->extraction;
1178:   PetscFunctionReturn(PETSC_SUCCESS);
1179: }

1181: /*@
1182:    EPSCISSGetExtraction - Gets the extraction technique used in the CISS solver.

1184:    Not Collective

1186:    Input Parameter:
1187: .  eps - the eigenproblem solver context

1189:    Output Parameters:
1190: .  extraction - extraction technique

1192:    Level: advanced

1194: .seealso: EPSCISSSetExtraction() EPSCISSExtraction
1195: @*/
1196: PetscErrorCode EPSCISSGetExtraction(EPS eps,EPSCISSExtraction *extraction)
1197: {
1198:   PetscFunctionBegin;
1200:   PetscAssertPointer(extraction,2);
1201:   PetscUseMethod(eps,"EPSCISSGetExtraction_C",(EPS,EPSCISSExtraction*),(eps,extraction));
1202:   PetscFunctionReturn(PETSC_SUCCESS);
1203: }

1205: static PetscErrorCode EPSCISSGetKSPs_CISS(EPS eps,PetscInt *nsolve,KSP **ksp)
1206: {
1207:   EPS_CISS         *ctx = (EPS_CISS*)eps->data;
1208:   SlepcContourData contour;
1209:   PetscInt         i,nsplit;
1210:   PC               pc;
1211:   MPI_Comm         child;

1213:   PetscFunctionBegin;
1214:   if (!ctx->contour) {  /* initialize contour data structure first */
1215:     PetscCall(RGCanUseConjugates(eps->rg,ctx->isreal,&ctx->useconj));
1216:     PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)eps,&ctx->contour));
1217:   }
1218:   contour = ctx->contour;
1219:   if (!contour->ksp) {
1220:     PetscCall(PetscMalloc1(contour->npoints,&contour->ksp));
1221:     PetscCall(EPSGetST(eps,&eps->st));
1222:     PetscCall(STGetSplitPreconditionerInfo(eps->st,&nsplit,NULL));
1223:     PetscCall(PetscSubcommGetChild(contour->subcomm,&child));
1224:     for (i=0;i<contour->npoints;i++) {
1225:       PetscCall(KSPCreate(child,&contour->ksp[i]));
1226:       PetscCall(PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)eps,1));
1227:       PetscCall(KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)eps)->prefix));
1228:       PetscCall(KSPAppendOptionsPrefix(contour->ksp[i],"eps_ciss_"));
1229:       PetscCall(PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)eps)->options));
1230:       PetscCall(KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE));
1231:       PetscCall(KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(eps->tol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
1232:       PetscCall(KSPGetPC(contour->ksp[i],&pc));
1233:       if (nsplit) {
1234:         PetscCall(KSPSetType(contour->ksp[i],KSPBCGS));
1235:         PetscCall(PCSetType(pc,PCBJACOBI));
1236:       } else {
1237:         PetscCall(KSPSetType(contour->ksp[i],KSPPREONLY));
1238:         PetscCall(PCSetType(pc,PCLU));
1239:       }
1240:     }
1241:   }
1242:   if (nsolve) *nsolve = contour->npoints;
1243:   if (ksp)    *ksp    = contour->ksp;
1244:   PetscFunctionReturn(PETSC_SUCCESS);
1245: }

1247: /*@C
1248:    EPSCISSGetKSPs - Retrieve the array of linear solver objects associated with
1249:    the CISS solver.

1251:    Not Collective

1253:    Input Parameter:
1254: .  eps - the eigenproblem solver solver

1256:    Output Parameters:
1257: +  nsolve - number of solver objects
1258: -  ksp - array of linear solver object

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.

1265:    Level: advanced

1267: .seealso: EPSCISSSetSizes()
1268: @*/
1269: PetscErrorCode EPSCISSGetKSPs(EPS eps,PetscInt *nsolve,KSP **ksp)
1270: {
1271:   PetscFunctionBegin;
1273:   PetscUseMethod(eps,"EPSCISSGetKSPs_C",(EPS,PetscInt*,KSP**),(eps,nsolve,ksp));
1274:   PetscFunctionReturn(PETSC_SUCCESS);
1275: }

1277: static PetscErrorCode EPSReset_CISS(EPS eps)
1278: {
1279:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

1281:   PetscFunctionBegin;
1282:   PetscCall(BVDestroy(&ctx->S));
1283:   PetscCall(BVDestroy(&ctx->V));
1284:   PetscCall(BVDestroy(&ctx->Y));
1285:   if (!ctx->usest) PetscCall(SlepcContourDataReset(ctx->contour));
1286:   PetscCall(BVDestroy(&ctx->pV));
1287:   PetscFunctionReturn(PETSC_SUCCESS);
1288: }

1290: static PetscErrorCode EPSSetFromOptions_CISS(EPS eps,PetscOptionItems *PetscOptionsObject)
1291: {
1292:   PetscReal         r3,r4;
1293:   PetscInt          i,i1,i2,i3,i4,i5,i6,i7;
1294:   PetscBool         b1,b2,flg,flg2,flg3,flg4,flg5,flg6;
1295:   EPS_CISS          *ctx = (EPS_CISS*)eps->data;
1296:   EPSCISSQuadRule   quad;
1297:   EPSCISSExtraction extraction;

1299:   PetscFunctionBegin;
1300:   PetscOptionsHeadBegin(PetscOptionsObject,"EPS CISS Options");

1302:     PetscCall(EPSCISSGetSizes(eps,&i1,&i2,&i3,&i4,&i5,&b1));
1303:     PetscCall(PetscOptionsInt("-eps_ciss_integration_points","Number of integration points","EPSCISSSetSizes",i1,&i1,&flg));
1304:     PetscCall(PetscOptionsInt("-eps_ciss_blocksize","Block size","EPSCISSSetSizes",i2,&i2,&flg2));
1305:     PetscCall(PetscOptionsInt("-eps_ciss_moments","Moment size","EPSCISSSetSizes",i3,&i3,&flg3));
1306:     PetscCall(PetscOptionsInt("-eps_ciss_partitions","Number of partitions","EPSCISSSetSizes",i4,&i4,&flg4));
1307:     PetscCall(PetscOptionsInt("-eps_ciss_maxblocksize","Maximum block size","EPSCISSSetSizes",i5,&i5,&flg5));
1308:     PetscCall(PetscOptionsBool("-eps_ciss_realmats","True if A and B are real","EPSCISSSetSizes",b1,&b1,&flg6));
1309:     if (flg || flg2 || flg3 || flg4 || flg5 || flg6) PetscCall(EPSCISSSetSizes(eps,i1,i2,i3,i4,i5,b1));

1311:     PetscCall(EPSCISSGetThreshold(eps,&r3,&r4));
1312:     PetscCall(PetscOptionsReal("-eps_ciss_delta","Threshold for numerical rank","EPSCISSSetThreshold",r3,&r3,&flg));
1313:     PetscCall(PetscOptionsReal("-eps_ciss_spurious_threshold","Threshold for the spurious eigenpairs","EPSCISSSetThreshold",r4,&r4,&flg2));
1314:     if (flg || flg2) PetscCall(EPSCISSSetThreshold(eps,r3,r4));

1316:     PetscCall(EPSCISSGetRefinement(eps,&i6,&i7));
1317:     PetscCall(PetscOptionsInt("-eps_ciss_refine_inner","Number of inner iterative refinement iterations","EPSCISSSetRefinement",i6,&i6,&flg));
1318:     PetscCall(PetscOptionsInt("-eps_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","EPSCISSSetRefinement",i7,&i7,&flg2));
1319:     if (flg || flg2) PetscCall(EPSCISSSetRefinement(eps,i6,i7));

1321:     PetscCall(EPSCISSGetUseST(eps,&b2));
1322:     PetscCall(PetscOptionsBool("-eps_ciss_usest","Use ST for linear solves","EPSCISSSetUseST",b2,&b2,&flg));
1323:     if (flg) PetscCall(EPSCISSSetUseST(eps,b2));

1325:     PetscCall(PetscOptionsEnum("-eps_ciss_quadrule","Quadrature rule","EPSCISSSetQuadRule",EPSCISSQuadRules,(PetscEnum)ctx->quad,(PetscEnum*)&quad,&flg));
1326:     if (flg) PetscCall(EPSCISSSetQuadRule(eps,quad));

1328:     PetscCall(PetscOptionsEnum("-eps_ciss_extraction","Extraction technique","EPSCISSSetExtraction",EPSCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg));
1329:     if (flg) PetscCall(EPSCISSSetExtraction(eps,extraction));

1331:   PetscOptionsHeadEnd();

1333:   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
1334:   PetscCall(RGSetFromOptions(eps->rg)); /* this is necessary here to set useconj */
1335:   if (!ctx->contour || !ctx->contour->ksp) PetscCall(EPSCISSGetKSPs(eps,NULL,NULL));
1336:   PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Something went wrong with EPSCISSGetKSPs()");
1337:   for (i=0;i<ctx->contour->npoints;i++) PetscCall(KSPSetFromOptions(ctx->contour->ksp[i]));
1338:   PetscCall(PetscSubcommSetFromOptions(ctx->contour->subcomm));
1339:   PetscFunctionReturn(PETSC_SUCCESS);
1340: }

1342: static PetscErrorCode EPSDestroy_CISS(EPS eps)
1343: {
1344:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

1346:   PetscFunctionBegin;
1347:   PetscCall(SlepcContourDataDestroy(&ctx->contour));
1348:   PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
1349:   PetscCall(PetscFree(eps->data));
1350:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",NULL));
1351:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",NULL));
1352:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",NULL));
1353:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",NULL));
1354:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",NULL));
1355:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",NULL));
1356:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetUseST_C",NULL));
1357:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetUseST_C",NULL));
1358:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetQuadRule_C",NULL));
1359:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetQuadRule_C",NULL));
1360:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetExtraction_C",NULL));
1361:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetExtraction_C",NULL));
1362:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetKSPs_C",NULL));
1363:   PetscFunctionReturn(PETSC_SUCCESS);
1364: }

1366: static PetscErrorCode EPSView_CISS(EPS eps,PetscViewer viewer)
1367: {
1368:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1369:   PetscBool      isascii;
1370:   PetscViewer    sviewer;

1372:   PetscFunctionBegin;
1373:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
1374:   if (isascii) {
1375:     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:     if (ctx->isreal) PetscCall(PetscViewerASCIIPrintf(viewer,"  exploiting symmetry of integration points\n"));
1377:     PetscCall(PetscViewerASCIIPrintf(viewer,"  threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold));
1378:     PetscCall(PetscViewerASCIIPrintf(viewer,"  iterative refinement { inner: %" PetscInt_FMT ", blocksize: %" PetscInt_FMT " }\n",ctx->refine_inner, ctx->refine_blocksize));
1379:     PetscCall(PetscViewerASCIIPrintf(viewer,"  extraction: %s\n",EPSCISSExtractions[ctx->extraction]));
1380:     PetscCall(PetscViewerASCIIPrintf(viewer,"  quadrature rule: %s\n",EPSCISSQuadRules[ctx->quad]));
1381:     if (ctx->usest) PetscCall(PetscViewerASCIIPrintf(viewer,"  using ST for linear solves\n"));
1382:     else {
1383:       if (!ctx->contour || !ctx->contour->ksp) PetscCall(EPSCISSGetKSPs(eps,NULL,NULL));
1384:       PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Something went wrong with EPSCISSGetKSPs()");
1385:       PetscCall(PetscViewerASCIIPushTab(viewer));
1386:       if (ctx->npart>1 && ctx->contour->subcomm) {
1387:         PetscCall(PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
1388:         if (!ctx->contour->subcomm->color) PetscCall(KSPView(ctx->contour->ksp[0],sviewer));
1389:         PetscCall(PetscViewerFlush(sviewer));
1390:         PetscCall(PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
1391:         /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1392:         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1393:       } else PetscCall(KSPView(ctx->contour->ksp[0],viewer));
1394:       PetscCall(PetscViewerASCIIPopTab(viewer));
1395:     }
1396:   }
1397:   PetscFunctionReturn(PETSC_SUCCESS);
1398: }

1400: static PetscErrorCode EPSSetDefaultST_CISS(EPS eps)
1401: {
1402:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1403:   PetscBool      usest = ctx->usest;
1404:   KSP            ksp;
1405:   PC             pc;

1407:   PetscFunctionBegin;
1408:   if (!((PetscObject)eps->st)->type_name) {
1409:     if (!ctx->usest_set) usest = (ctx->npart>1)? PETSC_FALSE: PETSC_TRUE;
1410:     if (usest) PetscCall(STSetType(eps->st,STSINVERT));
1411:     else {
1412:       /* we are not going to use ST, so avoid factorizing the matrix */
1413:       PetscCall(STSetType(eps->st,STSHIFT));
1414:       if (eps->isgeneralized) {
1415:         PetscCall(STGetKSP(eps->st,&ksp));
1416:         PetscCall(KSPGetPC(ksp,&pc));
1417:         PetscCall(PCSetType(pc,PCNONE));
1418:       }
1419:     }
1420:   }
1421:   PetscFunctionReturn(PETSC_SUCCESS);
1422: }

1424: SLEPC_EXTERN PetscErrorCode EPSCreate_CISS(EPS eps)
1425: {
1426:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

1428:   PetscFunctionBegin;
1429:   PetscCall(PetscNew(&ctx));
1430:   eps->data = ctx;

1432:   eps->useds = PETSC_TRUE;
1433:   eps->categ = EPS_CATEGORY_CONTOUR;

1435:   eps->ops->solve          = EPSSolve_CISS;
1436:   eps->ops->setup          = EPSSetUp_CISS;
1437:   eps->ops->setupsort      = EPSSetUpSort_CISS;
1438:   eps->ops->setfromoptions = EPSSetFromOptions_CISS;
1439:   eps->ops->destroy        = EPSDestroy_CISS;
1440:   eps->ops->reset          = EPSReset_CISS;
1441:   eps->ops->view           = EPSView_CISS;
1442:   eps->ops->computevectors = EPSComputeVectors_CISS;
1443:   eps->ops->setdefaultst   = EPSSetDefaultST_CISS;

1445:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",EPSCISSSetSizes_CISS));
1446:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",EPSCISSGetSizes_CISS));
1447:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",EPSCISSSetThreshold_CISS));
1448:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",EPSCISSGetThreshold_CISS));
1449:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",EPSCISSSetRefinement_CISS));
1450:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",EPSCISSGetRefinement_CISS));
1451:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetUseST_C",EPSCISSSetUseST_CISS));
1452:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetUseST_C",EPSCISSGetUseST_CISS));
1453:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetQuadRule_C",EPSCISSSetQuadRule_CISS));
1454:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetQuadRule_C",EPSCISSGetQuadRule_CISS));
1455:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetExtraction_C",EPSCISSSetExtraction_CISS));
1456:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetExtraction_C",EPSCISSGetExtraction_CISS));
1457:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetKSPs_C",EPSCISSGetKSPs_CISS));

1459:   /* set default values of parameters */
1460:   ctx->N                  = 32;
1461:   ctx->L                  = 16;
1462:   ctx->M                  = ctx->N/4;
1463:   ctx->delta              = SLEPC_DEFAULT_TOL*1e-4;
1464:   ctx->L_max              = 64;
1465:   ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1466:   ctx->usest              = PETSC_TRUE;
1467:   ctx->usest_set          = PETSC_FALSE;
1468:   ctx->isreal             = PETSC_FALSE;
1469:   ctx->refine_inner       = 0;
1470:   ctx->refine_blocksize   = 0;
1471:   ctx->npart              = 1;
1472:   ctx->quad               = (EPSCISSQuadRule)0;
1473:   ctx->extraction         = EPS_CISS_EXTRACTION_RITZ;
1474:   PetscFunctionReturn(PETSC_SUCCESS);
1475: }