Actual source code: nciss.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] J. Asakura, T. Sakurai, H. Tadano, T. Ikegami, K. Kimura, "A
 26:            numerical method for nonlinear eigenvalue problems using contour
 27:            integrals", JSIAM Lett. 1:52-55, 2009.

 29:        [2] S. Yokota and T. Sakurai, "A projection method for nonlinear
 30:            eigenvalue problems using contour integrals", JSIAM Lett.
 31:            5:41-44, 2013.
 32: */

 34: #include <slepc/private/nepimpl.h>
 35: #include <slepc/private/slepccontour.h>

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

 68: struct _n_nep_ciss_project {
 69:   NEP  nep;
 70:   BV   Q;
 71: };

 73: static PetscErrorCode NEPContourDSComputeMatrix(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat,void *ctx)
 74: {
 75:   NEP_CISS_PROJECT proj = (NEP_CISS_PROJECT)ctx;
 76:   Mat              M,fun;

 78:   PetscFunctionBegin;
 79:   if (!deriv) {
 80:     PetscCall(NEPComputeFunction(proj->nep,lambda,proj->nep->function,proj->nep->function));
 81:     fun = proj->nep->function;
 82:   } else {
 83:     PetscCall(NEPComputeJacobian(proj->nep,lambda,proj->nep->jacobian));
 84:     fun = proj->nep->jacobian;
 85:   }
 86:   PetscCall(DSGetMat(ds,mat,&M));
 87:   PetscCall(BVMatProject(proj->Q,fun,proj->Q,M));
 88:   PetscCall(DSRestoreMat(ds,mat,&M));
 89:   PetscFunctionReturn(PETSC_SUCCESS);
 90: }

 92: static PetscErrorCode NEPComputeFunctionSubcomm(NEP nep,PetscScalar lambda,Mat T,Mat P,PetscBool deriv)
 93: {
 94:   PetscInt       i;
 95:   PetscScalar    alpha;
 96:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;

 98:   PetscFunctionBegin;
 99:   PetscAssert(nep->fui!=NEP_USER_INTERFACE_CALLBACK,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONGSTATE,"Should not arrive here with callbacks");
100:   PetscCall(MatZeroEntries(T));
101:   if (!deriv && T != P) PetscCall(MatZeroEntries(P));
102:   for (i=0;i<nep->nt;i++) {
103:     if (!deriv) PetscCall(FNEvaluateFunction(nep->f[i],lambda,&alpha));
104:     else PetscCall(FNEvaluateDerivative(nep->f[i],lambda,&alpha));
105:     PetscCall(MatAXPY(T,alpha,ctx->contour->pA[i],nep->mstr));
106:     if (!deriv && T != P) PetscCall(MatAXPY(P,alpha,ctx->contour->pP[i],nep->mstrp));
107:   }
108:   PetscFunctionReturn(PETSC_SUCCESS);
109: }

111: /*
112:   Set up KSP solvers for every integration point
113: */
114: static PetscErrorCode NEPCISSSetUp(NEP nep,Mat T,Mat P)
115: {
116:   NEP_CISS         *ctx = (NEP_CISS*)nep->data;
117:   SlepcContourData contour;
118:   PetscInt         i,p_id;
119:   Mat              Amat,Pmat;

121:   PetscFunctionBegin;
122:   if (!ctx->contour || !ctx->contour->ksp) PetscCall(NEPCISSGetKSPs(nep,NULL,NULL));
123:   contour = ctx->contour;
124:   PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)nep),PETSC_ERR_PLIB,"Something went wrong with NEPCISSGetKSPs()");
125:   for (i=0;i<contour->npoints;i++) {
126:     p_id = i*contour->subcomm->n + contour->subcomm->color;
127:     PetscCall(MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&Amat));
128:     if (T != P) PetscCall(MatDuplicate(P,MAT_DO_NOT_COPY_VALUES,&Pmat)); else Pmat = Amat;
129:     if (contour->subcomm->n == 1 || nep->fui==NEP_USER_INTERFACE_CALLBACK) PetscCall(NEPComputeFunction(nep,ctx->omega[p_id],Amat,Pmat));
130:     else PetscCall(NEPComputeFunctionSubcomm(nep,ctx->omega[p_id],Amat,Pmat,PETSC_FALSE));
131:     PetscCall(NEP_KSPSetOperators(contour->ksp[i],Amat,Pmat));
132:     PetscCall(MatDestroy(&Amat));
133:     if (T != P) PetscCall(MatDestroy(&Pmat));
134:   }
135:   PetscFunctionReturn(PETSC_SUCCESS);
136: }

138: /*
139:   Y_i = F(z_i)^{-1}Fp(z_i)V for every integration point, Y=[Y_i] is in the context
140: */
141: static PetscErrorCode NEPCISSSolve(NEP nep,Mat dT,BV V,PetscInt L_start,PetscInt L_end)
142: {
143:   NEP_CISS         *ctx = (NEP_CISS*)nep->data;
144:   SlepcContourData contour = ctx->contour;
145:   PetscInt         i,p_id;
146:   Mat              MV,BMV=NULL,MC;

148:   PetscFunctionBegin;
149:   PetscCall(BVSetActiveColumns(V,L_start,L_end));
150:   PetscCall(BVGetMat(V,&MV));
151:   for (i=0;i<contour->npoints;i++) {
152:     p_id = i*contour->subcomm->n + contour->subcomm->color;
153:     if (contour->subcomm->n==1 || nep->fui==NEP_USER_INTERFACE_CALLBACK) PetscCall(NEPComputeJacobian(nep,ctx->omega[p_id],dT));
154:     else PetscCall(NEPComputeFunctionSubcomm(nep,ctx->omega[p_id],dT,NULL,PETSC_TRUE));
155:     PetscCall(BVSetActiveColumns(ctx->Y,i*ctx->L+L_start,i*ctx->L+L_end));
156:     PetscCall(BVGetMat(ctx->Y,&MC));
157:     if (!i) {
158:       PetscCall(MatProductCreate(dT,MV,NULL,&BMV));
159:       PetscCall(MatProductSetType(BMV,MATPRODUCT_AB));
160:       PetscCall(MatProductSetFromOptions(BMV));
161:       PetscCall(MatProductSymbolic(BMV));
162:     }
163:     PetscCall(MatProductNumeric(BMV));
164:     PetscCall(KSPMatSolve(contour->ksp[i],BMV,MC));
165:     PetscCall(BVRestoreMat(ctx->Y,&MC));
166:   }
167:   PetscCall(MatDestroy(&BMV));
168:   PetscCall(BVRestoreMat(V,&MV));
169:   PetscFunctionReturn(PETSC_SUCCESS);
170: }

