Actual source code: pciss.c

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

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

 13:    Method: Contour Integral Spectral Slicing

 15:    Algorithm:

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

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

 23:    References:

 25:        [1] J. Asakura, T. Sakurai, H. Tadano, T. Ikegami, K. Kimura, "A
 26:            numerical method for polynomial eigenvalue problems using contour
 27:            integral", Japan J. Indust. Appl. Math. 27:73-90, 2010.
 28: */

 30: #include <slepc/private/pepimpl.h>
 31: #include <slepc/private/slepccontour.h>

 33: typedef struct {
 34:   /* parameters */
 35:   PetscInt          N;             /* number of integration points (32) */
 36:   PetscInt          L;             /* block size (16) */
 37:   PetscInt          M;             /* moment degree (N/4 = 4) */
 38:   PetscReal         delta;         /* threshold of singular value (1e-12) */
 39:   PetscInt          L_max;         /* maximum number of columns of the source matrix V */
 40:   PetscReal         spurious_threshold; /* discard spurious eigenpairs */
 41:   PetscBool         isreal;        /* T(z) is real for real z */
 42:   PetscInt          npart;         /* number of partitions */
 43:   PetscInt          refine_inner;
 44:   PetscInt          refine_blocksize;
 45:   PEPCISSExtraction extraction;
 46:   /* private data */
 47:   SlepcContourData  contour;
 48:   PetscReal         *sigma;        /* threshold for numerical rank */
 49:   PetscScalar       *weight;
 50:   PetscScalar       *omega;
 51:   PetscScalar       *pp;
 52:   BV                V;
 53:   BV                S;
 54:   BV                Y;
 55:   PetscBool         useconj;
 56:   Mat               T,J;           /* auxiliary matrices */
 57:   BV                pV;
 58:   PetscObjectId     rgid;
 59:   PetscObjectState  rgstate;
 60: } PEP_CISS;

 62: static PetscErrorCode PEPComputeFunction(PEP pep,PetscScalar lambda,Mat T,PetscBool deriv)
 63: {
 64:   PetscErrorCode   ierr;
 65:   PetscInt         i;
 66:   PetscScalar      *coeff;
 67:   Mat              *A;
 68:   MatStructure     str;
 69:   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
 70:   SlepcContourData contour = ctx->contour;

 73:   A = (contour->pA)?contour->pA:pep->A;
 74:   PetscMalloc1(pep->nmat,&coeff);
 75:   if (deriv) {
 76:     PEPEvaluateBasisDerivative(pep,lambda,0,coeff,NULL);
 77:   } else {
 78:     PEPEvaluateBasis(pep,lambda,0,coeff,NULL);
 79:   }
 80:   STGetMatStructure(pep->st,&str);
 81:   MatZeroEntries(T);
 82:   i = deriv?1:0;
 83:   for (;i<pep->nmat;i++) {
 84:     MatAXPY(T,coeff[i],A[i],str);
 85:   }
 86:   PetscFree(coeff);
 87:   return(0);
 88: }

 90: /*
 91:   Y_i = F(z_i)^{-1}Fp(z_i)V for every integration point, Y=[Y_i] is in the context
 92: */
 93: static PetscErrorCode PEPCISSSolveSystem(PEP pep,Mat T,Mat dT,BV V,PetscInt L_start,PetscInt L_end,PetscBool initksp)
 94: {
 95:   PetscErrorCode   ierr;
 96:   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
 97:   SlepcContourData contour = ctx->contour;
 98:   PetscInt         i,p_id;
 99:   Mat              kspMat,MV,BMV=NULL,MC;

102:   if (!ctx->contour || !ctx->contour->ksp) { PEPCISSGetKSPs(pep,NULL,NULL); }
103:   BVSetActiveColumns(V,L_start,L_end);
104:   BVGetMat(V,&MV);
105:   for (i=0;i<contour->npoints;i++) {
106:     p_id = i*contour->subcomm->n + contour->subcomm->color;
107:     if (initksp) {
108:       PEPComputeFunction(pep,ctx->omega[p_id],T,PETSC_FALSE);
109:       MatDuplicate(T,MAT_COPY_VALUES,&kspMat);
110:       KSPSetOperators(contour->ksp[i],kspMat,kspMat);
111:       MatDestroy(&kspMat);
112:     }
113:     PEPComputeFunction(pep,ctx->omega[p_id],dT,PETSC_TRUE);
114:     BVSetActiveColumns(ctx->Y,i*ctx->L_max+L_start,i*ctx->L_max+L_end);
115:     BVGetMat(ctx->Y,&MC);
116:     if (!i) {
117:       MatProductCreate(dT,MV,NULL,&BMV);
118:       MatProductSetType(BMV,MATPRODUCT_AB);
119:       MatProductSetFromOptions(BMV);
120:       MatProductSymbolic(BMV);
121:     }
122:     MatProductNumeric(BMV);
123:     KSPMatSolve(contour->ksp[i],BMV,MC);
124:     BVRestoreMat(ctx->Y,&MC);
125:   }
126:   MatDestroy(&BMV);
127:   BVRestoreMat(V,&MV);
128:   return(0);
129: }

