Actual source code: nciss.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain

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

 13:    Method: Contour Integral Spectral Slicing

 15:    Algorithm:

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

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

 23:    References:

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

 28:        [2] T. Sakurai and H. Tadano, "CIRR: a Rayleigh-Ritz type method with
 29:            contour integral for generalized eigenvalue problems", Hokkaido
 30:            Math. J. 36:745-757, 2007.
 31: */

 33: #include <slepc/private/nepimpl.h>         /*I "slepcnep.h" I*/
 34: #include <slepcblaslapack.h>

 36: typedef struct {
 37:   /* parameters */
 38:   PetscInt     N;          /* number of integration points (32) */
 39:   PetscInt     L;          /* block size (16) */
 40:   PetscInt     M;          /* moment degree (N/4 = 4) */
 41:   PetscReal    delta;      /* threshold of singular value (1e-12) */
 42:   PetscInt     L_max;      /* maximum number of columns of the source matrix V */
 43:   PetscReal    spurious_threshold; /* discard spurious eigenpairs */
 44:   PetscBool    isreal;     /* T(z) is real for real z */
 45:   PetscInt     refine_inner;
 46:   PetscInt     refine_blocksize;
 47:   /* private data */
 48:   PetscReal    *sigma;     /* threshold for numerical rank */
 49:   PetscInt     num_subcomm;
 50:   PetscInt     subcomm_id;
 51:   PetscInt     num_solve_point;
 52:   PetscScalar  *weight;
 53:   PetscScalar  *omega;
 54:   PetscScalar  *pp;
 55:   BV           V;
 56:   BV           S;
 57:   BV           Y;
 58:   KSP          *ksp;
 59:   PetscBool    useconj;
 60:   PetscReal    est_eig;
 61:   PetscSubcomm subcomm;
 62:   PetscBool    usest;
 63: } NEP_CISS;

 65: static PetscErrorCode SetSolverComm(NEP nep)
 66: {
 68:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
 69:   PetscInt       N = ctx->N;

 72:   if (ctx->useconj) N = N/2;
 73:   PetscSubcommDestroy(&ctx->subcomm);
 74:   PetscSubcommCreate(PetscObjectComm((PetscObject)nep),&ctx->subcomm);
 75:   PetscSubcommSetNumber(ctx->subcomm,ctx->num_subcomm);
 76:   PetscSubcommSetType(ctx->subcomm,PETSC_SUBCOMM_INTERLACED);
 77:   PetscLogObjectMemory((PetscObject)nep,sizeof(PetscSubcomm));
 78:   PetscSubcommSetFromOptions(ctx->subcomm);
 79:   ctx->subcomm_id = ctx->subcomm->color;
 80:   ctx->num_solve_point = N / ctx->num_subcomm;
 81:   if ((N%ctx->num_subcomm) > ctx->subcomm_id) ctx->num_solve_point+=1;
 82:   return(0);
 83: }

 85: static PetscErrorCode SetPathParameter(NEP nep)
 86: {
 88:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
 89:   PetscInt       i;
 90:   PetscScalar    center;
 91:   PetscReal      theta,radius,vscale,rgscale;
 92:   PetscBool      isellipse=PETSC_FALSE;

 95:   PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse);
 96:   RGGetScale(nep->rg,&rgscale);
 97:   if (isellipse) {
 98:     RGEllipseGetParameters(nep->rg,&center,&radius,&vscale);
 99:   } else SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Region must be Ellipse");
100:   for (i=0;i<ctx->N;i++) {
101:     theta = ((2*PETSC_PI)/ctx->N)*(i+0.5);
102:     ctx->pp[i] = PetscCosReal(theta) + PETSC_i*vscale*PetscSinReal(theta);
103:     ctx->weight[i] = radius*(vscale*PetscCosReal(theta) + PETSC_i*PetscSinReal(theta))/(PetscReal)ctx->N;
104:     ctx->omega[i] = rgscale*(center + radius*ctx->pp[i]);
105:   }
106:   return(0);
107: }

109: static PetscErrorCode CISSVecSetRandom(BV V,PetscInt i0,PetscInt i1)
110: {
112:   PetscInt       i,j,nlocal;
113:   PetscScalar    *vdata;
114:   Vec            x;

117:   BVGetSizes(V,&nlocal,NULL,NULL);
118:   for (i=i0;i<i1;i++) {
119:     BVSetRandomColumn(V,i);
120:     BVGetColumn(V,i,&x);
121:     VecGetArray(x,&vdata);
122:     for (j=0;j<nlocal;j++) {
123:       vdata[j] = PetscRealPart(vdata[j]);
124:       if (PetscRealPart(vdata[j]) < 0.5) vdata[j] = -1.0;
125:       else vdata[j] = 1.0;
126:     }
127:     VecRestoreArray(x,&vdata);
128:     BVRestoreColumn(V,i,&x);
129:   }
130:   return(0);
131: }