172: static PetscErrorCode NEPSetUp_CISS(NEP nep)
173: {
174:   NEP_CISS         *ctx = (NEP_CISS*)nep->data;
175:   SlepcContourData contour;
176:   PetscInt         nwork;
177:   PetscBool        istrivial,isellipse,flg;
178:   NEP_CISS_PROJECT dsctxf;
179:   PetscObjectId    id;
180:   PetscObjectState state;
181:   Vec              v0;

183:   PetscFunctionBegin;
184:   if (nep->ncv==PETSC_DETERMINE) nep->ncv = ctx->L_max*ctx->M;
185:   else {
186:     ctx->L_max = nep->ncv/ctx->M;
187:     if (!ctx->L_max) {
188:       ctx->L_max = 1;
189:       nep->ncv = ctx->L_max*ctx->M;
190:     }
191:   }
192:   ctx->L = PetscMin(ctx->L,ctx->L_max);
193:   if (nep->max_it==PETSC_DETERMINE) nep->max_it = 5;
194:   if (nep->mpd==PETSC_DETERMINE) nep->mpd = nep->ncv;
195:   if (!nep->which) nep->which = NEP_ALL;
196:   PetscCheck(nep->which==NEP_ALL,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver supports only computing all eigenvalues");
197:   NEPCheckUnsupported(nep,NEP_FEATURE_STOPPING | NEP_FEATURE_TWOSIDED);

199:   /* check region */
200:   PetscCall(RGIsTrivial(nep->rg,&istrivial));
201:   PetscCheck(!istrivial,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
202:   PetscCall(RGGetComplement(nep->rg,&flg));
203:   PetscCheck(!flg,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
204:   PetscCall(PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse));
205:   PetscCheck(isellipse,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Currently only implemented for elliptic regions");

207:   /* if the region has changed, then reset contour data */
208:   PetscCall(PetscObjectGetId((PetscObject)nep->rg,&id));
209:   PetscCall(PetscObjectStateGet((PetscObject)nep->rg,&state));
210:   if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
211:     PetscCall(SlepcContourDataDestroy(&ctx->contour));
212:     PetscCall(PetscInfo(nep,"Resetting the contour data structure due to a change of region\n"));
213:     ctx->rgid = id; ctx->rgstate = state;
214:   }

216:   /* create contour data structure */
217:   if (!ctx->contour) {
218:     PetscCall(RGCanUseConjugates(nep->rg,ctx->isreal,&ctx->useconj));
219:     PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)nep,&ctx->contour));
220:   }

222:   PetscCall(NEPAllocateSolution(nep,0));
223:   if (ctx->weight) PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
224:   PetscCall(PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma));

226:   /* allocate basis vectors */
227:   PetscCall(BVDestroy(&ctx->S));
228:   PetscCall(BVDuplicateResize(nep->V,ctx->L*ctx->M,&ctx->S));
229:   PetscCall(BVDestroy(&ctx->V));
230:   PetscCall(BVDuplicateResize(nep->V,ctx->L,&ctx->V));

232:   contour = ctx->contour;
233:   if (contour->subcomm && contour->subcomm->n != 1 && nep->fui==NEP_USER_INTERFACE_CALLBACK) {
234:     PetscCall(NEPComputeFunction(nep,0,nep->function,nep->function_pre));
235:     PetscCall(SlepcContourRedundantMat(contour,1,&nep->function,(nep->function!=nep->function_pre)?&nep->function_pre:NULL));
236:   } else PetscCall(SlepcContourRedundantMat(contour,nep->nt,nep->A,nep->P));
237:   if (contour->pA) {
238:     if (!ctx->J) PetscCall(MatDuplicate(contour->pA[0],MAT_DO_NOT_COPY_VALUES,&ctx->J));
239:     PetscCall(BVGetColumn(ctx->V,0,&v0));
240:     PetscCall(SlepcContourScatterCreate(contour,v0));
241:     PetscCall(BVRestoreColumn(ctx->V,0,&v0));
242:     PetscCall(BVDestroy(&ctx->pV));
243:     PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV));
244:     PetscCall(BVSetSizesFromVec(ctx->pV,contour->xsub,nep->n));
245:     PetscCall(BVSetFromOptions(ctx->pV));
246:     PetscCall(BVResize(ctx->pV,ctx->L,PETSC_FALSE));
247:   }

249:   PetscCall(BVDestroy(&ctx->Y));
250:   if (contour->pA) {
251:     PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y));
252:     PetscCall(BVSetSizesFromVec(ctx->Y,contour->xsub,nep->n));
253:     PetscCall(BVSetFromOptions(ctx->Y));
254:     PetscCall(BVResize(ctx->Y,contour->npoints*ctx->L,PETSC_FALSE));
255:   } else PetscCall(BVDuplicateResize(nep->V,contour->npoints*ctx->L,&ctx->Y));

257:   if (ctx->extraction == NEP_CISS_EXTRACTION_RITZ)  {
258:     PetscCall(DSSetMethod(nep->ds,1));
259:     PetscCall(DSNEPSetRG(nep->ds,nep->rg));
260:     if (nep->fui==NEP_USER_INTERFACE_SPLIT) PetscCall(DSNEPSetFN(nep->ds,nep->nt,nep->f));
261:     else {
262:       PetscCall(PetscNew(&dsctxf));
263:       PetscCall(DSNEPSetComputeMatrixFunction(nep->ds,NEPContourDSComputeMatrix,dsctxf));
264:       dsctxf->nep = nep;
265:       ctx->dsctxf = dsctxf;
266:     }
267:   }
268:   PetscCall(DSAllocate(nep->ds,nep->ncv));
269:   nwork = (nep->fui==NEP_USER_INTERFACE_SPLIT)? 2: 1;
270:   PetscCall(NEPSetWorkVecs(nep,nwork));
271:   PetscFunctionReturn(PETSC_SUCCESS);
272: }