131: PetscErrorCode PEPSetUp_CISS(PEP pep)
132: {
133:   PetscErrorCode   ierr;
134:   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
135:   SlepcContourData contour;
136:   PetscInt         nwork;
137:   PetscBool        istrivial,isellipse,flg;
138:   PetscObjectId    id;
139:   PetscObjectState state;
140:   Vec              v0;

143:   if (pep->ncv==PETSC_DEFAULT) pep->ncv = ctx->L_max*ctx->M;
144:   else {
145:     ctx->L_max = pep->ncv/ctx->M;
146:     if (!ctx->L_max) {
147:       ctx->L_max = 1;
148:       pep->ncv = ctx->L_max*ctx->M;
149:     }
150:   }
151:   ctx->L = PetscMin(ctx->L,ctx->L_max);
152:   if (pep->max_it==PETSC_DEFAULT) pep->max_it = 5;
153:   if (pep->mpd==PETSC_DEFAULT) pep->mpd = pep->ncv;
154:   if (!pep->which) pep->which = PEP_ALL;
155:   if (pep->which!=PEP_ALL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"This solver supports only computing all eigenvalues");
156:   PEPCheckUnsupported(pep,PEP_FEATURE_STOPPING);
157:   PEPCheckIgnored(pep,PEP_FEATURE_SCALE);

159:   /* check region */
160:   RGIsTrivial(pep->rg,&istrivial);
161:   if (istrivial) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
162:   RGGetComplement(pep->rg,&flg);
163:   if (flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
164:   PetscObjectTypeCompare((PetscObject)pep->rg,RGELLIPSE,&isellipse);
165:   if (!isellipse) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Currently only implemented for elliptic regions");

167:   /* if the region has changed, then reset contour data */
168:   PetscObjectGetId((PetscObject)pep->rg,&id);
169:   PetscObjectStateGet((PetscObject)pep->rg,&state);
170:   if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
171:     SlepcContourDataDestroy(&ctx->contour);
172:     PetscInfo(pep,"Resetting the contour data structure due to a change of region\n");
173:     ctx->rgid = id; ctx->rgstate = state;
174:   }

176:   /* create contour data structure */
177:   if (!ctx->contour) {
178:     RGCanUseConjugates(pep->rg,ctx->isreal,&ctx->useconj);
179:     SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)pep,&ctx->contour);
180:   }

182:   PEPAllocateSolution(pep,0);
183:   if (ctx->weight) { PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma); }
184:   PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma);
185:   PetscLogObjectMemory((PetscObject)pep,3*ctx->N*sizeof(PetscScalar)+ctx->L_max*ctx->N*sizeof(PetscReal));

187:   /* allocate basis vectors */
188:   BVDestroy(&ctx->S);
189:   BVDuplicateResize(pep->V,ctx->L_max*ctx->M,&ctx->S);
190:   PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->S);
191:   BVDestroy(&ctx->V);
192:   BVDuplicateResize(pep->V,ctx->L_max,&ctx->V);
193:   PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->V);

195:   contour = ctx->contour;
196:   SlepcContourRedundantMat(contour,pep->nmat,pep->A);
197:   if (!ctx->T) {
198:     MatDuplicate(contour->pA?contour->pA[0]:pep->A[0],MAT_DO_NOT_COPY_VALUES,&ctx->T);
199:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->T);
200:   }
201:   if (!ctx->J) {
202:     MatDuplicate(contour->pA?contour->pA[0]:pep->A[0],MAT_DO_NOT_COPY_VALUES,&ctx->J);
203:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->J);
204:   }
205:   if (contour->pA) {
206:     BVGetColumn(ctx->V,0,&v0);
207:     SlepcContourScatterCreate(contour,v0);
208:     BVRestoreColumn(ctx->V,0,&v0);
209:     BVDestroy(&ctx->pV);
210:     BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV);
211:     BVSetSizesFromVec(ctx->pV,contour->xsub,pep->n);
212:     BVSetFromOptions(ctx->pV);
213:     BVResize(ctx->pV,ctx->L_max,PETSC_FALSE);
214:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->pV);
215:   }

217:   BVDestroy(&ctx->Y);
218:   if (contour->pA) {
219:     BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y);
220:     BVSetSizesFromVec(ctx->Y,contour->xsub,pep->n);
221:     BVSetFromOptions(ctx->Y);
222:     BVResize(ctx->Y,contour->npoints*ctx->L_max,PETSC_FALSE);
223:   } else {
224:     BVDuplicateResize(pep->V,contour->npoints*ctx->L_max,&ctx->Y);
225:   }

227:   if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
228:     DSSetType(pep->ds,DSGNHEP);
229:   } else if (ctx->extraction == PEP_CISS_EXTRACTION_CAA) {
230:     DSSetType(pep->ds,DSNHEP);
231:   } else {
232:     DSSetType(pep->ds,DSPEP);
233:     DSPEPSetDegree(pep->ds,pep->nmat-1);
234:     DSPEPSetCoefficients(pep->ds,pep->pbc);
235:   }
236:   DSAllocate(pep->ds,pep->ncv);
237:   nwork = 2;
238:   PEPSetWorkVecs(pep,nwork);
239:   return(0);
240: }