133: static PetscErrorCode SolveLinearSystem(NEP nep,Mat T,Mat dT,BV V,PetscInt L_start,PetscInt L_end,PetscBool initksp)
134: {
136:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
137:   PetscInt       i,j,p_id;
138:   Mat            Fz,kspMat;
139:   PC             pc;
140:   Vec            Bvj,vj,yj;
141:   KSP            ksp;

144:   if (ctx->usest) {
145:     NEPComputeFunction(nep,0,T,T);
146:     MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&Fz);
147:     KSPCreate(PetscObjectComm((PetscObject)nep),&ksp);
148:   }
149:   BVCreateVec(V,&Bvj);
150:   for (i=0;i<ctx->num_solve_point;i++) {
151:     p_id = i*ctx->subcomm->n + ctx->subcomm_id;
152:     NEPComputeFunction(nep,ctx->omega[p_id],T,T);
153:     NEPComputeJacobian(nep,ctx->omega[p_id],dT);
154:     if (!ctx->usest && initksp == PETSC_TRUE) {
155:       MatDuplicate(T,MAT_COPY_VALUES,&kspMat);
156:       KSPSetOperators(ctx->ksp[i],kspMat,kspMat);
157:       MatDestroy(&kspMat);
158:       KSPSetType(ctx->ksp[i],KSPPREONLY);
159:       KSPGetPC(ctx->ksp[i],&pc);
160:       PCSetType(pc,PCLU);
161:       KSPSetFromOptions(ctx->ksp[i]);
162:     } else if (ctx->usest) {
163:       MatCopy(T,Fz,DIFFERENT_NONZERO_PATTERN);
164:       KSPSetOperators(ksp,Fz,Fz);
165:       KSPSetType(ksp,KSPPREONLY);
166:       KSPGetPC(ksp,&pc);
167:       PCSetType(pc,PCLU);
168:       KSPSetTolerances(ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
169:       KSPSetFromOptions(ksp);
170:     }
171:     for (j=L_start;j<L_end;j++) {
172:       BVGetColumn(V,j,&vj);
173:       BVGetColumn(ctx->Y,i*ctx->L_max+j,&yj);
174:       MatMult(dT,vj,Bvj);
175:       if (ctx->usest) {
176:         KSPSolve(ksp,Bvj,yj);
177:       } else {
178:         KSPSolve(ctx->ksp[i],Bvj,yj);
179:       }
180:       BVRestoreColumn(V,j,&vj);
181:       BVRestoreColumn(ctx->Y,i*ctx->L_max+j,&yj);
182:     }
183:     if (ctx->usest && i<ctx->num_solve_point-1) { KSPReset(ksp); }
184:   }
185:   if (ctx->usest) {
186:     MatDestroy(&Fz);
187:     KSPDestroy(&ksp);
188:   }
189:   VecDestroy(&Bvj);
190:   return(0);
191: }

193: static PetscErrorCode EstimateNumberEigs(NEP nep,PetscInt *L_add)
194: {
196:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
197:   PetscInt       i,j,p_id;
198:   PetscScalar    tmp,m = 1,sum = 0.0;
199:   PetscReal      eta;
200:   Vec            v,vtemp,vj,yj;

203:   BVGetColumn(ctx->Y,0,&yj);
204:   VecDuplicate(yj,&v);
205:   BVRestoreColumn(ctx->Y,0,&yj);
206:   BVCreateVec(ctx->V,&vtemp);
207:   for (j=0;j<ctx->L;j++) {
208:     VecSet(v,0);
209:     for (i=0;i<ctx->num_solve_point; i++) {
210:       p_id = i*ctx->subcomm->n + ctx->subcomm_id;
211:       BVSetActiveColumns(ctx->Y,i*ctx->L_max+j,i*ctx->L_max+j+1);
212:       BVMultVec(ctx->Y,ctx->weight[p_id],1,v,&m);
213:     }
214:     BVGetColumn(ctx->V,j,&vj);
215:     VecDot(vj,v,&tmp);
216:     BVRestoreColumn(ctx->V,j,&vj);
217:     if (ctx->useconj) sum += PetscRealPart(tmp)*2;
218:     else sum += tmp;
219:   }
220:   ctx->est_eig = PetscAbsScalar(sum/(PetscReal)ctx->L);
221:   eta = PetscPowReal(10,-PetscLog10Real(nep->tol)/ctx->N);
222:   PetscInfo1(nep,"Estimation_#Eig %f\n",(double)ctx->est_eig);
223:   *L_add = (PetscInt)PetscCeilReal((ctx->est_eig*eta)/ctx->M) - ctx->L;
224:   if (*L_add < 0) *L_add = 0;
225:   if (*L_add>ctx->L_max-ctx->L) {
226:     PetscInfo(nep,"Number of eigenvalues around the contour path may be too large\n");
227:     *L_add = ctx->L_max-ctx->L;
228:   }
229:   VecDestroy(&v);
230:   VecDestroy(&vtemp);
231:   return(0);
232: }