274: static PetscErrorCode NEPSolve_CISS(NEP nep)
275: {
276:   NEP_CISS         *ctx = (NEP_CISS*)nep->data;
277:   SlepcContourData contour = ctx->contour;
278:   Mat              X,M,E,T,P,J;
279:   BV               V;
280:   PetscInt         i,j,ld,L_add=0,nv=0,L_base=ctx->L,inner,*inside;
281:   PetscScalar      *Mu,*H0,*H1,*rr,*temp,center;
282:   PetscReal        error,max_error,radius,rgscale,est_eig,eta;
283:   PetscBool        isellipse,*fl1;
284:   Vec              si;
285:   SlepcSC          sc;
286:   PetscRandom      rand;

288:   PetscFunctionBegin;
289:   PetscCall(DSGetSlepcSC(nep->ds,&sc));
290:   sc->comparison    = SlepcCompareLargestMagnitude;
291:   sc->comparisonctx = NULL;
292:   sc->map           = NULL;
293:   sc->mapobj        = NULL;
294:   PetscCall(DSGetLeadingDimension(nep->ds,&ld));
295:   PetscCall(RGComputeQuadrature(nep->rg,RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight));
296:   if (contour->pA) {
297:     T = contour->pA[0];
298:     P = ((nep->fui==NEP_USER_INTERFACE_SPLIT && nep->P) || (nep->fui==NEP_USER_INTERFACE_CALLBACK && contour->pP))? contour->pP[0]: T;
299:   } else {
300:     T = nep->function;
301:     P = nep->function_pre? nep->function_pre: nep->function;
302:   }
303:   PetscCall(NEPCISSSetUp(nep,T,P));
304:   PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
305:   PetscCall(BVSetRandomSign(ctx->V));
306:   PetscCall(BVGetRandomContext(ctx->V,&rand));
307:   if (contour->pA) {
308:     J = ctx->J;
309:     V = ctx->pV;
310:     PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
311:   } else {
312:     J = nep->jacobian;
313:     V = ctx->V;
314:   }
315:   PetscCall(NEPCISSSolve(nep,J,V,0,ctx->L));
316:   PetscCall(PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse));
317:   if (isellipse) {
318:     PetscCall(BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig));
319:     PetscCall(PetscInfo(nep,"Estimated eigenvalue count: %f\n",(double)est_eig));
320:     eta = PetscPowReal(10.0,-PetscLog10Real(nep->tol)/ctx->N);
321:     L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
322:     if (L_add>ctx->L_max-ctx->L) {
323:       PetscCall(PetscInfo(nep,"Number of eigenvalues inside the contour path may be too large\n"));
324:       L_add = ctx->L_max-ctx->L;
325:     }
326:   }
327:   /* Updates L after estimate the number of eigenvalue */
328:   if (L_add>0) {
329:     PetscCall(PetscInfo(nep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by Estimate #Eig\n",ctx->L,ctx->L+L_add));
330:     PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
331:     PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
332:     PetscCall(BVSetRandomSign(ctx->V));
333:     if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
334:     ctx->L += L_add;
335:     PetscCall(NEPCISSSolve(nep,J,V,ctx->L-L_add,ctx->L));
336:   }

338:   PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
339:   for (i=0;i<ctx->refine_blocksize;i++) {
340:     PetscCall(BVDotQuadrature(ctx->Y,V,Mu,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj));
341:     PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
342:     PetscCall(PetscLogEventBegin(NEP_CISS_SVD,nep,0,0,0));
343:     PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
344:     PetscCall(PetscLogEventEnd(NEP_CISS_SVD,nep,0,0,0));
345:     if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
346:     L_add = L_base;
347:     if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
348:     PetscCall(PetscInfo(nep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by SVD(H0)\n",ctx->L,ctx->L+L_add));
349:     PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
350:     PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
351:     PetscCall(BVSetRandomSign(ctx->V));
352:     if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
353:     ctx->L += L_add;
354:     PetscCall(NEPCISSSolve(nep,J,V,ctx->L-L_add,ctx->L));
355:     if (L_add) {
356:       PetscCall(PetscFree2(Mu,H0));
357:       PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
358:     }
359:   }

361:   PetscCall(RGGetScale(nep->rg,&rgscale));
362:   PetscCall(RGEllipseGetParameters(nep->rg,&center,&radius,NULL));

364:   if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) PetscCall(PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1));

366:   while (nep->reason == NEP_CONVERGED_ITERATING) {
367:     nep->its++;
368:     for (inner=0;inner<=ctx->refine_inner;inner++) {
369:       if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
370:         PetscCall(BVDotQuadrature(ctx->Y,V,Mu,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj));
371:         PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
372:         PetscCall(PetscLogEventBegin(NEP_CISS_SVD,nep,0,0,0));
373:         PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
374:         PetscCall(PetscLogEventEnd(NEP_CISS_SVD,nep,0,0,0));
375:       } else {
376:         PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
377:         /* compute SVD of S */
378:         PetscCall(BVSVDAndRank(ctx->S,ctx->M,ctx->L,ctx->delta,(ctx->extraction==NEP_CISS_EXTRACTION_CAA)?BV_SVD_METHOD_QR_CAA:BV_SVD_METHOD_QR,H0,ctx->sigma,&nv));
379:       }
380:       PetscCall(PetscInfo(nep,"Estimated rank: nv = %" PetscInt_FMT "\n",nv));
381:       if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
382:         PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
383:         PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
384:         PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
385:         PetscCall(BVCopy(ctx->S,ctx->V));
386:         if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
387:         PetscCall(NEPCISSSolve(nep,J,V,0,ctx->L));
388:       } else break;
389:     }
390:     nep->nconv = 0;
391:     if (nv == 0) { nep->reason = NEP_CONVERGED_TOL; break; }
392:     else {
393:       /* Extracting eigenpairs */
394:       PetscCall(DSSetDimensions(nep->ds,nv,0,0));
395:       PetscCall(DSSetState(nep->ds,DS_STATE_RAW));
396:       if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
397:         PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
398:         PetscCall(CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1));
399:         PetscCall(DSGetArray(nep->ds,DS_MAT_A,&temp));
400:         for (j=0;j<nv;j++)
401:           for (i=0;i<nv;i++)
402:             temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
403:         PetscCall(DSRestoreArray(nep->ds,DS_MAT_A,&temp));
404:         PetscCall(DSGetArray(nep->ds,DS_MAT_B,&temp));
405:         for (j=0;j<nv;j++)
406:           for (i=0;i<nv;i++)
407:             temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
408:         PetscCall(DSRestoreArray(nep->ds,DS_MAT_B,&temp));
409:       } else if (ctx->extraction == NEP_CISS_EXTRACTION_CAA) {
410:         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
411:         PetscCall(DSGetArray(nep->ds,DS_MAT_A,&temp));
412:         for (i=0;i<nv;i++) PetscCall(PetscArraycpy(temp+i*ld,H0+i*nv,nv));
413:         PetscCall(DSRestoreArray(nep->ds,DS_MAT_A,&temp));
414:       } else {
415:         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
416:         if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
417:           for (i=0;i<nep->nt;i++) {
418:             PetscCall(DSGetMat(nep->ds,DSMatExtra[i],&E));
419:             PetscCall(BVMatProject(ctx->S,nep->A[i],ctx->S,E));
420:             PetscCall(DSRestoreMat(nep->ds,DSMatExtra[i],&E));
421:           }
422:         } else { ctx->dsctxf->Q = ctx->S; }
423:       }
424:       PetscCall(DSSolve(nep->ds,nep->eigr,nep->eigi));
425:       PetscCall(DSSynchronize(nep->ds,nep->eigr,nep->eigi));
426:       PetscCall(DSGetDimensions(nep->ds,NULL,NULL,NULL,&nv));
427:       if (ctx->extraction == NEP_CISS_EXTRACTION_CAA || ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
428:         for (i=0;i<nv;i++) {
429:           nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
430:         }
431:       }
432:       PetscCall(PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr));
433:       PetscCall(DSVectors(nep->ds,DS_MAT_X,NULL,NULL));
434:       PetscCall(DSGetMat(nep->ds,DS_MAT_X,&X));
435:       PetscCall(SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1));
436:       PetscCall(DSRestoreMat(nep->ds,DS_MAT_X,&X));
437:       PetscCall(RGCheckInside(nep->rg,nv,nep->eigr,nep->eigi,inside));
438:       for (i=0;i<nv;i++) {
439:         if (fl1[i] && inside[i]>=0) {
440:           rr[i] = 1.0;
441:           nep->nconv++;
442:         } else rr[i] = 0.0;
443:       }
444:       PetscCall(DSSort(nep->ds,nep->eigr,nep->eigi,rr,NULL,&nep->nconv));
445:       PetscCall(DSSynchronize(nep->ds,nep->eigr,nep->eigi));
446:       if (ctx->extraction == NEP_CISS_EXTRACTION_CAA || ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
447:         for (i=0;i<nv;i++) nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
448:       }
449:       PetscCall(PetscFree3(fl1,inside,rr));
450:       PetscCall(BVSetActiveColumns(nep->V,0,nv));
451:       PetscCall(DSVectors(nep->ds,DS_MAT_X,NULL,NULL));
452:       if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) {
453:         PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
454:         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
455:         PetscCall(BVCopy(ctx->S,nep->V));
456:         PetscCall(DSGetMat(nep->ds,DS_MAT_X,&X));
457:         PetscCall(BVMultInPlace(ctx->S,X,0,nep->nconv));
458:         PetscCall(BVMultInPlace(nep->V,X,0,nep->nconv));
459:         PetscCall(DSRestoreMat(nep->ds,DS_MAT_X,&X));
460:       } else {
461:         PetscCall(DSGetMat(nep->ds,DS_MAT_X,&X));
462:         PetscCall(BVMultInPlace(ctx->S,X,0,nep->nconv));
463:         PetscCall(DSRestoreMat(nep->ds,DS_MAT_X,&X));
464:         PetscCall(BVSetActiveColumns(ctx->S,0,nep->nconv));
465:         PetscCall(BVCopy(ctx->S,nep->V));
466:       }
467:       max_error = 0.0;
468:       for (i=0;i<nep->nconv;i++) {
469:         PetscCall(BVGetColumn(nep->V,i,&si));
470:         PetscCall(VecNormalize(si,NULL));
471:         PetscCall(NEPComputeResidualNorm_Private(nep,PETSC_FALSE,nep->eigr[i],si,nep->work,&error));
472:         PetscCall((*nep->converged)(nep,nep->eigr[i],0,error,&error,nep->convergedctx));
473:         PetscCall(BVRestoreColumn(nep->V,i,&si));
474:         max_error = PetscMax(max_error,error);
475:       }
476:       if (max_error <= nep->tol) nep->reason = NEP_CONVERGED_TOL;
477:       else if (nep->its > nep->max_it) nep->reason = NEP_DIVERGED_ITS;
478:       else {
479:         if (nep->nconv > ctx->L) nv = nep->nconv;
480:         else if (ctx->L > nv) nv = ctx->L;
481:         nv = PetscMin(nv,ctx->L*ctx->M);
482:         PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M));
483:         PetscCall(MatSetRandom(M,rand));
484:         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
485:         PetscCall(BVMultInPlace(ctx->S,M,0,ctx->L));
486:         PetscCall(MatDestroy(&M));
487:         PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
488:         PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
489:         PetscCall(BVCopy(ctx->S,ctx->V));
490:         if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
491:         PetscCall(NEPCISSSolve(nep,J,V,0,ctx->L));
492:       }
493:     }
494:   }
495:   PetscCall(PetscFree2(Mu,H0));
496:   if (ctx->extraction == NEP_CISS_EXTRACTION_HANKEL) PetscCall(PetscFree(H1));
497:   PetscFunctionReturn(PETSC_SUCCESS);
498: }