242: PetscErrorCode PEPSolve_CISS(PEP pep)
243: {
244:   PetscErrorCode   ierr;
245:   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
246:   SlepcContourData contour = ctx->contour;
247:   Mat              X,M,E;
248:   PetscInt         i,j,ld,L_add=0,nv=0,L_base=ctx->L,inner,*inside;
249:   PetscScalar      *Mu,*H0,*H1,*rr,*temp,center;
250:   PetscReal        error,max_error,radius,rgscale,est_eig,eta;
251:   PetscBool        isellipse,*fl1;
252:   Vec              si;
253:   SlepcSC          sc;
254:   PetscRandom      rand;

257:   DSSetFromOptions(pep->ds);
258:   DSGetSlepcSC(pep->ds,&sc);
259:   sc->comparison    = SlepcCompareLargestMagnitude;
260:   sc->comparisonctx = NULL;
261:   sc->map           = NULL;
262:   sc->mapobj        = NULL;
263:   DSGetLeadingDimension(pep->ds,&ld);
264:   RGComputeQuadrature(pep->rg,RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight);
265:   BVSetActiveColumns(ctx->V,0,ctx->L);
266:   BVSetRandomSign(ctx->V);
267:   BVGetRandomContext(ctx->V,&rand);
268:   if (contour->pA) {
269:     BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
270:   }
271:   PEPCISSSolveSystem(pep,ctx->T,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L,PETSC_TRUE);
272:   PetscObjectTypeCompare((PetscObject)pep->rg,RGELLIPSE,&isellipse);
273:   if (isellipse) {
274:     BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L_max,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig);
275:     PetscInfo1(pep,"Estimated eigenvalue count: %f\n",(double)est_eig);
276:     eta = PetscPowReal(10.0,-PetscLog10Real(pep->tol)/ctx->N);
277:     L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
278:     if (L_add>ctx->L_max-ctx->L) {
279:       PetscInfo(pep,"Number of eigenvalues inside the contour path may be too large\n");
280:       L_add = ctx->L_max-ctx->L;
281:     }
282:   }
283:   /* Updates L after estimate the number of eigenvalue */
284:   if (L_add>0) {
285:     PetscInfo2(pep,"Changing L %D -> %D by Estimate #Eig\n",ctx->L,ctx->L+L_add);
286:     BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
287:     BVSetRandomSign(ctx->V);
288:     if (contour->pA) {
289:       BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
290:     }
291:     PEPCISSSolveSystem(pep,ctx->T,ctx->J,(contour->pA)?ctx->pV:ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
292:     ctx->L += L_add;
293:   }

295:   PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
296:   for (i=0;i<ctx->refine_blocksize;i++) {
297:     BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj);
298:     CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
299:     PetscLogEventBegin(PEP_CISS_SVD,pep,0,0,0);
300:     SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
301:     PetscLogEventEnd(PEP_CISS_SVD,pep,0,0,0);
302:     if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
303:     L_add = L_base;
304:     if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
305:     PetscInfo2(pep,"Changing L %D -> %D by SVD(H0)\n",ctx->L,ctx->L+L_add);
306:     BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
307:     BVSetRandomSign(ctx->V);
308:     if (contour->pA) {
309:       BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
310:     }
311:     PEPCISSSolveSystem(pep,ctx->T,ctx->J,(contour->pA)?ctx->pV:ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
312:     ctx->L += L_add;
313:     if (L_add) {
314:       PetscFree2(Mu,H0);
315:       PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
316:     }
317:   }

319:   RGGetScale(pep->rg,&rgscale);
320:   RGEllipseGetParameters(pep->rg,&center,&radius,NULL);

322:   if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
323:     PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1);
324:   }