234: static PetscErrorCode CalcMu(NEP nep, PetscScalar *Mu)
235: {
237:   PetscMPIInt    sub_size,len;
238:   PetscInt       i,j,k,s;
239:   PetscScalar    *m,*temp,*temp2,*ppk,alp;
240:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
241:   Mat            M;

244:   MPI_Comm_size(PetscSubcommChild(ctx->subcomm),&sub_size);
245:   PetscMalloc3(ctx->num_solve_point*ctx->L*(ctx->L+1),&temp,2*ctx->M*ctx->L*ctx->L,&temp2,ctx->num_solve_point,&ppk);
246:   MatCreateSeqDense(PETSC_COMM_SELF,ctx->L,ctx->L_max*ctx->num_solve_point,NULL,&M);
247:   for (i=0;i<2*ctx->M*ctx->L*ctx->L;i++) temp2[i] = 0;
248:   BVSetActiveColumns(ctx->Y,0,ctx->L_max*ctx->num_solve_point);
249:   BVSetActiveColumns(ctx->V,0,ctx->L);
250:   BVDot(ctx->Y,ctx->V,M);
251:   MatDenseGetArray(M,&m);
252:   for (i=0;i<ctx->num_solve_point;i++) {
253:     for (j=0;j<ctx->L;j++) {
254:       for (k=0;k<ctx->L;k++) {
255:         temp[k+j*ctx->L+i*ctx->L*ctx->L]=m[k+j*ctx->L+i*ctx->L*ctx->L_max];
256:       }
257:     }
258:   }
259:   MatDenseRestoreArray(M,&m);
260:   for (i=0;i<ctx->num_solve_point;i++) ppk[i] = 1;
261:   for (k=0;k<2*ctx->M;k++) {
262:     for (j=0;j<ctx->L;j++) {
263:       for (i=0;i<ctx->num_solve_point;i++) {
264:         alp = ppk[i]*ctx->weight[i*ctx->subcomm->n + ctx->subcomm_id];
265:         for (s=0;s<ctx->L;s++) {
266:           if (ctx->useconj) temp2[s+(j+k*ctx->L)*ctx->L] += PetscRealPart(alp*temp[s+(j+i*ctx->L)*ctx->L])*2;
267:           else temp2[s+(j+k*ctx->L)*ctx->L] += alp*temp[s+(j+i*ctx->L)*ctx->L];
268:         }
269:       }
270:     }
271:     for (i=0;i<ctx->num_solve_point;i++)
272:       ppk[i] *= ctx->pp[i*ctx->subcomm->n + ctx->subcomm_id];
273:   }
274:   for (i=0;i<2*ctx->M*ctx->L*ctx->L;i++) temp2[i] /= sub_size;
275:   PetscMPIIntCast(2*ctx->M*ctx->L*ctx->L,&len);
276:   MPI_Allreduce(temp2,Mu,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)nep));
277:   PetscFree3(temp,temp2,ppk);
278:   MatDestroy(&M);
279:   return(0);
280: }

282: static PetscErrorCode BlockHankel(NEP nep,PetscScalar *Mu,PetscInt s,PetscScalar *H)
283: {
284:   NEP_CISS *ctx = (NEP_CISS*)nep->data;
285:   PetscInt  i,j,k,L=ctx->L,M=ctx->M;

288:   for (k=0;k<L*M;k++)
289:     for (j=0;j<M;j++)
290:       for (i=0;i<L;i++)
291:         H[j*L+i+k*L*M] = Mu[i+k*L+(j+s)*L*L];
292:   return(0);
293: }

295: static PetscErrorCode SVD_H0(NEP nep,PetscScalar *S,PetscInt *K)
296: {
297: #if defined(PETSC_MISSING_LAPACK_GESVD)
299:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable");
300: #else
302:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
303:   PetscInt       i,ml=ctx->L*ctx->M;
304:   PetscBLASInt   m,n,lda,ldu,ldvt,lwork,info;
305:   PetscScalar    *work;
306: #if defined(PETSC_USE_COMPLEX)
307:   PetscReal      *rwork;
308: #endif

311:   PetscMalloc1(5*ml,&work);
312: #if defined(PETSC_USE_COMPLEX)
313:   PetscMalloc1(5*ml,&rwork);
314: #endif
315:   PetscBLASIntCast(ml,&m);
316:   n = m; lda = m; ldu = m; ldvt = m; lwork = 5*m;
317:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
318: #if defined(PETSC_USE_COMPLEX)
319:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,S,&lda,ctx->sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,rwork,&info));
320: #else
321:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,S,&lda,ctx->sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,&info));
322: #endif
323:   SlepcCheckLapackInfo("gesvd",info);
324:   PetscFPTrapPop();
325:   (*K) = 0;
326:   for (i=0;i<ml;i++) {
327:     if (ctx->sigma[i]/PetscMax(ctx->sigma[0],1)>ctx->delta) (*K)++;
328:   }
329:   PetscFree(work);
330: #if defined(PETSC_USE_COMPLEX)
331:   PetscFree(rwork);
332: #endif
333:   return(0);
334: #endif
335: }

337: static PetscErrorCode ConstructS(NEP nep)
338: {
340:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
341:   PetscInt       i,j,k,vec_local_size,p_id;
342:   Vec            v,sj,yj;
343:   PetscScalar    *ppk, *v_data, m = 1;

346:   BVGetSizes(ctx->Y,&vec_local_size,NULL,NULL);
347:   PetscMalloc1(ctx->num_solve_point,&ppk);
348:   for (i=0;i<ctx->num_solve_point;i++) ppk[i] = 1;
349:   BVGetColumn(ctx->Y,0,&yj);
350:   VecDuplicate(yj,&v);
351:   BVRestoreColumn(ctx->Y,0,&yj);
352:   for (k=0;k<ctx->M;k++) {
353:     for (j=0;j<ctx->L;j++) {
354:       VecSet(v,0);
355:       for (i=0;i<ctx->num_solve_point;i++) {
356:         p_id = i*ctx->subcomm->n + ctx->subcomm_id;
357:         BVSetActiveColumns(ctx->Y,i*ctx->L_max+j,i*ctx->L_max+j+1);
358:         BVMultVec(ctx->Y,ppk[i]*ctx->weight[p_id],1,v,&m);
359:       }
360:       if (ctx->useconj) {
361:         VecGetArray(v,&v_data);
362:         for (i=0;i<vec_local_size;i++) v_data[i] = PetscRealPart(v_data[i])*2;
363:         VecRestoreArray(v,&v_data);
364:       }
365:       BVGetColumn(ctx->S,k*ctx->L+j,&sj);
366:       VecCopy(v,sj);
367:       BVRestoreColumn(ctx->S,k*ctx->L+j,&sj);
368:     }
369:     for (i=0;i<ctx->num_solve_point;i++) {
370:       p_id = i*ctx->subcomm->n + ctx->subcomm_id;
371:       ppk[i] *= ctx->pp[p_id];
372:     }
373:   }
374:   PetscFree(ppk);
375:   VecDestroy(&v);
376:   return(0);
377: }