500: static PetscErrorCode NEPCISSSetSizes_CISS(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
501: {
502:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
503:   PetscInt       oN,oL,oM,oLmax,onpart;
504:   PetscMPIInt    size;

506:   PetscFunctionBegin;
507:   oN = ctx->N;
508:   if (ip == PETSC_DETERMINE) {
509:     if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
510:   } else if (ip != PETSC_CURRENT) {
511:     PetscCheck(ip>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
512:     PetscCheck(ip%2==0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
513:     if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
514:   }
515:   oL = ctx->L;
516:   if (bs == PETSC_DETERMINE) {
517:     ctx->L = 16;
518:   } else if (bs != PETSC_CURRENT) {
519:     PetscCheck(bs>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
520:     ctx->L = bs;
521:   }
522:   oM = ctx->M;
523:   if (ms == PETSC_DETERMINE) {
524:     ctx->M = ctx->N/4;
525:   } else if (ms != PETSC_CURRENT) {
526:     PetscCheck(ms>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
527:     PetscCheck(ms<=ctx->N,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be less than or equal to the number of integration points");
528:     ctx->M = ms;
529:   }
530:   onpart = ctx->npart;
531:   if (npart == PETSC_DETERMINE) {
532:     ctx->npart = 1;
533:   } else if (npart != PETSC_CURRENT) {
534:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)nep),&size));
535:     PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
536:     ctx->npart = npart;
537:   }
538:   oLmax = ctx->L_max;
539:   if (bsmax == PETSC_DETERMINE) {
540:     ctx->L_max = 64;
541:   } else if (bsmax != PETSC_CURRENT) {
542:     PetscCheck(bsmax>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
543:     ctx->L_max = PetscMax(bsmax,ctx->L);
544:   }
545:   if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
546:     PetscCall(SlepcContourDataDestroy(&ctx->contour));
547:     PetscCall(PetscInfo(nep,"Resetting the contour data structure due to a change of parameters\n"));
548:     nep->state = NEP_STATE_INITIAL;
549:   }
550:   ctx->isreal = realmats;
551:   if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) nep->state = NEP_STATE_INITIAL;
552:   PetscFunctionReturn(PETSC_SUCCESS);
553: }

555: /*@
556:    NEPCISSSetSizes - Sets the values of various size parameters in the CISS solver.

558:    Logically Collective

560:    Input Parameters:
561: +  nep   - the nonlinear eigensolver context
562: .  ip    - number of integration points
563: .  bs    - block size
564: .  ms    - moment size
565: .  npart - number of partitions when splitting the communicator
566: .  bsmax - max block size
567: -  realmats - T(z) is real for real z

569:    Options Database Keys:
570: +  -nep_ciss_integration_points - Sets the number of integration points
571: .  -nep_ciss_blocksize - Sets the block size
572: .  -nep_ciss_moments - Sets the moment size
573: .  -nep_ciss_partitions - Sets the number of partitions
574: .  -nep_ciss_maxblocksize - Sets the maximum block size
575: -  -nep_ciss_realmats - T(z) is real for real z

577:    Notes:
578:    For all integer arguments, you can use PETSC_CURRENT to keep the current value, and
579:    PETSC_DETERMINE to set them to a default value.

581:    The default number of partitions is 1. This means the internal KSP object is shared
582:    among all processes of the NEP communicator. Otherwise, the communicator is split
583:    into npart communicators, so that npart KSP solves proceed simultaneously.

585:    The realmats flag can be set to true when T(.) is guaranteed to be real
586:    when the argument is a real value, for example, when all matrices in
587:    the split form are real. When set to true, the solver avoids some computations.

589:    Level: advanced

591: .seealso: NEPCISSGetSizes()
592: @*/
593: PetscErrorCode NEPCISSSetSizes(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
594: {
595:   PetscFunctionBegin;
603:   PetscTryMethod(nep,"NEPCISSSetSizes_C",(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(nep,ip,bs,ms,npart,bsmax,realmats));
604:   PetscFunctionReturn(PETSC_SUCCESS);
605: }

607: static PetscErrorCode NEPCISSGetSizes_CISS(NEP nep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
608: {
609:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

611:   PetscFunctionBegin;
612:   if (ip) *ip = ctx->N;
613:   if (bs) *bs = ctx->L;
614:   if (ms) *ms = ctx->M;
615:   if (npart) *npart = ctx->npart;
616:   if (bsmax) *bsmax = ctx->L_max;
617:   if (realmats) *realmats = ctx->isreal;
618:   PetscFunctionReturn(PETSC_SUCCESS);
619: }

621: /*@
622:    NEPCISSGetSizes - Gets the values of various size parameters in the CISS solver.

624:    Not Collective

626:    Input Parameter:
627: .  nep - the nonlinear eigensolver context

629:    Output Parameters:
630: +  ip    - number of integration points
631: .  bs    - block size
632: .  ms    - moment size
633: .  npart - number of partitions when splitting the communicator
634: .  bsmax - max block size
635: -  realmats - T(z) is real for real z

637:    Level: advanced

639: .seealso: NEPCISSSetSizes()
640: @*/
641: PetscErrorCode NEPCISSGetSizes(NEP nep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
642: {
643:   PetscFunctionBegin;
645:   PetscUseMethod(nep,"NEPCISSGetSizes_C",(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(nep,ip,bs,ms,npart,bsmax,realmats));
646:   PetscFunctionReturn(PETSC_SUCCESS);
647: }

649: static PetscErrorCode NEPCISSSetThreshold_CISS(NEP nep,PetscReal delta,PetscReal spur)
650: {
651:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

653:   PetscFunctionBegin;
654:   if (delta == (PetscReal)PETSC_DETERMINE) {
655:     ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
656:   } else if (delta != (PetscReal)PETSC_CURRENT) {
657:     PetscCheck(delta>0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
658:     ctx->delta = delta;
659:   }
660:   if (spur == (PetscReal)PETSC_DETERMINE) {
661:     ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
662:   } else if (spur != (PetscReal)PETSC_CURRENT) {
663:     PetscCheck(spur>0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
664:     ctx->spurious_threshold = spur;
665:   }
666:   PetscFunctionReturn(PETSC_SUCCESS);
667: }

669: /*@
670:    NEPCISSSetThreshold - Sets the values of various threshold parameters in
671:    the CISS solver.

673:    Logically Collective

675:    Input Parameters:
676: +  nep   - the nonlinear eigensolver context
677: .  delta - threshold for numerical rank
678: -  spur  - spurious threshold (to discard spurious eigenpairs)

680:    Options Database Keys:
681: +  -nep_ciss_delta - Sets the delta
682: -  -nep_ciss_spurious_threshold - Sets the spurious threshold

684:    Note:
685:    PETSC_CURRENT can be used to preserve the current value of any of the
686:    arguments, and PETSC_DETERMINE to set them to a default value.

688:    Level: advanced

690: .seealso: NEPCISSGetThreshold()
691: @*/
692: PetscErrorCode NEPCISSSetThreshold(NEP nep,PetscReal delta,PetscReal spur)
693: {
694:   PetscFunctionBegin;
698:   PetscTryMethod(nep,"NEPCISSSetThreshold_C",(NEP,PetscReal,PetscReal),(nep,delta,spur));
699:   PetscFunctionReturn(PETSC_SUCCESS);
700: }

702: static PetscErrorCode NEPCISSGetThreshold_CISS(NEP nep,PetscReal *delta,PetscReal *spur)
703: {
704:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

706:   PetscFunctionBegin;
707:   if (delta) *delta = ctx->delta;
708:   if (spur)  *spur = ctx->spurious_threshold;
709:   PetscFunctionReturn(PETSC_SUCCESS);
710: }

712: /*@
713:    NEPCISSGetThreshold - Gets the values of various threshold parameters in
714:    the CISS solver.

716:    Not Collective

718:    Input Parameter:
719: .  nep - the nonlinear eigensolver context

721:    Output Parameters:
722: +  delta - threshold for numerical rank
723: -  spur  - spurious threshold (to discard spurious eigenpairs)

725:    Level: advanced

727: .seealso: NEPCISSSetThreshold()
728: @*/
729: PetscErrorCode NEPCISSGetThreshold(NEP nep,PetscReal *delta,PetscReal *spur)
730: {
731:   PetscFunctionBegin;
733:   PetscUseMethod(nep,"NEPCISSGetThreshold_C",(NEP,PetscReal*,PetscReal*),(nep,delta,spur));
734:   PetscFunctionReturn(PETSC_SUCCESS);
735: }

737: static PetscErrorCode NEPCISSSetRefinement_CISS(NEP nep,PetscInt inner,PetscInt blsize)
738: {
739:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

741:   PetscFunctionBegin;
742:   if (inner == PETSC_DETERMINE) {
743:     ctx->refine_inner = 0;
744:   } else if (inner != PETSC_CURRENT) {
745:     PetscCheck(inner>=0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
746:     ctx->refine_inner = inner;
747:   }
748:   if (blsize == PETSC_DETERMINE) {
749:     ctx->refine_blocksize = 0;
750:   } else if (blsize != PETSC_CURRENT) {
751:     PetscCheck(blsize>=0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
752:     ctx->refine_blocksize = blsize;
753:   }
754:   PetscFunctionReturn(PETSC_SUCCESS);
755: }

757: /*@
758:    NEPCISSSetRefinement - Sets the values of various refinement parameters
759:    in the CISS solver.

761:    Logically Collective

763:    Input Parameters:
764: +  nep    - the nonlinear eigensolver context
765: .  inner  - number of iterative refinement iterations (inner loop)
766: -  blsize - number of iterative refinement iterations (blocksize loop)

768:    Options Database Keys:
769: +  -nep_ciss_refine_inner - Sets number of inner iterations
770: -  -nep_ciss_refine_blocksize - Sets number of blocksize iterations

772:    Note:
773:    PETSC_CURRENT can be used to preserve the current value of any of the
774:    arguments, and PETSC_DETERMINE to set them to a default value.

776:    Level: advanced

778: .seealso: NEPCISSGetRefinement()
779: @*/
780: PetscErrorCode NEPCISSSetRefinement(NEP nep,PetscInt inner,PetscInt blsize)
781: {
782:   PetscFunctionBegin;
786:   PetscTryMethod(nep,"NEPCISSSetRefinement_C",(NEP,PetscInt,PetscInt),(nep,inner,blsize));
787:   PetscFunctionReturn(PETSC_SUCCESS);
788: }

790: static PetscErrorCode NEPCISSGetRefinement_CISS(NEP nep,PetscInt *inner,PetscInt *blsize)
791: {
792:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

794:   PetscFunctionBegin;
795:   if (inner)  *inner = ctx->refine_inner;
796:   if (blsize) *blsize = ctx->refine_blocksize;
797:   PetscFunctionReturn(PETSC_SUCCESS);
798: }

800: /*@
801:    NEPCISSGetRefinement - Gets the values of various refinement parameters
802:    in the CISS solver.

804:    Not Collective

806:    Input Parameter:
807: .  nep - the nonlinear eigensolver context

809:    Output Parameters:
810: +  inner  - number of iterative refinement iterations (inner loop)
811: -  blsize - number of iterative refinement iterations (blocksize loop)

813:    Level: advanced

815: .seealso: NEPCISSSetRefinement()
816: @*/
817: PetscErrorCode NEPCISSGetRefinement(NEP nep, PetscInt *inner, PetscInt *blsize)
818: {
819:   PetscFunctionBegin;
821:   PetscUseMethod(nep,"NEPCISSGetRefinement_C",(NEP,PetscInt*,PetscInt*),(nep,inner,blsize));
822:   PetscFunctionReturn(PETSC_SUCCESS);
823: }

825: static PetscErrorCode NEPCISSSetExtraction_CISS(NEP nep,NEPCISSExtraction extraction)
826: {
827:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

829:   PetscFunctionBegin;
830:   if (ctx->extraction != extraction) {
831:     ctx->extraction = extraction;
832:     nep->state      = NEP_STATE_INITIAL;
833:   }
834:   PetscFunctionReturn(PETSC_SUCCESS);
835: }

837: /*@
838:    NEPCISSSetExtraction - Sets the extraction technique used in the CISS solver.

840:    Logically Collective

842:    Input Parameters:
843: +  nep        - the nonlinear eigensolver context
844: -  extraction - the extraction technique

846:    Options Database Key:
847: .  -nep_ciss_extraction - Sets the extraction technique (either 'ritz', 'hankel' or 'caa')

849:    Notes:
850:    By default, the Rayleigh-Ritz extraction is used (NEP_CISS_EXTRACTION_RITZ).

852:    If the 'hankel' or the 'caa' option is specified (NEP_CISS_EXTRACTION_HANKEL or
853:    NEP_CISS_EXTRACTION_CAA), then the Block Hankel method, or the Communication-avoiding
854:    Arnoldi method, respectively, is used for extracting eigenpairs.

856:    Level: advanced

858: .seealso: NEPCISSGetExtraction(), NEPCISSExtraction
859: @*/
860: PetscErrorCode NEPCISSSetExtraction(NEP nep,NEPCISSExtraction extraction)
861: {
862:   PetscFunctionBegin;
865:   PetscTryMethod(nep,"NEPCISSSetExtraction_C",(NEP,NEPCISSExtraction),(nep,extraction));
866:   PetscFunctionReturn(PETSC_SUCCESS);
867: }

869: static PetscErrorCode NEPCISSGetExtraction_CISS(NEP nep,NEPCISSExtraction *extraction)
870: {
871:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

873:   PetscFunctionBegin;
874:   *extraction = ctx->extraction;
875:   PetscFunctionReturn(PETSC_SUCCESS);
876: }

878: /*@
879:    NEPCISSGetExtraction - Gets the extraction technique used in the CISS solver.

881:    Not Collective

883:    Input Parameter:
884: .  nep - the nonlinear eigensolver context

886:    Output Parameters:
887: .  extraction - extraction technique

889:    Level: advanced

891: .seealso: NEPCISSSetExtraction() NEPCISSExtraction
892: @*/
893: PetscErrorCode NEPCISSGetExtraction(NEP nep,NEPCISSExtraction *extraction)
894: {
895:   PetscFunctionBegin;
897:   PetscAssertPointer(extraction,2);
898:   PetscUseMethod(nep,"NEPCISSGetExtraction_C",(NEP,NEPCISSExtraction*),(nep,extraction));
899:   PetscFunctionReturn(PETSC_SUCCESS);
900: }

902: static PetscErrorCode NEPCISSGetKSPs_CISS(NEP nep,PetscInt *nsolve,KSP **ksp)
903: {
904:   NEP_CISS         *ctx = (NEP_CISS*)nep->data;
905:   SlepcContourData contour;
906:   PetscInt         i;
907:   PC               pc;
908:   MPI_Comm         child;

910:   PetscFunctionBegin;
911:   if (!ctx->contour) {  /* initialize contour data structure first */
912:     PetscCall(RGCanUseConjugates(nep->rg,ctx->isreal,&ctx->useconj));
913:     PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)nep,&ctx->contour));
914:   }
915:   contour = ctx->contour;
916:   if (!contour->ksp) {
917:     PetscCall(PetscMalloc1(contour->npoints,&contour->ksp));
918:     PetscCall(PetscSubcommGetChild(contour->subcomm,&child));
919:     for (i=0;i<contour->npoints;i++) {
920:       PetscCall(KSPCreate(child,&contour->ksp[i]));
921:       PetscCall(PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)nep,1));
922:       PetscCall(KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)nep)->prefix));
923:       PetscCall(KSPAppendOptionsPrefix(contour->ksp[i],"nep_ciss_"));
924:       PetscCall(PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)nep)->options));
925:       PetscCall(KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE));
926:       PetscCall(KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(nep->tol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
927:       PetscCall(KSPGetPC(contour->ksp[i],&pc));
928:       if ((nep->fui==NEP_USER_INTERFACE_SPLIT && nep->P) || (nep->fui==NEP_USER_INTERFACE_CALLBACK && nep->function_pre!=nep->function)) {
929:         PetscCall(KSPSetType(contour->ksp[i],KSPBCGS));
930:         PetscCall(PCSetType(pc,PCBJACOBI));
931:       } else {
932:         PetscCall(KSPSetType(contour->ksp[i],KSPPREONLY));
933:         PetscCall(PCSetType(pc,PCLU));
934:       }
935:     }
936:   }
937:   if (nsolve) *nsolve = contour->npoints;
938:   if (ksp)    *ksp    = contour->ksp;
939:   PetscFunctionReturn(PETSC_SUCCESS);
940: }