326:   while (pep->reason == PEP_CONVERGED_ITERATING) {
327:     pep->its++;
328:     for (inner=0;inner<=ctx->refine_inner;inner++) {
329:       if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
330:         BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj);
331:         CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
332:         PetscLogEventBegin(PEP_CISS_SVD,pep,0,0,0);
333:         SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
334:         PetscLogEventEnd(PEP_CISS_SVD,pep,0,0,0);
335:       } else {
336:         BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
337:         /* compute SVD of S */
338:         BVSVDAndRank(ctx->S,ctx->M,ctx->L,ctx->delta,(ctx->extraction==PEP_CISS_EXTRACTION_CAA)?BV_SVD_METHOD_QR_CAA:BV_SVD_METHOD_QR,H0,ctx->sigma,&nv);
339:       }
340:       if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
341:         BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
342:         BVSetActiveColumns(ctx->S,0,ctx->L);
343:         BVSetActiveColumns(ctx->V,0,ctx->L);
344:         BVCopy(ctx->S,ctx->V);
345:         if (contour->pA) {
346:           BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
347:         }
348:         PEPCISSSolveSystem(pep,ctx->T,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L,PETSC_FALSE);
349:       } else break;
350:     }
351:     pep->nconv = 0;
352:     if (nv == 0) { pep->reason = PEP_CONVERGED_TOL; break; }
353:     else {
354:       /* Extracting eigenpairs */
355:       DSSetDimensions(pep->ds,nv,0,0);
356:       DSSetState(pep->ds,DS_STATE_RAW);
357:       if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
358:         CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
359:         CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1);
360:         DSGetArray(pep->ds,DS_MAT_A,&temp);
361:         for (j=0;j<nv;j++)
362:           for (i=0;i<nv;i++)
363:             temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
364:         DSRestoreArray(pep->ds,DS_MAT_A,&temp);
365:         DSGetArray(pep->ds,DS_MAT_B,&temp);
366:         for (j=0;j<nv;j++)
367:           for (i=0;i<nv;i++)
368:             temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
369:         DSRestoreArray(pep->ds,DS_MAT_B,&temp);
370:       } else if (ctx->extraction == PEP_CISS_EXTRACTION_CAA) {
371:         BVSetActiveColumns(ctx->S,0,nv);
372:         DSGetArray(pep->ds,DS_MAT_A,&temp);
373:         for (i=0;i<nv;i++) {
374:           PetscArraycpy(temp+i*ld,H0+i*nv,nv);
375:         }
376:         DSRestoreArray(pep->ds,DS_MAT_A,&temp);
377:       } else {
378:         BVSetActiveColumns(ctx->S,0,nv);
379:         for (i=0;i<pep->nmat;i++) {
380:           DSGetMat(pep->ds,DSMatExtra[i],&E);
381:           BVMatProject(ctx->S,pep->A[i],ctx->S,E);
382:           DSRestoreMat(pep->ds,DSMatExtra[i],&E);
383:         }
384:         nv = (pep->nmat-1)*nv;
385:       }
386:       DSSolve(pep->ds,pep->eigr,pep->eigi);
387:       DSSynchronize(pep->ds,pep->eigr,pep->eigi);
388:       if (ctx->extraction == PEP_CISS_EXTRACTION_CAA || ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
389:         for (i=0;i<nv;i++) {
390:           pep->eigr[i] = (pep->eigr[i]*radius+center)*rgscale;
391:         }
392:       }
393:       PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr);
394:       DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
395:       DSGetMat(pep->ds,DS_MAT_X,&X);
396:       SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1);
397:       MatDestroy(&X);
398:       RGCheckInside(pep->rg,nv,pep->eigr,pep->eigi,inside);
399:       for (i=0;i<nv;i++) {
400:         if (fl1[i] && inside[i]>=0) {
401:           rr[i] = 1.0;
402:           pep->nconv++;
403:         } else rr[i] = 0.0;
404:       }
405:       DSSort(pep->ds,pep->eigr,pep->eigi,rr,NULL,&pep->nconv);
406:       DSSynchronize(pep->ds,pep->eigr,pep->eigi);
407:       if (ctx->extraction == PEP_CISS_EXTRACTION_CAA || ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
408:         for (i=0;i<nv;i++) pep->eigr[i] = (pep->eigr[i]*radius+center)*rgscale;
409:       }
410:       PetscFree3(fl1,inside,rr);
411:       BVSetActiveColumns(pep->V,0,nv);
412:       DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
413:       if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
414:         BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
415:         BVSetActiveColumns(ctx->S,0,nv);
416:         BVCopy(ctx->S,pep->V);
417:         DSGetMat(pep->ds,DS_MAT_X,&X);
418:         BVMultInPlace(ctx->S,X,0,pep->nconv);
419:         BVMultInPlace(pep->V,X,0,pep->nconv);
420:         MatDestroy(&X);
421:       } else {
422:         DSGetMat(pep->ds,DS_MAT_X,&X);
423:         BVMultInPlace(ctx->S,X,0,pep->nconv);
424:         MatDestroy(&X);
425:         BVSetActiveColumns(ctx->S,0,nv);
426:         BVCopy(ctx->S,pep->V);
427:       }
428:       max_error = 0.0;
429:       for (i=0;i<pep->nconv;i++) {
430:         BVGetColumn(pep->V,i,&si);
431:         VecNormalize(si,NULL);
432:         PEPComputeResidualNorm_Private(pep,pep->eigr[i],0,si,NULL,pep->work,&error);
433:         (*pep->converged)(pep,pep->eigr[i],0,error,&error,pep->convergedctx);
434:         BVRestoreColumn(pep->V,i,&si);
435:         max_error = PetscMax(max_error,error);
436:       }
437:       if (max_error <= pep->tol) pep->reason = PEP_CONVERGED_TOL;
438:       else if (pep->its > pep->max_it) pep->reason = PEP_DIVERGED_ITS;
439:       else {
440:         if (pep->nconv > ctx->L) nv = pep->nconv;
441:         else if (ctx->L > nv) nv = ctx->L;
442:         MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M);
443:         MatSetRandom(M,rand);
444:         BVSetActiveColumns(ctx->S,0,nv);
445:         BVMultInPlace(ctx->S,M,0,ctx->L);
446:         MatDestroy(&M);
447:         BVSetActiveColumns(ctx->S,0,ctx->L);
448:         BVSetActiveColumns(ctx->V,0,ctx->L);
449:         BVCopy(ctx->S,ctx->V);
450:         if (contour->pA) {
451:           BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
452:         }
453:         PEPCISSSolveSystem(pep,ctx->T,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L,PETSC_FALSE);
454:       }
455:     }
456:   }
457:   PetscFree2(Mu,H0);
458:   if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
459:     PetscFree(H1);
460:   }
461:   return(0);
462: }

464: static PetscErrorCode PEPCISSSetSizes_CISS(PEP pep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
465: {
467:   PEP_CISS       *ctx = (PEP_CISS*)pep->data;
468:   PetscInt       oN,oL,oM,oLmax,onpart;

471:   oN = ctx->N;
472:   if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
473:     if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
474:   } else {
475:     if (ip<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
476:     if (ip%2) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
477:     if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
478:   }
479:   oL = ctx->L;
480:   if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
481:     ctx->L = 16;
482:   } else {
483:     if (bs<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
484:     ctx->L = bs;
485:   }
486:   oM = ctx->M;
487:   if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
488:     ctx->M = ctx->N/4;
489:   } else {
490:     if (ms<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
491:     if (ms>ctx->N) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be less than or equal to the number of integration points");
492:     ctx->M = PetscMax(ms,2);
493:   }
494:   onpart = ctx->npart;
495:   if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
496:     ctx->npart = 1;
497:   } else {
498:     if (npart<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The npart argument must be > 0");
499:     ctx->npart = npart;
500:   }
501:   oLmax = ctx->L_max;
502:   if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
503:     ctx->L_max = 64;
504:   } else {
505:     if (bsmax<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
506:     ctx->L_max = PetscMax(bsmax,ctx->L);
507:   }
508:   if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
509:     SlepcContourDataDestroy(&ctx->contour);
510:     PetscInfo(pep,"Resetting the contour data structure due to a change of parameters\n");
511:     pep->state = PEP_STATE_INITIAL;
512:   }
513:   ctx->isreal = realmats;
514:   if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) pep->state = PEP_STATE_INITIAL;
515:   return(0);
516: }