379: static PetscErrorCode isGhost(NEP nep,PetscInt ld,PetscInt nv,PetscBool *fl)
380: {
382:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
383:   PetscInt       i,j;
384:   PetscScalar    *pX;
385:   PetscReal      *tau,s1,s2,tau_max=0.0;

388:   PetscMalloc1(nv,&tau);
389:   DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
390:   DSGetArray(nep->ds,DS_MAT_X,&pX);

392:   for (i=0;i<nv;i++) {
393:     s1 = 0;
394:     s2 = 0;
395:     for (j=0;j<nv;j++) {
396:       s1 += PetscAbsScalar(PetscPowScalar(pX[i*ld+j],2));
397:       s2 += PetscPowReal(PetscAbsScalar(pX[i*ld+j]),2)/ctx->sigma[j];
398:     }
399:     tau[i] = s1/s2;
400:     tau_max = PetscMax(tau_max,tau[i]);
401:   }
402:   DSRestoreArray(nep->ds,DS_MAT_X,&pX);
403:   for (i=0;i<nv;i++) {
404:     tau[i] /= tau_max;
405:   }
406:   for (i=0;i<nv;i++) {
407:     if (tau[i]>=ctx->spurious_threshold) fl[i] = PETSC_TRUE;
408:     else fl[i] = PETSC_FALSE;
409:   }
410:   PetscFree(tau);
411:   return(0);
412: }

414: PetscErrorCode NEPSetUp_CISS(NEP nep)
415: {
417:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
418:   PetscInt       i,nwork;
419:   PetscBool      istrivial,isellipse,flg;
420:   PetscScalar    center;

423:   if (!nep->ncv) nep->ncv = ctx->L_max*ctx->M;
424:   else {
425:     ctx->L_max = nep->ncv/ctx->M;
426:     if (ctx->L_max == 0) {
427:       ctx->L_max = 1;
428:       nep->ncv = ctx->L_max*ctx->M;
429:     }
430:     if (ctx->L > ctx->L_max) ctx->L = ctx->L_max;
431:   }
432:   if (!nep->max_it) nep->max_it = 1;
433:   if (!nep->mpd) nep->mpd = nep->ncv;
434:   if (!nep->which) nep->which = NEP_ALL;
435:   if (nep->stopping!=NEPStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support user-defined stopping test");

437:   /* check region */
438:   RGIsTrivial(nep->rg,&istrivial);
439:   if (istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
440:   RGGetComplement(nep->rg,&flg);
441:   if (flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
442:   PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse);
443:   if (!isellipse) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Currently only implemented for elliptic or arc regions");
444:   RGEllipseGetParameters(nep->rg,&center,NULL,NULL);
445:   if (ctx->isreal && PetscImaginaryPart(center) == 0.0) ctx->useconj = PETSC_TRUE;
446:   else ctx->useconj = PETSC_FALSE;

448:   /* create split comm */
449:   ctx->num_subcomm = 1;
450:   SetSolverComm(nep);

452:   NEPAllocateSolution(nep,0);
453:   if (ctx->weight) { PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma); }
454:   PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma);
455:   PetscLogObjectMemory((PetscObject)nep,3*ctx->N*sizeof(PetscScalar)+ctx->L_max*ctx->N*sizeof(PetscReal));

457:   /* allocate basis vectors */
458:   BVDestroy(&ctx->S);
459:   BVDuplicateResize(nep->V,ctx->L_max*ctx->M,&ctx->S);
460:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->S);
461:   BVDestroy(&ctx->V);
462:   BVDuplicateResize(nep->V,ctx->L_max,&ctx->V);
463:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->V);

465:   if (!ctx->usest) {
466:     if (ctx->ksp) {
467:       for (i=0;i<ctx->num_solve_point;i++) {
468:         KSPDestroy(&ctx->ksp[i]);
469:       }
470:       PetscFree(ctx->ksp);
471:     }
472:     PetscMalloc1(ctx->num_solve_point,&ctx->ksp);
473:     PetscLogObjectMemory((PetscObject)nep,ctx->num_solve_point*sizeof(KSP));
474:     for (i=0;i<ctx->num_solve_point;i++) {
475:       KSPCreate(PetscSubcommChild(ctx->subcomm),&ctx->ksp[i]);
476:       PetscObjectIncrementTabLevel((PetscObject)ctx->ksp[i],(PetscObject)nep,1);
477:       PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp[i]);
478:       KSPSetOptionsPrefix(ctx->ksp[i],((PetscObject)nep)->prefix);
479:       KSPAppendOptionsPrefix(ctx->ksp[i],"nep_ciss_");
480:       KSPSetErrorIfNotConverged(ctx->ksp[i],PETSC_TRUE);
481:       KSPSetTolerances(ctx->ksp[i],SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
482:     }
483:   }