942: /*@C
943:    NEPCISSGetKSPs - Retrieve the array of linear solver objects associated with
944:    the CISS solver.

946:    Collective

948:    Input Parameter:
949: .  nep - nonlinear eigenvalue solver

951:    Output Parameters:
952: +  nsolve - number of solver objects
953: -  ksp - array of linear solver object

955:    Notes:
956:    The number of KSP solvers is equal to the number of integration points divided by
957:    the number of partitions. This value is halved in the case of real matrices with
958:    a region centered at the real axis.

960:    Level: advanced

962: .seealso: NEPCISSSetSizes()
963: @*/
964: PetscErrorCode NEPCISSGetKSPs(NEP nep,PetscInt *nsolve,KSP **ksp)
965: {
966:   PetscFunctionBegin;
968:   PetscUseMethod(nep,"NEPCISSGetKSPs_C",(NEP,PetscInt*,KSP**),(nep,nsolve,ksp));
969:   PetscFunctionReturn(PETSC_SUCCESS);
970: }

972: static PetscErrorCode NEPReset_CISS(NEP nep)
973: {
974:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;

976:   PetscFunctionBegin;
977:   PetscCall(BVDestroy(&ctx->S));
978:   PetscCall(BVDestroy(&ctx->V));
979:   PetscCall(BVDestroy(&ctx->Y));
980:   PetscCall(SlepcContourDataReset(ctx->contour));
981:   PetscCall(MatDestroy(&ctx->J));
982:   PetscCall(BVDestroy(&ctx->pV));
983:   if (ctx->extraction == NEP_CISS_EXTRACTION_RITZ && nep->fui==NEP_USER_INTERFACE_CALLBACK) PetscCall(PetscFree(ctx->dsctxf));
984:   PetscFunctionReturn(PETSC_SUCCESS);
985: }