518: /*@
519:    PEPCISSSetSizes - Sets the values of various size parameters in the CISS solver.

521:    Logically Collective on pep

523:    Input Parameters:
524: +  pep   - the polynomial eigensolver context
525: .  ip    - number of integration points
526: .  bs    - block size
527: .  ms    - moment size
528: .  npart - number of partitions when splitting the communicator
529: .  bsmax - max block size
530: -  realmats - all coefficient matrices of P(.) are real

532:    Options Database Keys:
533: +  -pep_ciss_integration_points - Sets the number of integration points
534: .  -pep_ciss_blocksize - Sets the block size
535: .  -pep_ciss_moments - Sets the moment size
536: .  -pep_ciss_partitions - Sets the number of partitions
537: .  -pep_ciss_maxblocksize - Sets the maximum block size
538: -  -pep_ciss_realmats - all coefficient matrices of P(.) are real

540:    Notes:
541:    The default number of partitions is 1. This means the internal KSP object is shared
542:    among all processes of the PEP communicator. Otherwise, the communicator is split
543:    into npart communicators, so that npart KSP solves proceed simultaneously.

545:    Level: advanced

547: .seealso: PEPCISSGetSizes()
548: @*/
549: PetscErrorCode PEPCISSSetSizes(PEP pep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
550: {

561:   PetscTryMethod(pep,"PEPCISSSetSizes_C",(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(pep,ip,bs,ms,npart,bsmax,realmats));
562:   return(0);
563: }

565: static PetscErrorCode PEPCISSGetSizes_CISS(PEP pep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
566: {
567:   PEP_CISS *ctx = (PEP_CISS*)pep->data;

570:   if (ip) *ip = ctx->N;
571:   if (bs) *bs = ctx->L;
572:   if (ms) *ms = ctx->M;
573:   if (npart) *npart = ctx->npart;
574:   if (bsmax) *bsmax = ctx->L_max;
575:   if (realmats) *realmats = ctx->isreal;
576:   return(0);
577: }

579: /*@
580:    PEPCISSGetSizes - Gets the values of various size parameters in the CISS solver.

582:    Not Collective

584:    Input Parameter:
585: .  pep - the polynomial eigensolver context

587:    Output Parameters:
588: +  ip    - number of integration points
589: .  bs    - block size
590: .  ms    - moment size
591: .  npart - number of partitions when splitting the communicator
592: .  bsmax - max block size
593: -  realmats - all coefficient matrices of P(.) are real

595:    Level: advanced

597: .seealso: PEPCISSSetSizes()
598: @*/
599: PetscErrorCode PEPCISSGetSizes(PEP pep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
600: {

605:   PetscUseMethod(pep,"PEPCISSGetSizes_C",(PEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(pep,ip,bs,ms,npart,bsmax,realmats));
606:   return(0);
607: }

609: static PetscErrorCode PEPCISSSetThreshold_CISS(PEP pep,PetscReal delta,PetscReal spur)
610: {
611:   PEP_CISS *ctx = (PEP_CISS*)pep->data;

614:   if (delta == PETSC_DEFAULT) {
615:     ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
616:   } else {
617:     if (delta<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
618:     ctx->delta = delta;
619:   }
620:   if (spur == PETSC_DEFAULT) {
621:     ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
622:   } else {
623:     if (spur<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
624:     ctx->spurious_threshold = spur;
625:   }
626:   return(0);
627: }

629: /*@
630:    PEPCISSSetThreshold - Sets the values of various threshold parameters in
631:    the CISS solver.

633:    Logically Collective on pep

635:    Input Parameters:
636: +  pep   - the polynomial eigensolver context
637: .  delta - threshold for numerical rank
638: -  spur  - spurious threshold (to discard spurious eigenpairs)

640:    Options Database Keys:
641: +  -pep_ciss_delta - Sets the delta
642: -  -pep_ciss_spurious_threshold - Sets the spurious threshold

644:    Level: advanced

646: .seealso: PEPCISSGetThreshold()
647: @*/
648: PetscErrorCode PEPCISSSetThreshold(PEP pep,PetscReal delta,PetscReal spur)
649: {

656:   PetscTryMethod(pep,"PEPCISSSetThreshold_C",(PEP,PetscReal,PetscReal),(pep,delta,spur));
657:   return(0);
658: }

660: static PetscErrorCode PEPCISSGetThreshold_CISS(PEP pep,PetscReal *delta,PetscReal *spur)
661: {
662:   PEP_CISS *ctx = (PEP_CISS*)pep->data;

665:   if (delta) *delta = ctx->delta;
666:   if (spur)  *spur = ctx->spurious_threshold;
667:   return(0);
668: }

670: /*@
671:    PEPCISSGetThreshold - Gets the values of various threshold parameters in
672:    the CISS solver.

674:    Not Collective

676:    Input Parameter:
677: .  pep - the polynomial eigensolver context

679:    Output Parameters:
680: +  delta - threshold for numerical rank
681: -  spur  - spurious threshold (to discard spurious eigenpairs)

683:    Level: advanced

685: .seealso: PEPCISSSetThreshold()
686: @*/
687: PetscErrorCode PEPCISSGetThreshold(PEP pep,PetscReal *delta,PetscReal *spur)
688: {

693:   PetscUseMethod(pep,"PEPCISSGetThreshold_C",(PEP,PetscReal*,PetscReal*),(pep,delta,spur));
694:   return(0);
695: }

697: static PetscErrorCode PEPCISSSetRefinement_CISS(PEP pep,PetscInt inner,PetscInt blsize)
698: {
699:   PEP_CISS *ctx = (PEP_CISS*)pep->data;

702:   if (inner == PETSC_DEFAULT) {
703:     ctx->refine_inner = 0;
704:   } else {
705:     if (inner<0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
706:     ctx->refine_inner = inner;
707:   }
708:   if (blsize == PETSC_DEFAULT) {
709:     ctx->refine_blocksize = 0;
710:   } else {
711:     if (blsize<0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
712:     ctx->refine_blocksize = blsize;
713:   }
714:   return(0);
715: }

717: /*@
718:    PEPCISSSetRefinement - Sets the values of various refinement parameters
719:    in the CISS solver.

721:    Logically Collective on pep

723:    Input Parameters:
724: +  pep    - the polynomial eigensolver context
725: .  inner  - number of iterative refinement iterations (inner loop)
726: -  blsize - number of iterative refinement iterations (blocksize loop)

728:    Options Database Keys:
729: +  -pep_ciss_refine_inner - Sets number of inner iterations
730: -  -pep_ciss_refine_blocksize - Sets number of blocksize iterations

732:    Level: advanced

734: .seealso: PEPCISSGetRefinement()
735: @*/
736: PetscErrorCode PEPCISSSetRefinement(PEP pep,PetscInt inner,PetscInt blsize)
737: {

744:   PetscTryMethod(pep,"PEPCISSSetRefinement_C",(PEP,PetscInt,PetscInt),(pep,inner,blsize));
745:   return(0);
746: }

748: static PetscErrorCode PEPCISSGetRefinement_CISS(PEP pep,PetscInt *inner,PetscInt *blsize)
749: {
750:   PEP_CISS *ctx = (PEP_CISS*)pep->data;

753:   if (inner)  *inner = ctx->refine_inner;
754:   if (blsize) *blsize = ctx->refine_blocksize;
755:   return(0);
756: }

758: /*@
759:    PEPCISSGetRefinement - Gets the values of various refinement parameters
760:    in the CISS solver.

762:    Not Collective

764:    Input Parameter:
765: .  pep - the polynomial eigensolver context

767:    Output Parameters:
768: +  inner  - number of iterative refinement iterations (inner loop)
769: -  blsize - number of iterative refinement iterations (blocksize loop)

771:    Level: advanced

773: .seealso: PEPCISSSetRefinement()
774: @*/
775: PetscErrorCode PEPCISSGetRefinement(PEP pep, PetscInt *inner, PetscInt *blsize)
776: {

781:   PetscUseMethod(pep,"PEPCISSGetRefinement_C",(PEP,PetscInt*,PetscInt*),(pep,inner,blsize));
782:   return(0);
783: }

785: static PetscErrorCode PEPCISSSetExtraction_CISS(PEP pep,PEPCISSExtraction extraction)
786: {
787:   PEP_CISS *ctx = (PEP_CISS*)pep->data;

790:   if (ctx->extraction != extraction) {
791:     ctx->extraction = extraction;
792:     pep->state      = PEP_STATE_INITIAL;
793:   }
794:   return(0);
795: }

797: /*@
798:    PEPCISSSetExtraction - Sets the extraction technique used in the CISS solver.

800:    Logically Collective on pep

802:    Input Parameters:
803: +  pep        - the polynomial eigensolver context
804: -  extraction - the extraction technique

806:    Options Database Key:
807: .  -pep_ciss_extraction - Sets the extraction technique (either 'ritz', 'hankel' or 'caa')

809:    Notes:
810:    By default, the Rayleigh-Ritz extraction is used (PEP_CISS_EXTRACTION_RITZ).

812:    If the 'hankel' or the 'caa' option is specified (PEP_CISS_EXTRACTION_HANKEL or
813:    PEP_CISS_EXTRACTION_CAA), then the Block Hankel method, or the Communication-avoiding
814:    Arnoldi method, respectively, is used for extracting eigenpairs.

816:    Level: advanced

818: .seealso: PEPCISSGetExtraction(), PEPCISSExtraction
819: @*/
820: PetscErrorCode PEPCISSSetExtraction(PEP pep,PEPCISSExtraction extraction)
821: {

827:   PetscTryMethod(pep,"PEPCISSSetExtraction_C",(PEP,PEPCISSExtraction),(pep,extraction));
828:   return(0);
829: }

831: static PetscErrorCode PEPCISSGetExtraction_CISS(PEP pep,PEPCISSExtraction *extraction)
832: {
833:   PEP_CISS *ctx = (PEP_CISS*)pep->data;

836:   *extraction = ctx->extraction;
837:   return(0);
838: }

840: /*@
841:    PEPCISSGetExtraction - Gets the extraction technique used in the CISS solver.

843:    Not Collective

845:    Input Parameter:
846: .  pep - the polynomial eigensolver context

848:    Output Parameters:
849: .  extraction - extraction technique

851:    Level: advanced

853: .seealso: PEPCISSSetExtraction() PEPCISSExtraction
854: @*/
855: PetscErrorCode PEPCISSGetExtraction(PEP pep,PEPCISSExtraction *extraction)
856: {

862:   PetscUseMethod(pep,"PEPCISSGetExtraction_C",(PEP,PEPCISSExtraction*),(pep,extraction));
863:   return(0);
864: }

866: static PetscErrorCode PEPCISSGetKSPs_CISS(PEP pep,PetscInt *nsolve,KSP **ksp)
867: {
868:   PetscErrorCode   ierr;
869:   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
870:   SlepcContourData contour;
871:   PetscInt         i;
872:   PC               pc;

875:   if (!ctx->contour) {  /* initialize contour data structure first */
876:     RGCanUseConjugates(pep->rg,ctx->isreal,&ctx->useconj);
877:     SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)pep,&ctx->contour);
878:   }
879:   contour = ctx->contour;
880:   if (!contour->ksp) {
881:     PetscMalloc1(contour->npoints,&contour->ksp);
882:     for (i=0;i<contour->npoints;i++) {
883:       KSPCreate(PetscSubcommChild(contour->subcomm),&contour->ksp[i]);
884:       PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)pep,1);
885:       KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)pep)->prefix);
886:       KSPAppendOptionsPrefix(contour->ksp[i],"pep_ciss_");
887:       PetscLogObjectParent((PetscObject)pep,(PetscObject)contour->ksp[i]);
888:       PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)pep)->options);
889:       KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE);
890:       KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(pep->tol),PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
891:       KSPGetPC(contour->ksp[i],&pc);
892:       KSPSetType(contour->ksp[i],KSPPREONLY);
893:       PCSetType(pc,PCLU);
894:     }
895:   }
896:   if (nsolve) *nsolve = contour->npoints;
897:   if (ksp)    *ksp    = contour->ksp;
898:   return(0);
899: }