485:   BVDestroy(&ctx->Y);
486:   BVDuplicateResize(nep->V,ctx->num_solve_point*ctx->L_max,&ctx->Y);

488:   DSSetType(nep->ds,DSGNHEP);
489:   DSAllocate(nep->ds,nep->ncv);
490:   nwork = (nep->fui==NEP_USER_INTERFACE_SPLIT)? 2: 1;
491:   NEPSetWorkVecs(nep,nwork);
492:   return(0);
493: }

495: PetscErrorCode NEPSolve_CISS(NEP nep)
496: {
498:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
499:   Mat            X,M;
500:   PetscInt       i,j,ld,L_add=0,nv=0,L_base=ctx->L,inner,*inside;
501:   PetscScalar    *Mu,*H0,*H1,*rr,*temp,center;
502:   PetscReal      error,max_error,radius,rgscale;
503:   PetscBool      *fl1;
504:   Vec            si;
505:   SlepcSC        sc;
506:   PetscRandom    rand;

509:   DSGetSlepcSC(nep->ds,&sc);
510:   sc->comparison    = SlepcCompareLargestMagnitude;
511:   sc->comparisonctx = NULL;
512:   sc->map           = NULL;
513:   sc->mapobj        = NULL;
514:   DSGetLeadingDimension(nep->ds,&ld);
515:   SetPathParameter(nep);
516:   CISSVecSetRandom(ctx->V,0,ctx->L);
517:   BVGetRandomContext(ctx->V,&rand);

519:   SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,0,ctx->L,PETSC_TRUE);
520:   EstimateNumberEigs(nep,&L_add);
521:   if (L_add>0) {
522:     PetscInfo2(nep,"Changing L %D -> %D by Estimate #Eig\n",ctx->L,ctx->L+L_add);
523:     CISSVecSetRandom(ctx->V,ctx->L,ctx->L+L_add);
524:     SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
525:     ctx->L += L_add;
526:   }
527:   PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
528:   for (i=0;i<ctx->refine_blocksize;i++) {
529:     CalcMu(nep,Mu);
530:     BlockHankel(nep,Mu,0,H0);
531:     SVD_H0(nep,H0,&nv);
532:     if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
533:     L_add = L_base;
534:     if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
535:     PetscInfo2(nep,"Changing L %D -> %D by SVD(H0)\n",ctx->L,ctx->L+L_add);
536:     CISSVecSetRandom(ctx->V,ctx->L,ctx->L+L_add);
537:     SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
538:     ctx->L += L_add;
539:   }
540:   PetscFree2(Mu,H0);

542:   RGGetScale(nep->rg,&rgscale);
543:   RGEllipseGetParameters(nep->rg,&center,&radius,NULL);

545:   PetscMalloc3(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0,ctx->L*ctx->M*ctx->L*ctx->M,&H1);
546:   while (nep->reason == NEP_CONVERGED_ITERATING) {
547:     nep->its++;
548:     for (inner=0;inner<=ctx->refine_inner;inner++) {
549:       CalcMu(nep,Mu);
550:       BlockHankel(nep,Mu,0,H0);
551:       SVD_H0(nep,H0,&nv);
552:       if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
553:         ConstructS(nep);
554:         BVSetActiveColumns(ctx->S,0,ctx->L);
555:         BVCopy(ctx->S,ctx->V);
556:         SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,0,ctx->L,PETSC_FALSE);
557:       } else break;
558:     }

560:     nep->nconv = 0;
561:     if (nv == 0) break;
562:     BlockHankel(nep,Mu,0,H0);
563:     BlockHankel(nep,Mu,1,H1);
564:     DSSetDimensions(nep->ds,nv,0,0,0);
565:     DSSetState(nep->ds,DS_STATE_RAW);
566:     DSGetArray(nep->ds,DS_MAT_A,&temp);
567:     for (j=0;j<nv;j++)
568:       for (i=0;i<nv;i++)
569:         temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
570:     DSRestoreArray(nep->ds,DS_MAT_A,&temp);
571:     DSGetArray(nep->ds,DS_MAT_B,&temp);
572:     for (j=0;j<nv;j++)
573:       for (i=0;i<nv;i++)
574:         temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
575:     DSRestoreArray(nep->ds,DS_MAT_B,&temp);
576:     DSSolve(nep->ds,nep->eigr,nep->eigi);
577:     DSSynchronize(nep->ds,nep->eigr,nep->eigi);
578:     DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
579:     for (i=0;i<nv;i++){
580:       nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
581: #if !defined(PETSC_USE_COMPLEX)
582:       nep->eigi[i] = nep->eigi[i]*radius*rgscale;
583: #endif
584:     }
585:     PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr);
586:     isGhost(nep,ld,nv,fl1);
587:     RGCheckInside(nep->rg,nv,nep->eigr,nep->eigi,inside);
588:     for (i=0;i<nv;i++) {
589:       if (fl1[i] && inside[i]>0) {
590:         rr[i] = 1.0;
591:         nep->nconv++;
592:       } else rr[i] = 0.0;
593:     }
594:     DSSort(nep->ds,nep->eigr,nep->eigi,rr,NULL,&nep->nconv);
595:     DSSynchronize(nep->ds,nep->eigr,nep->eigi);
596:     for (i=0;i<nv;i++){
597:       nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
598: #if !defined(PETSC_USE_COMPLEX)
599:       nep->eigi[i] = nep->eigi[i]*radius*rgscale;
600: #endif
601:     }
602:     PetscFree3(fl1,inside,rr);
603:     BVSetActiveColumns(nep->V,0,nv);
604:     ConstructS(nep);
605:     BVSetActiveColumns(ctx->S,0,nv);
606:     BVCopy(ctx->S,nep->V);