987: static PetscErrorCode NEPSetFromOptions_CISS(NEP nep,PetscOptionItems *PetscOptionsObject)
988: {
989:   NEP_CISS          *ctx = (NEP_CISS*)nep->data;
990:   PetscReal         r1,r2;
991:   PetscInt          i,i1,i2,i3,i4,i5,i6,i7;
992:   PetscBool         b1,flg,flg2,flg3,flg4,flg5,flg6;
993:   NEPCISSExtraction extraction;

995:   PetscFunctionBegin;
996:   PetscOptionsHeadBegin(PetscOptionsObject,"NEP CISS Options");

998:     PetscCall(NEPCISSGetSizes(nep,&i1,&i2,&i3,&i4,&i5,&b1));
999:     PetscCall(PetscOptionsInt("-nep_ciss_integration_points","Number of integration points","NEPCISSSetSizes",i1,&i1,&flg));
1000:     PetscCall(PetscOptionsInt("-nep_ciss_blocksize","Block size","NEPCISSSetSizes",i2,&i2,&flg2));
1001:     PetscCall(PetscOptionsInt("-nep_ciss_moments","Moment size","NEPCISSSetSizes",i3,&i3,&flg3));
1002:     PetscCall(PetscOptionsInt("-nep_ciss_partitions","Number of partitions","NEPCISSSetSizes",i4,&i4,&flg4));
1003:     PetscCall(PetscOptionsInt("-nep_ciss_maxblocksize","Maximum block size","NEPCISSSetSizes",i5,&i5,&flg5));
1004:     PetscCall(PetscOptionsBool("-nep_ciss_realmats","True if T(z) is real for real z","NEPCISSSetSizes",b1,&b1,&flg6));
1005:     if (flg || flg2 || flg3 || flg4 || flg5 || flg6) PetscCall(NEPCISSSetSizes(nep,i1,i2,i3,i4,i5,b1));

1007:     PetscCall(NEPCISSGetThreshold(nep,&r1,&r2));
1008:     PetscCall(PetscOptionsReal("-nep_ciss_delta","Threshold for numerical rank","NEPCISSSetThreshold",r1,&r1,&flg));
1009:     PetscCall(PetscOptionsReal("-nep_ciss_spurious_threshold","Threshold for the spurious eigenpairs","NEPCISSSetThreshold",r2,&r2,&flg2));
1010:     if (flg || flg2) PetscCall(NEPCISSSetThreshold(nep,r1,r2));

1012:     PetscCall(NEPCISSGetRefinement(nep,&i6,&i7));
1013:     PetscCall(PetscOptionsInt("-nep_ciss_refine_inner","Number of inner iterative refinement iterations","NEPCISSSetRefinement",i6,&i6,&flg));
1014:     PetscCall(PetscOptionsInt("-nep_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","NEPCISSSetRefinement",i7,&i7,&flg2));
1015:     if (flg || flg2) PetscCall(NEPCISSSetRefinement(nep,i6,i7));

1017:     PetscCall(PetscOptionsEnum("-nep_ciss_extraction","Extraction technique","NEPCISSSetExtraction",NEPCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg));
1018:     if (flg) PetscCall(NEPCISSSetExtraction(nep,extraction));

1020:   PetscOptionsHeadEnd();

1022:   if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg));
1023:   PetscCall(RGSetFromOptions(nep->rg)); /* this is necessary here to set useconj */
1024:   if (!ctx->contour || !ctx->contour->ksp) PetscCall(NEPCISSGetKSPs(nep,NULL,NULL));
1025:   PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)nep),PETSC_ERR_PLIB,"Something went wrong with NEPCISSGetKSPs()");
1026:   for (i=0;i<ctx->contour->npoints;i++) PetscCall(KSPSetFromOptions(ctx->contour->ksp[i]));
1027:   PetscCall(PetscSubcommSetFromOptions(ctx->contour->subcomm));
1028:   PetscFunctionReturn(PETSC_SUCCESS);
1029: }