901: /*@C
902:    PEPCISSGetKSPs - Retrieve the array of linear solver objects associated with
903:    the CISS solver.

905:    Not Collective

907:    Input Parameter:
908: .  pep - polynomial eigenvalue solver

910:    Output Parameters:
911: +  nsolve - number of solver objects
912: -  ksp - array of linear solver object

914:    Notes:
915:    The number of KSP solvers is equal to the number of integration points divided by
916:    the number of partitions. This value is halved in the case of real matrices with
917:    a region centered at the real axis.

919:    Level: advanced

921: .seealso: PEPCISSSetSizes()
922: @*/
923: PetscErrorCode PEPCISSGetKSPs(PEP pep,PetscInt *nsolve,KSP **ksp)
924: {

929:   PetscUseMethod(pep,"PEPCISSGetKSPs_C",(PEP,PetscInt*,KSP**),(pep,nsolve,ksp));
930:   return(0);
931: }

933: PetscErrorCode PEPReset_CISS(PEP pep)
934: {
936:   PEP_CISS       *ctx = (PEP_CISS*)pep->data;

939:   BVDestroy(&ctx->S);
940:   BVDestroy(&ctx->V);
941:   BVDestroy(&ctx->Y);
942:   SlepcContourDataReset(ctx->contour);
943:   MatDestroy(&ctx->T);
944:   MatDestroy(&ctx->J);
945:   BVDestroy(&ctx->pV);
946:   return(0);
947: }