608:     DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
609:     DSGetMat(nep->ds,DS_MAT_X,&X);
610:     BVMultInPlace(ctx->S,X,0,nep->nconv);
611:     BVMultInPlace(nep->V,X,0,nep->nconv);
612:     MatDestroy(&X);
613:     max_error = 0.0;
614:     for (i=0;i<nep->nconv;i++) {
615:       BVGetColumn(nep->V,i,&si);
616:       VecNormalize(si,NULL);
617:       NEPComputeResidualNorm_Private(nep,nep->eigr[i],si,nep->work,&error);
618:       (*nep->converged)(nep,nep->eigr[i],0,error,&error,nep->convergedctx);
619:       BVRestoreColumn(nep->V,i,&si);
620:       max_error = PetscMax(max_error,error);
621:     }
622:     if (max_error <= nep->tol) nep->reason = NEP_CONVERGED_TOL;
623:     else if (nep->its > nep->max_it) nep->reason = NEP_DIVERGED_ITS;
624:     else {
625:       if (nep->nconv > ctx->L) nv = nep->nconv;
626:       else if (ctx->L > nv) nv = ctx->L;
627:       MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M);
628:       MatDenseGetArray(M,&temp);
629:       for (i=0;i<ctx->L*nv;i++) {
630:         PetscRandomGetValue(rand,&temp[i]);
631:         temp[i] = PetscRealPart(temp[i]);
632:       }
633:       MatDenseRestoreArray(M,&temp);
634:       BVSetActiveColumns(ctx->S,0,nv);
635:       BVMultInPlace(ctx->S,M,0,ctx->L);
636:       MatDestroy(&M);
637:       BVSetActiveColumns(ctx->S,0,ctx->L);
638:       BVCopy(ctx->S,ctx->V);
639:       SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,0,ctx->L,PETSC_FALSE);
640:     }
641:   }
642:   PetscFree3(Mu,H0,H1);
643:   return(0);
644: }