1031: static PetscErrorCode NEPDestroy_CISS(NEP nep)
1032: {
1033:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;

1035:   PetscFunctionBegin;
1036:   PetscCall(SlepcContourDataDestroy(&ctx->contour));
1037:   PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
1038:   PetscCall(PetscFree(nep->data));
1039:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NULL));
1040:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NULL));
1041:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NULL));
1042:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NULL));
1043:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NULL));
1044:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NULL));
1045:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetExtraction_C",NULL));
1046:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetExtraction_C",NULL));
1047:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetKSPs_C",NULL));
1048:   PetscFunctionReturn(PETSC_SUCCESS);
1049: }

1051: static PetscErrorCode NEPView_CISS(NEP nep,PetscViewer viewer)
1052: {
1053:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
1054:   PetscBool      isascii;
1055:   PetscViewer    sviewer;

1057:   PetscFunctionBegin;
1058:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
1059:   if (isascii) {
1060:     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));
1061:     if (ctx->isreal) PetscCall(PetscViewerASCIIPrintf(viewer,"  exploiting symmetry of integration points\n"));
1062:     PetscCall(PetscViewerASCIIPrintf(viewer,"  threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold));
1063:     PetscCall(PetscViewerASCIIPrintf(viewer,"  iterative refinement  { inner: %" PetscInt_FMT ", blocksize: %" PetscInt_FMT " }\n",ctx->refine_inner, ctx->refine_blocksize));
1064:     PetscCall(PetscViewerASCIIPrintf(viewer,"  extraction: %s\n",NEPCISSExtractions[ctx->extraction]));
1065:     if (!ctx->contour || !ctx->contour->ksp) PetscCall(NEPCISSGetKSPs(nep,NULL,NULL));
1066:     PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)nep),PETSC_ERR_PLIB,"Something went wrong with NEPCISSGetKSPs()");
1067:     PetscCall(PetscViewerASCIIPushTab(viewer));
1068:     if (ctx->npart>1 && ctx->contour->subcomm) {
1069:       PetscCall(PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
1070:       if (!ctx->contour->subcomm->color) PetscCall(KSPView(ctx->contour->ksp[0],sviewer));
1071:       PetscCall(PetscViewerFlush(sviewer));
1072:       PetscCall(PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
1073:       /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1074:       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1075:     } else PetscCall(KSPView(ctx->contour->ksp[0],viewer));
1076:     PetscCall(PetscViewerASCIIPopTab(viewer));
1077:   }
1078:   PetscFunctionReturn(PETSC_SUCCESS);
1079: }