949: PetscErrorCode PEPSetFromOptions_CISS(PetscOptionItems *PetscOptionsObject,PEP pep)
950: {
951:   PetscErrorCode    ierr;
952:   PEP_CISS          *ctx = (PEP_CISS*)pep->data;
953:   PetscReal         r1,r2;
954:   PetscInt          i,i1,i2,i3,i4,i5,i6,i7;
955:   PetscBool         b1,flg,flg2,flg3,flg4,flg5,flg6;
956:   PEPCISSExtraction extraction;

959:   PetscOptionsHead(PetscOptionsObject,"PEP CISS Options");

961:     PEPCISSGetSizes(pep,&i1,&i2,&i3,&i4,&i5,&b1);
962:     PetscOptionsInt("-pep_ciss_integration_points","Number of integration points","PEPCISSSetSizes",i1,&i1,&flg);
963:     PetscOptionsInt("-pep_ciss_blocksize","Block size","PEPCISSSetSizes",i2,&i2,&flg2);
964:     PetscOptionsInt("-pep_ciss_moments","Moment size","PEPCISSSetSizes",i3,&i3,&flg3);
965:     PetscOptionsInt("-pep_ciss_partitions","Number of partitions","PEPCISSSetSizes",i4,&i4,&flg4);
966:     PetscOptionsInt("-pep_ciss_maxblocksize","Maximum block size","PEPCISSSetSizes",i5,&i5,&flg5);
967:     PetscOptionsBool("-pep_ciss_realmats","True if all coefficient matrices of P(.) are real","PEPCISSSetSizes",b1,&b1,&flg6);
968:     if (flg || flg2 || flg3 || flg4 || flg5 || flg6) { PEPCISSSetSizes(pep,i1,i2,i3,i4,i5,b1); }

970:     PEPCISSGetThreshold(pep,&r1,&r2);
971:     PetscOptionsReal("-pep_ciss_delta","Threshold for numerical rank","PEPCISSSetThreshold",r1,&r1,&flg);
972:     PetscOptionsReal("-pep_ciss_spurious_threshold","Threshold for the spurious eigenpairs","PEPCISSSetThreshold",r2,&r2,&flg2);
973:     if (flg || flg2) { PEPCISSSetThreshold(pep,r1,r2); }

975:     PEPCISSGetRefinement(pep,&i6,&i7);
976:     PetscOptionsInt("-pep_ciss_refine_inner","Number of inner iterative refinement iterations","PEPCISSSetRefinement",i6,&i6,&flg);
977:     PetscOptionsInt("-pep_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","PEPCISSSetRefinement",i7,&i7,&flg2);
978:     if (flg || flg2) { PEPCISSSetRefinement(pep,i6,i7); }

980:     PetscOptionsEnum("-pep_ciss_extraction","Extraction technique","PEPCISSSetExtraction",PEPCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg);
981:     if (flg) { PEPCISSSetExtraction(pep,extraction); }

983:   PetscOptionsTail();

985:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
986:   RGSetFromOptions(pep->rg); /* this is necessary here to set useconj */
987:   if (!ctx->contour || !ctx->contour->ksp) { PEPCISSGetKSPs(pep,NULL,NULL); }
988:   for (i=0;i<ctx->contour->npoints;i++) {
989:     KSPSetFromOptions(ctx->contour->ksp[i]);
990:   }
991:   PetscSubcommSetFromOptions(ctx->contour->subcomm);
992:   return(0);
993: }