646: static PetscErrorCode NEPCISSSetSizes_CISS(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
647: {
648:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

651:   if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
652:     if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
653:   } else {
654:     if (ip<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
655:     if (ip%2) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
656:     if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
657:   }
658:   if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
659:     ctx->L = 16;
660:   } else {
661:     if (bs<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
662:     if (bs>ctx->L_max) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be less than or equal to the maximum number of block size");
663:     ctx->L = bs;
664:   }
665:   if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
666:     ctx->M = ctx->N/4;
667:   } else {
668:     if (ms<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
669:     if (ms>ctx->N) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be less than or equal to the number of integration points");
670:     ctx->M = ms;
671:   }
672:   if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
673:     ctx->num_subcomm = 1;
674:   } else {
675:     if (npart<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The npart argument must be > 0");
676:     ctx->num_subcomm = npart;
677:   }
678:   if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
679:     ctx->L = 256;
680:   } else {
681:     if (bsmax<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
682:     if (bsmax<ctx->L) ctx->L_max = ctx->L;
683:     else ctx->L_max = bsmax;
684:   }
685:   ctx->isreal = realmats;
686:   nep->state = NEP_STATE_INITIAL;
687:   return(0);
688: }

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

693:    Logically Collective on NEP

695:    Input Parameters:
696: +  nep   - the eigenproblem solver context
697: .  ip    - number of integration points
698: .  bs    - block size
699: .  ms    - moment size
700: .  npart - number of partitions when splitting the communicator
701: .  bsmax - max block size
702: -  realmats - T(z) is real for real z

704:    Options Database Keys:
705: +  -nep_ciss_integration_points - Sets the number of integration points
706: .  -nep_ciss_blocksize - Sets the block size
707: .  -nep_ciss_moments - Sets the moment size
708: .  -nep_ciss_partitions - Sets the number of partitions
709: .  -nep_ciss_maxblocksize - Sets the maximum block size
710: -  -nep_ciss_realmats - T(z) is real for real z

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

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

721:    Level: advanced

723: .seealso: NEPCISSGetSizes()
724: @*/
725: PetscErrorCode NEPCISSSetSizes(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
726: {

737:   PetscTryMethod(nep,"NEPCISSSetSizes_C",(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(nep,ip,bs,ms,npart,bsmax,realmats));
738:   return(0);
739: }

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

746:   if (ip) *ip = ctx->N;
747:   if (bs) *bs = ctx->L;
748:   if (ms) *ms = ctx->M;
749:   if (npart) *npart = ctx->num_subcomm;
750:   if (bsmax) *bsmax = ctx->L_max;
751:   if (realmats) *realmats = ctx->isreal;
752:   return(0);
753: }

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

758:    Not Collective

760:    Input Parameter:
761: .  nep - the eigenproblem solver context

763:    Output Parameters:
764: +  ip    - number of integration points
765: .  bs    - block size
766: .  ms    - moment size
767: .  npart - number of partitions when splitting the communicator
768: .  bsmax - max block size
769: -  realmats - T(z) is real for real z

771:    Level: advanced

773: .seealso: NEPCISSSetSizes()
774: @*/
775: PetscErrorCode NEPCISSGetSizes(NEP nep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
776: {

781:   PetscUseMethod(nep,"NEPCISSGetSizes_C",(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(nep,ip,bs,ms,npart,bsmax,realmats));
782:   return(0);
783: }

785: static PetscErrorCode NEPCISSSetThreshold_CISS(NEP nep,PetscReal delta,PetscReal spur)
786: {
787:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

790:   if (delta == PETSC_DEFAULT) {
791:     ctx->delta = 1e-12;
792:   } else {
793:     if (delta<=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
794:     ctx->delta = delta;
795:   }
796:   if (spur == PETSC_DEFAULT) {
797:     ctx->spurious_threshold = 1e-4;
798:   } else {
799:     if (spur<=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
800:     ctx->spurious_threshold = spur;
801:   }
802:   return(0);
803: }

805: /*@
806:    NEPCISSSetThreshold - Sets the values of various threshold parameters in
807:    the CISS solver.

809:    Logically Collective on NEP

811:    Input Parameters:
812: +  nep   - the eigenproblem solver context
813: .  delta - threshold for numerical rank
814: -  spur  - spurious threshold (to discard spurious eigenpairs)

816:    Options Database Keys:
817: +  -nep_ciss_delta - Sets the delta
818: -  -nep_ciss_spurious_threshold - Sets the spurious threshold

820:    Level: advanced

822: .seealso: NEPCISSGetThreshold()
823: @*/
824: PetscErrorCode NEPCISSSetThreshold(NEP nep,PetscReal delta,PetscReal spur)
825: {

832:   PetscTryMethod(nep,"NEPCISSSetThreshold_C",(NEP,PetscReal,PetscReal),(nep,delta,spur));
833:   return(0);
834: }

836: static PetscErrorCode NEPCISSGetThreshold_CISS(NEP nep,PetscReal *delta,PetscReal *spur)
837: {
838:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

841:   if (delta) *delta = ctx->delta;
842:   if (spur)  *spur = ctx->spurious_threshold;
843:   return(0);
844: }

846: /*@
847:    NEPCISSGetThreshold - Gets the values of various threshold parameters in
848:    the CISS solver.

850:    Not Collective

852:    Input Parameter:
853: .  nep - the eigenproblem solver context

855:    Output Parameters:
856: +  delta - threshold for numerical rank
857: -  spur  - spurious threshold (to discard spurious eigenpairs)

859:    Level: advanced

861: .seealso: NEPCISSSetThreshold()
862: @*/
863: PetscErrorCode NEPCISSGetThreshold(NEP nep,PetscReal *delta,PetscReal *spur)
864: {

869:   PetscUseMethod(nep,"NEPCISSGetThreshold_C",(NEP,PetscReal*,PetscReal*),(nep,delta,spur));
870:   return(0);
871: }

873: static PetscErrorCode NEPCISSSetRefinement_CISS(NEP nep,PetscInt inner,PetscInt blsize)
874: {
875:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

878:   if (inner == PETSC_DEFAULT) {
879:     ctx->refine_inner = 0;
880:   } else {
881:     if (inner<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
882:     ctx->refine_inner = inner;
883:   }
884:   if (blsize == PETSC_DEFAULT) {
885:     ctx->refine_blocksize = 0;
886:   } else {
887:     if (blsize<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
888:     ctx->refine_blocksize = blsize;
889:   }
890:   return(0);
891: }

893: /*@
894:    NEPCISSSetRefinement - Sets the values of various refinement parameters
895:    in the CISS solver.

897:    Logically Collective on NEP

899:    Input Parameters:
900: +  nep    - the eigenproblem solver context
901: .  inner  - number of iterative refinement iterations (inner loop)
902: -  blsize - number of iterative refinement iterations (blocksize loop)

904:    Options Database Keys:
905: +  -nep_ciss_refine_inner - Sets number of inner iterations
906: -  -nep_ciss_refine_blocksize - Sets number of blocksize iterations

908:    Level: advanced

910: .seealso: NEPCISSGetRefinement()
911: @*/
912: PetscErrorCode NEPCISSSetRefinement(NEP nep,PetscInt inner,PetscInt blsize)
913: {

920:   PetscTryMethod(nep,"NEPCISSSetRefinement_C",(NEP,PetscInt,PetscInt),(nep,inner,blsize));
921:   return(0);
922: }

924: static PetscErrorCode NEPCISSGetRefinement_CISS(NEP nep,PetscInt *inner,PetscInt *blsize)
925: {
926:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

929:   if (inner)  *inner = ctx->refine_inner;
930:   if (blsize) *blsize = ctx->refine_blocksize;
931:   return(0);
932: }

934: /*@
935:    NEPCISSGetRefinement - Gets the values of various refinement parameters
936:    in the CISS solver.

938:    Not Collective

940:    Input Parameter:
941: .  nep - the eigenproblem solver context

943:    Output Parameters:
944: +  inner  - number of iterative refinement iterations (inner loop)
945: -  blsize - number of iterative refinement iterations (blocksize loop)

947:    Level: advanced

949: .seealso: NEPCISSSetRefinement()
950: @*/
951: PetscErrorCode NEPCISSGetRefinement(NEP nep, PetscInt *inner, PetscInt *blsize)
952: {

957:   PetscUseMethod(nep,"NEPCISSGetRefinement_C",(NEP,PetscInt*,PetscInt*),(nep,inner,blsize));
958:   return(0);
959: }

961: PetscErrorCode NEPReset_CISS(NEP nep)
962: {
964:   PetscInt       i;
965:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;

968:   BVDestroy(&ctx->S);
969:   BVDestroy(&ctx->V);
970:   BVDestroy(&ctx->Y);
971:   if (!ctx->usest) {
972:     for (i=0;i<ctx->num_solve_point;i++) {
973:       KSPDestroy(&ctx->ksp[i]);
974:     }
975:     PetscFree(ctx->ksp);
976:   }
977:   return(0);
978: }

980: PetscErrorCode NEPSetFromOptions_CISS(PetscOptionItems *PetscOptionsObject,NEP nep)
981: {
983:   PetscReal      r1,r2;
984:   PetscInt       i1,i2,i3,i4,i5,i6,i7;
985:   PetscBool      b1;

988:   PetscOptionsHead(PetscOptionsObject,"NEP CISS Options");

990:     NEPCISSGetSizes(nep,&i1,&i2,&i3,&i4,&i5,&b1);
991:     PetscOptionsInt("-nep_ciss_integration_points","Number of integration points","NEPCISSSetSizes",i1,&i1,NULL);
992:     PetscOptionsInt("-nep_ciss_blocksize","Block size","NEPCISSSetSizes",i2,&i2,NULL);
993:     PetscOptionsInt("-nep_ciss_moments","Moment size","NEPCISSSetSizes",i3,&i3,NULL);
994:     PetscOptionsInt("-nep_ciss_partitions","Number of partitions","NEPCISSSetSizes",i4,&i4,NULL);
995:     PetscOptionsInt("-nep_ciss_maxblocksize","Maximum block size","NEPCISSSetSizes",i5,&i5,NULL);
996:     PetscOptionsBool("-nep_ciss_realmats","True if T(z) is real for real z","NEPCISSSetSizes",b1,&b1,NULL);
997:     NEPCISSSetSizes(nep,i1,i2,i3,i4,i5,b1);

999:     NEPCISSGetThreshold(nep,&r1,&r2);
1000:     PetscOptionsReal("-nep_ciss_delta","Threshold for numerical rank","NEPCISSSetThreshold",r1,&r1,NULL);
1001:     PetscOptionsReal("-nep_ciss_spurious_threshold","Threshold for the spurious eigenpairs","NEPCISSSetThreshold",r2,&r2,NULL);
1002:     NEPCISSSetThreshold(nep,r1,r2);

1004:     NEPCISSGetRefinement(nep,&i6,&i7);
1005:     PetscOptionsInt("-nep_ciss_refine_inner","Number of inner iterative refinement iterations","NEPCISSSetRefinement",i6,&i6,NULL);
1006:     PetscOptionsInt("-nep_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","NEPCISSSetRefinement",i7,&i7,NULL);
1007:     NEPCISSSetRefinement(nep,i6,i7);

1009:   PetscOptionsTail();
1010:   return(0);
1011: }

1013: PetscErrorCode NEPDestroy_CISS(NEP nep)
1014: {
1016:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;

1019:   PetscSubcommDestroy(&ctx->subcomm);
1020:   PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
1021:   PetscFree(nep->data);
1022:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NULL);
1023:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NULL);
1024:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NULL);
1025:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NULL);
1026:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NULL);
1027:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NULL);
1028:   return(0);
1029: }

1031: PetscErrorCode NEPView_CISS(NEP nep,PetscViewer viewer)
1032: {
1034:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
1035:   PetscBool      isascii;

1038:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1039:   if (isascii) {
1040:     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->num_subcomm,ctx->L_max);
1041:     if (ctx->isreal) {
1042:       PetscViewerASCIIPrintf(viewer,"  exploiting symmetry of integration points\n");
1043:     }
1044:     PetscViewerASCIIPrintf(viewer,"  threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold);
1045:     PetscViewerASCIIPrintf(viewer,"  iterative refinement  { inner: %D, blocksize: %D }\n",ctx->refine_inner, ctx->refine_blocksize);
1046:     PetscViewerASCIIPushTab(viewer);
1047:     if (!ctx->usest && ctx->ksp[0]) { KSPView(ctx->ksp[0],viewer); }
1048:     PetscViewerASCIIPopTab(viewer);
1049:   }
1050:   return(0);
1051: }

1053: PETSC_EXTERN PetscErrorCode NEPCreate_CISS(NEP nep)
1054: {
1056:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;

1059:   PetscNewLog(nep,&ctx);
1060:   nep->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              = 1e-12;
1066:   ctx->L_max              = 64;
1067:   ctx->spurious_threshold = 1e-4;
1068:   ctx->usest              = PETSC_FALSE;
1069:   ctx->isreal             = PETSC_FALSE;
1070:   ctx->num_subcomm        = 1;

1072:   nep->ops->solve          = NEPSolve_CISS;
1073:   nep->ops->setup          = NEPSetUp_CISS;
1074:   nep->ops->setfromoptions = NEPSetFromOptions_CISS;
1075:   nep->ops->reset          = NEPReset_CISS;
1076:   nep->ops->destroy        = NEPDestroy_CISS;
1077:   nep->ops->view           = NEPView_CISS;

1079:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NEPCISSSetSizes_CISS);
1080:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NEPCISSGetSizes_CISS);
1081:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NEPCISSSetThreshold_CISS);
1082:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NEPCISSGetThreshold_CISS);
1083:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NEPCISSSetRefinement_CISS);
1084:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NEPCISSGetRefinement_CISS);
1085:   return(0);
1086: }