1081: static PetscErrorCode NEPSetDSType_CISS(NEP nep)
1082: {
1083:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;

1085:   PetscFunctionBegin;
1086:   switch (ctx->extraction) {
1087:     case NEP_CISS_EXTRACTION_RITZ:   PetscCall(DSSetType(nep->ds,DSNEP)); break;
1088:     case NEP_CISS_EXTRACTION_HANKEL: PetscCall(DSSetType(nep->ds,DSGNHEP)); break;
1089:     case NEP_CISS_EXTRACTION_CAA:    PetscCall(DSSetType(nep->ds,DSNHEP)); break;
1090:   }
1091:   PetscFunctionReturn(PETSC_SUCCESS);
1092: }

1094: SLEPC_EXTERN PetscErrorCode NEPCreate_CISS(NEP nep)
1095: {
1096:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;

1098:   PetscFunctionBegin;
1099:   PetscCall(PetscNew(&ctx));
1100:   nep->data = ctx;
1101:   /* set default values of parameters */
1102:   ctx->N                  = 32;
1103:   ctx->L                  = 16;
1104:   ctx->M                  = ctx->N/4;
1105:   ctx->delta              = SLEPC_DEFAULT_TOL*1e-4;
1106:   ctx->L_max              = 64;
1107:   ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1108:   ctx->isreal             = PETSC_FALSE;
1109:   ctx->npart              = 1;

1111:   nep->useds = PETSC_TRUE;

1113:   nep->ops->solve          = NEPSolve_CISS;
1114:   nep->ops->setup          = NEPSetUp_CISS;
1115:   nep->ops->setfromoptions = NEPSetFromOptions_CISS;
1116:   nep->ops->reset          = NEPReset_CISS;
1117:   nep->ops->destroy        = NEPDestroy_CISS;
1118:   nep->ops->view           = NEPView_CISS;
1119:   nep->ops->setdstype      = NEPSetDSType_CISS;

1121:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NEPCISSSetSizes_CISS));
1122:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NEPCISSGetSizes_CISS));
1123:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NEPCISSSetThreshold_CISS));
1124:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NEPCISSGetThreshold_CISS));
1125:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NEPCISSSetRefinement_CISS));
1126:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NEPCISSGetRefinement_CISS));
1127:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetExtraction_C",NEPCISSSetExtraction_CISS));
1128:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetExtraction_C",NEPCISSGetExtraction_CISS));
1129:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetKSPs_C",NEPCISSGetKSPs_CISS));
1130:   PetscFunctionReturn(PETSC_SUCCESS);
1131: }