995: PetscErrorCode PEPDestroy_CISS(PEP pep)
996: {
998:   PEP_CISS       *ctx = (PEP_CISS*)pep->data;

1001:   SlepcContourDataDestroy(&ctx->contour);
1002:   PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
1003:   PetscFree(pep->data);
1004:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetSizes_C",NULL);
1005:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetSizes_C",NULL);
1006:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetThreshold_C",NULL);
1007:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetThreshold_C",NULL);
1008:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetRefinement_C",NULL);
1009:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetRefinement_C",NULL);
1010:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetExtraction_C",NULL);
1011:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetExtraction_C",NULL);
1012:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetKSPs_C",NULL);
1013:   return(0);
1014: }

1016: PetscErrorCode PEPView_CISS(PEP pep,PetscViewer viewer)
1017: {
1019:   PEP_CISS       *ctx = (PEP_CISS*)pep->data;
1020:   PetscBool      isascii;
1021:   PetscViewer    sviewer;

1024:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1025:   if (isascii) {
1026:     PetscViewerASCIIPrintf(viewer,"  sizes { integration points: %D, block size: %D, moment size: %D, partitions: %D, maximum block size: %D }\n",ctx->N,ctx->L,ctx->M,ctx->npart,ctx->L_max);
1027:     if (ctx->isreal) {
1028:       PetscViewerASCIIPrintf(viewer,"  exploiting symmetry of integration points\n");
1029:     }
1030:     PetscViewerASCIIPrintf(viewer,"  threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold);
1031:     PetscViewerASCIIPrintf(viewer,"  iterative refinement  { inner: %D, blocksize: %D }\n",ctx->refine_inner, ctx->refine_blocksize);
1032:     PetscViewerASCIIPrintf(viewer,"  extraction: %s\n",PEPCISSExtractions[ctx->extraction]);
1033:     if (!ctx->contour || !ctx->contour->ksp) { PEPCISSGetKSPs(pep,NULL,NULL); }
1034:     PetscViewerASCIIPushTab(viewer);
1035:     if (ctx->npart>1 && ctx->contour->subcomm) {
1036:       PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1037:       if (!ctx->contour->subcomm->color) {
1038:         KSPView(ctx->contour->ksp[0],sviewer);
1039:       }
1040:       PetscViewerFlush(sviewer);
1041:       PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1042:       PetscViewerFlush(viewer);
1043:       /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1044:       PetscViewerASCIIPopSynchronized(viewer);
1045:     } else {
1046:       KSPView(ctx->contour->ksp[0],viewer);
1047:     }
1048:     PetscViewerASCIIPopTab(viewer);
1049:   }
1050:   return(0);
1051: }

1053: SLEPC_EXTERN PetscErrorCode PEPCreate_CISS(PEP pep)
1054: {
1056:   PEP_CISS       *ctx = (PEP_CISS*)pep->data;

1059:   PetscNewLog(pep,&ctx);
1060:   pep->data = ctx;
1061:   /* set default values of parameters */
1062:   ctx->N                  = 32;
1063:   ctx->L                  = 16;
1064:   ctx->M                  = ctx->N/4;
1065:   ctx->delta              = SLEPC_DEFAULT_TOL*1e-4;
1066:   ctx->L_max              = 64;
1067:   ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1068:   ctx->isreal             = PETSC_FALSE;
1069:   ctx->npart              = 1;

1071:   pep->ops->solve          = PEPSolve_CISS;
1072:   pep->ops->setup          = PEPSetUp_CISS;
1073:   pep->ops->setfromoptions = PEPSetFromOptions_CISS;
1074:   pep->ops->reset          = PEPReset_CISS;
1075:   pep->ops->destroy        = PEPDestroy_CISS;
1076:   pep->ops->view           = PEPView_CISS;

1078:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetSizes_C",PEPCISSSetSizes_CISS);
1079:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetSizes_C",PEPCISSGetSizes_CISS);
1080:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetThreshold_C",PEPCISSSetThreshold_CISS);
1081:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetThreshold_C",PEPCISSGetThreshold_CISS);
1082:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetRefinement_C",PEPCISSSetRefinement_CISS);
1083:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetRefinement_C",PEPCISSGetRefinement_CISS);
1084:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetExtraction_C",PEPCISSSetExtraction_CISS);
1085:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetExtraction_C",PEPCISSGetExtraction_CISS);
1086:   PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetKSPs_C",PEPCISSGetKSPs_CISS);
1087:   return(0);
1088: }