Actual source code: nciss.c
slepc-3.9.1 2018-05-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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 npart; /* number of partitions */
46: PetscInt refine_inner;
47: PetscInt refine_blocksize;
48: /* private data */
49: PetscReal *sigma; /* threshold for numerical rank */
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; /* ksp array for storing factorizations at integration points */
59: PetscBool useconj;
60: PetscReal est_eig;
61: PetscSubcomm subcomm;
62: } NEP_CISS;
64: /* destroy KSP objects when the number of solve points changes */
65: PETSC_STATIC_INLINE PetscErrorCode NEPCISSResetSolvers(NEP nep)
66: {
68: PetscInt i;
69: NEP_CISS *ctx = (NEP_CISS*)nep->data;
72: if (ctx->ksp) {
73: for (i=0;i<ctx->num_solve_point;i++) {
74: KSPDestroy(&ctx->ksp[i]);
75: }
76: PetscFree(ctx->ksp);
77: }
78: return(0);
79: }
81: /* clean PetscSubcomm object when the number of partitions changes */
82: PETSC_STATIC_INLINE PetscErrorCode NEPCISSResetSubcomm(NEP nep)
83: {
85: NEP_CISS *ctx = (NEP_CISS*)nep->data;
88: NEPCISSResetSolvers(nep);
89: PetscSubcommDestroy(&ctx->subcomm);
90: return(0);
91: }
93: /* determine whether half of integration points can be avoided (use its conjugates);
94: depends on isreal and the center of the region */
95: PETSC_STATIC_INLINE PetscErrorCode NEPCISSSetUseConj(NEP nep,PetscBool *useconj)
96: {
98: NEP_CISS *ctx = (NEP_CISS*)nep->data;
99: PetscScalar center;
102: *useconj = PETSC_FALSE;
103: RGEllipseGetParameters(nep->rg,¢er,NULL,NULL);
104: *useconj = (ctx->isreal && PetscImaginaryPart(center) == 0.0)? PETSC_TRUE: PETSC_FALSE;
105: return(0);
106: }
108: /* create PetscSubcomm object and determine num_solve_point (depends on npart, N, useconj) */
109: PETSC_STATIC_INLINE PetscErrorCode NEPCISSSetUpSubComm(NEP nep,PetscInt *num_solve_point)
110: {
112: NEP_CISS *ctx = (NEP_CISS*)nep->data;
113: PetscInt N = ctx->N;
116: PetscSubcommCreate(PetscObjectComm((PetscObject)nep),&ctx->subcomm);
117: PetscSubcommSetNumber(ctx->subcomm,ctx->npart);
118: PetscSubcommSetType(ctx->subcomm,PETSC_SUBCOMM_INTERLACED);
119: PetscLogObjectMemory((PetscObject)nep,sizeof(PetscSubcomm));
120: ctx->subcomm_id = ctx->subcomm->color;
121: NEPCISSSetUseConj(nep,&ctx->useconj);
122: if (ctx->useconj) N = N/2;
123: *num_solve_point = N / ctx->npart;
124: if (N%ctx->npart > ctx->subcomm_id) (*num_solve_point)++;
125: return(0);
126: }
128: static PetscErrorCode SetPathParameter(NEP nep)
129: {
131: NEP_CISS *ctx = (NEP_CISS*)nep->data;
132: PetscInt i;
133: PetscScalar center;
134: PetscReal theta,radius,vscale,rgscale;
135: PetscBool isellipse=PETSC_FALSE;
138: PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse);
139: RGGetScale(nep->rg,&rgscale);
140: if (isellipse) {
141: RGEllipseGetParameters(nep->rg,¢er,&radius,&vscale);
142: } else SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Region must be Ellipse");
143: for (i=0;i<ctx->N;i++) {
144: theta = ((2*PETSC_PI)/ctx->N)*(i+0.5);
145: ctx->pp[i] = PetscCosReal(theta) + PETSC_i*vscale*PetscSinReal(theta);
146: ctx->weight[i] = radius*(vscale*PetscCosReal(theta) + PETSC_i*PetscSinReal(theta))/(PetscReal)ctx->N;
147: ctx->omega[i] = rgscale*(center + radius*ctx->pp[i]);
148: }
149: return(0);
150: }
152: static PetscErrorCode CISSVecSetRandom(BV V,PetscInt i0,PetscInt i1)
153: {
155: PetscInt i,j,nlocal;
156: PetscScalar *vdata;
157: Vec x;
160: BVGetSizes(V,&nlocal,NULL,NULL);
161: for (i=i0;i<i1;i++) {
162: BVSetRandomColumn(V,i);
163: BVGetColumn(V,i,&x);
164: VecGetArray(x,&vdata);
165: for (j=0;j<nlocal;j++) {
166: vdata[j] = PetscRealPart(vdata[j]);
167: if (PetscRealPart(vdata[j]) < 0.5) vdata[j] = -1.0;
168: else vdata[j] = 1.0;
169: }
170: VecRestoreArray(x,&vdata);
171: BVRestoreColumn(V,i,&x);
172: }
173: return(0);
174: }
176: static PetscErrorCode SolveLinearSystem(NEP nep,Mat T,Mat dT,BV V,PetscInt L_start,PetscInt L_end,PetscBool initksp)
177: {
179: NEP_CISS *ctx = (NEP_CISS*)nep->data;
180: PetscInt i,j,p_id;
181: Mat kspMat;
182: Vec Bvj,vj,yj;
185: if (!ctx->ksp) { NEPCISSGetKSPs(nep,&ctx->num_solve_point,&ctx->ksp); }
186: BVCreateVec(V,&Bvj);
187: for (i=0;i<ctx->num_solve_point;i++) {
188: p_id = i*ctx->subcomm->n + ctx->subcomm_id;
189: if (initksp) {
190: NEPComputeFunction(nep,ctx->omega[p_id],T,T);
191: MatDuplicate(T,MAT_COPY_VALUES,&kspMat);
192: KSPSetOperators(ctx->ksp[i],kspMat,kspMat);
193: MatDestroy(&kspMat);
194: }
195: NEPComputeJacobian(nep,ctx->omega[p_id],dT);
196: for (j=L_start;j<L_end;j++) {
197: BVGetColumn(V,j,&vj);
198: BVGetColumn(ctx->Y,i*ctx->L_max+j,&yj);
199: MatMult(dT,vj,Bvj);
200: KSPSolve(ctx->ksp[i],Bvj,yj);
201: BVRestoreColumn(V,j,&vj);
202: BVRestoreColumn(ctx->Y,i*ctx->L_max+j,&yj);
203: }
204: }
205: VecDestroy(&Bvj);
206: return(0);
207: }
209: static PetscErrorCode EstimateNumberEigs(NEP nep,PetscInt *L_add)
210: {
212: NEP_CISS *ctx = (NEP_CISS*)nep->data;
213: PetscInt i,j,p_id;
214: PetscScalar tmp,m = 1,sum = 0.0;
215: PetscReal eta;
216: Vec v,vtemp,vj,yj;
219: BVGetColumn(ctx->Y,0,&yj);
220: VecDuplicate(yj,&v);
221: BVRestoreColumn(ctx->Y,0,&yj);
222: BVCreateVec(ctx->V,&vtemp);
223: for (j=0;j<ctx->L;j++) {
224: VecSet(v,0);
225: for (i=0;i<ctx->num_solve_point; i++) {
226: p_id = i*ctx->subcomm->n + ctx->subcomm_id;
227: BVSetActiveColumns(ctx->Y,i*ctx->L_max+j,i*ctx->L_max+j+1);
228: BVMultVec(ctx->Y,ctx->weight[p_id],1,v,&m);
229: }
230: BVGetColumn(ctx->V,j,&vj);
231: VecDot(vj,v,&tmp);
232: BVRestoreColumn(ctx->V,j,&vj);
233: if (ctx->useconj) sum += PetscRealPart(tmp)*2;
234: else sum += tmp;
235: }
236: ctx->est_eig = PetscAbsScalar(sum/(PetscReal)ctx->L);
237: eta = PetscPowReal(10,-PetscLog10Real(nep->tol)/ctx->N);
238: PetscInfo1(nep,"Estimation_#Eig %f\n",(double)ctx->est_eig);
239: *L_add = (PetscInt)PetscCeilReal((ctx->est_eig*eta)/ctx->M) - ctx->L;
240: if (*L_add < 0) *L_add = 0;
241: if (*L_add>ctx->L_max-ctx->L) {
242: PetscInfo(nep,"Number of eigenvalues around the contour path may be too large\n");
243: *L_add = ctx->L_max-ctx->L;
244: }
245: VecDestroy(&v);
246: VecDestroy(&vtemp);
247: return(0);
248: }
250: static PetscErrorCode CalcMu(NEP nep, PetscScalar *Mu)
251: {
253: PetscMPIInt sub_size,len;
254: PetscInt i,j,k,s;
255: PetscScalar *m,*temp,*temp2,*ppk,alp;
256: NEP_CISS *ctx = (NEP_CISS*)nep->data;
257: Mat M;
260: MPI_Comm_size(PetscSubcommChild(ctx->subcomm),&sub_size);
261: PetscMalloc3(ctx->num_solve_point*ctx->L*(ctx->L+1),&temp,2*ctx->M*ctx->L*ctx->L,&temp2,ctx->num_solve_point,&ppk);
262: MatCreateSeqDense(PETSC_COMM_SELF,ctx->L,ctx->L_max*ctx->num_solve_point,NULL,&M);
263: for (i=0;i<2*ctx->M*ctx->L*ctx->L;i++) temp2[i] = 0;
264: BVSetActiveColumns(ctx->Y,0,ctx->L_max*ctx->num_solve_point);
265: BVSetActiveColumns(ctx->V,0,ctx->L);
266: BVDot(ctx->Y,ctx->V,M);
267: MatDenseGetArray(M,&m);
268: for (i=0;i<ctx->num_solve_point;i++) {
269: for (j=0;j<ctx->L;j++) {
270: for (k=0;k<ctx->L;k++) {
271: temp[k+j*ctx->L+i*ctx->L*ctx->L]=m[k+j*ctx->L+i*ctx->L*ctx->L_max];
272: }
273: }
274: }
275: MatDenseRestoreArray(M,&m);
276: for (i=0;i<ctx->num_solve_point;i++) ppk[i] = 1;
277: for (k=0;k<2*ctx->M;k++) {
278: for (j=0;j<ctx->L;j++) {
279: for (i=0;i<ctx->num_solve_point;i++) {
280: alp = ppk[i]*ctx->weight[i*ctx->subcomm->n + ctx->subcomm_id];
281: for (s=0;s<ctx->L;s++) {
282: if (ctx->useconj) temp2[s+(j+k*ctx->L)*ctx->L] += PetscRealPart(alp*temp[s+(j+i*ctx->L)*ctx->L])*2;
283: else temp2[s+(j+k*ctx->L)*ctx->L] += alp*temp[s+(j+i*ctx->L)*ctx->L];
284: }
285: }
286: }
287: for (i=0;i<ctx->num_solve_point;i++)
288: ppk[i] *= ctx->pp[i*ctx->subcomm->n + ctx->subcomm_id];
289: }
290: for (i=0;i<2*ctx->M*ctx->L*ctx->L;i++) temp2[i] /= sub_size;
291: PetscMPIIntCast(2*ctx->M*ctx->L*ctx->L,&len);
292: MPI_Allreduce(temp2,Mu,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)nep));
293: PetscFree3(temp,temp2,ppk);
294: MatDestroy(&M);
295: return(0);
296: }
298: static PetscErrorCode BlockHankel(NEP nep,PetscScalar *Mu,PetscInt s,PetscScalar *H)
299: {
300: NEP_CISS *ctx = (NEP_CISS*)nep->data;
301: PetscInt i,j,k,L=ctx->L,M=ctx->M;
304: for (k=0;k<L*M;k++)
305: for (j=0;j<M;j++)
306: for (i=0;i<L;i++)
307: H[j*L+i+k*L*M] = Mu[i+k*L+(j+s)*L*L];
308: return(0);
309: }
311: static PetscErrorCode SVD_H0(NEP nep,PetscScalar *S,PetscInt *K)
312: {
313: #if defined(PETSC_MISSING_LAPACK_GESVD)
315: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable");
316: #else
318: NEP_CISS *ctx = (NEP_CISS*)nep->data;
319: PetscInt i,ml=ctx->L*ctx->M;
320: PetscBLASInt m,n,lda,ldu,ldvt,lwork,info;
321: PetscScalar *work;
322: #if defined(PETSC_USE_COMPLEX)
323: PetscReal *rwork;
324: #endif
327: PetscMalloc1(5*ml,&work);
328: #if defined(PETSC_USE_COMPLEX)
329: PetscMalloc1(5*ml,&rwork);
330: #endif
331: PetscBLASIntCast(ml,&m);
332: n = m; lda = m; ldu = m; ldvt = m; lwork = 5*m;
333: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
334: #if defined(PETSC_USE_COMPLEX)
335: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,S,&lda,ctx->sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,rwork,&info));
336: #else
337: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,S,&lda,ctx->sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,&info));
338: #endif
339: SlepcCheckLapackInfo("gesvd",info);
340: PetscFPTrapPop();
341: (*K) = 0;
342: for (i=0;i<ml;i++) {
343: if (ctx->sigma[i]/PetscMax(ctx->sigma[0],1)>ctx->delta) (*K)++;
344: }
345: PetscFree(work);
346: #if defined(PETSC_USE_COMPLEX)
347: PetscFree(rwork);
348: #endif
349: return(0);
350: #endif
351: }
353: static PetscErrorCode ConstructS(NEP nep)
354: {
356: NEP_CISS *ctx = (NEP_CISS*)nep->data;
357: PetscInt i,j,k,vec_local_size,p_id;
358: Vec v,sj,yj;
359: PetscScalar *ppk, *v_data, m = 1;
362: BVGetSizes(ctx->Y,&vec_local_size,NULL,NULL);
363: PetscMalloc1(ctx->num_solve_point,&ppk);
364: for (i=0;i<ctx->num_solve_point;i++) ppk[i] = 1;
365: BVGetColumn(ctx->Y,0,&yj);
366: VecDuplicate(yj,&v);
367: BVRestoreColumn(ctx->Y,0,&yj);
368: for (k=0;k<ctx->M;k++) {
369: for (j=0;j<ctx->L;j++) {
370: VecSet(v,0);
371: for (i=0;i<ctx->num_solve_point;i++) {
372: p_id = i*ctx->subcomm->n + ctx->subcomm_id;
373: BVSetActiveColumns(ctx->Y,i*ctx->L_max+j,i*ctx->L_max+j+1);
374: BVMultVec(ctx->Y,ppk[i]*ctx->weight[p_id],1,v,&m);
375: }
376: if (ctx->useconj) {
377: VecGetArray(v,&v_data);
378: for (i=0;i<vec_local_size;i++) v_data[i] = PetscRealPart(v_data[i])*2;
379: VecRestoreArray(v,&v_data);
380: }
381: BVGetColumn(ctx->S,k*ctx->L+j,&sj);
382: VecCopy(v,sj);
383: BVRestoreColumn(ctx->S,k*ctx->L+j,&sj);
384: }
385: for (i=0;i<ctx->num_solve_point;i++) {
386: p_id = i*ctx->subcomm->n + ctx->subcomm_id;
387: ppk[i] *= ctx->pp[p_id];
388: }
389: }
390: PetscFree(ppk);
391: VecDestroy(&v);
392: return(0);
393: }
395: static PetscErrorCode isGhost(NEP nep,PetscInt ld,PetscInt nv,PetscBool *fl)
396: {
398: NEP_CISS *ctx = (NEP_CISS*)nep->data;
399: PetscInt i,j;
400: PetscScalar *pX;
401: PetscReal *tau,s1,s2,tau_max=0.0;
404: PetscMalloc1(nv,&tau);
405: DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
406: DSGetArray(nep->ds,DS_MAT_X,&pX);
408: for (i=0;i<nv;i++) {
409: s1 = 0;
410: s2 = 0;
411: for (j=0;j<nv;j++) {
412: s1 += PetscAbsScalar(PetscPowScalar(pX[i*ld+j],2));
413: s2 += PetscPowReal(PetscAbsScalar(pX[i*ld+j]),2)/ctx->sigma[j];
414: }
415: tau[i] = s1/s2;
416: tau_max = PetscMax(tau_max,tau[i]);
417: }
418: DSRestoreArray(nep->ds,DS_MAT_X,&pX);
419: for (i=0;i<nv;i++) {
420: tau[i] /= tau_max;
421: }
422: for (i=0;i<nv;i++) {
423: if (tau[i]>=ctx->spurious_threshold) fl[i] = PETSC_TRUE;
424: else fl[i] = PETSC_FALSE;
425: }
426: PetscFree(tau);
427: return(0);
428: }
430: PetscErrorCode NEPSetUp_CISS(NEP nep)
431: {
433: NEP_CISS *ctx = (NEP_CISS*)nep->data;
434: PetscInt nwork;
435: PetscBool istrivial,isellipse,flg,useconj;
438: if (!nep->ncv) {
439: nep->ncv = ctx->L_max*ctx->M;
440: if (nep->ncv>nep->n) {
441: nep->ncv = nep->n;
442: ctx->L_max = nep->ncv/ctx->M;
443: if (!ctx->L_max) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Cannot adjust solver parameters, try setting a smaller value of M (moment size)");
444: }
445: } else {
446: ctx->L_max = nep->ncv/ctx->M;
447: if (!ctx->L_max) {
448: ctx->L_max = 1;
449: nep->ncv = ctx->L_max*ctx->M;
450: }
451: }
452: ctx->L = PetscMin(ctx->L,ctx->L_max);
453: if (!nep->max_it) nep->max_it = 1;
454: if (!nep->mpd) nep->mpd = nep->ncv;
455: if (!nep->which) nep->which = NEP_ALL;
456: if (nep->stopping!=NEPStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support user-defined stopping test");
458: /* check region */
459: RGIsTrivial(nep->rg,&istrivial);
460: if (istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
461: RGGetComplement(nep->rg,&flg);
462: if (flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
463: PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse);
464: if (!isellipse) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Currently only implemented for elliptic regions");
465: /* if useconj has changed, then reset subcomm data */
466: NEPCISSSetUseConj(nep,&useconj);
467: if (useconj!=ctx->useconj) { NEPCISSResetSubcomm(nep); }
469: /* create split comm */
470: if (!ctx->subcomm) { NEPCISSSetUpSubComm(nep,&ctx->num_solve_point); }
472: NEPAllocateSolution(nep,0);
473: if (ctx->weight) { PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma); }
474: PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma);
475: PetscLogObjectMemory((PetscObject)nep,3*ctx->N*sizeof(PetscScalar)+ctx->L_max*ctx->N*sizeof(PetscReal));
477: /* allocate basis vectors */
478: BVDestroy(&ctx->S);
479: BVDuplicateResize(nep->V,ctx->L_max*ctx->M,&ctx->S);
480: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->S);
481: BVDestroy(&ctx->V);
482: BVDuplicateResize(nep->V,ctx->L_max,&ctx->V);
483: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->V);
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,¢er,&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: for (i=0;i<nv;i++){
579: nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
580: #if !defined(PETSC_USE_COMPLEX)
581: nep->eigi[i] = nep->eigi[i]*radius*rgscale;
582: #endif
583: }
584: PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr);
585: isGhost(nep,ld,nv,fl1);
586: RGCheckInside(nep->rg,nv,nep->eigr,nep->eigi,inside);
587: for (i=0;i<nv;i++) {
588: if (fl1[i] && inside[i]>0) {
589: rr[i] = 1.0;
590: nep->nconv++;
591: } else rr[i] = 0.0;
592: }
593: DSSort(nep->ds,nep->eigr,nep->eigi,rr,NULL,&nep->nconv);
594: DSSynchronize(nep->ds,nep->eigr,nep->eigi);
595: for (i=0;i<nv;i++){
596: nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
597: #if !defined(PETSC_USE_COMPLEX)
598: nep->eigi[i] = nep->eigi[i]*radius*rgscale;
599: #endif
600: }
601: PetscFree3(fl1,inside,rr);
602: BVSetActiveColumns(nep->V,0,nv);
603: ConstructS(nep);
604: BVSetActiveColumns(ctx->S,0,nv);
605: BVCopy(ctx->S,nep->V);
607: DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
608: DSGetMat(nep->ds,DS_MAT_X,&X);
609: BVMultInPlace(ctx->S,X,0,nep->nconv);
610: BVMultInPlace(nep->V,X,0,nep->nconv);
611: MatDestroy(&X);
612: max_error = 0.0;
613: for (i=0;i<nep->nconv;i++) {
614: BVGetColumn(nep->V,i,&si);
615: VecNormalize(si,NULL);
616: NEPComputeResidualNorm_Private(nep,nep->eigr[i],si,nep->work,&error);
617: (*nep->converged)(nep,nep->eigr[i],0,error,&error,nep->convergedctx);
618: BVRestoreColumn(nep->V,i,&si);
619: max_error = PetscMax(max_error,error);
620: }
621: if (max_error <= nep->tol) nep->reason = NEP_CONVERGED_TOL;
622: else if (nep->its > nep->max_it) nep->reason = NEP_DIVERGED_ITS;
623: else {
624: if (nep->nconv > ctx->L) nv = nep->nconv;
625: else if (ctx->L > nv) nv = ctx->L;
626: MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M);
627: MatDenseGetArray(M,&temp);
628: for (i=0;i<ctx->L*nv;i++) {
629: PetscRandomGetValue(rand,&temp[i]);
630: temp[i] = PetscRealPart(temp[i]);
631: }
632: MatDenseRestoreArray(M,&temp);
633: BVSetActiveColumns(ctx->S,0,nv);
634: BVMultInPlace(ctx->S,M,0,ctx->L);
635: MatDestroy(&M);
636: BVSetActiveColumns(ctx->S,0,ctx->L);
637: BVCopy(ctx->S,ctx->V);
638: SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,0,ctx->L,PETSC_FALSE);
639: }
640: }
641: PetscFree3(Mu,H0,H1);
642: return(0);
643: }
645: static PetscErrorCode NEPCISSSetSizes_CISS(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
646: {
648: NEP_CISS *ctx = (NEP_CISS*)nep->data;
649: PetscInt oN,onpart;
652: oN = ctx->N;
653: if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
654: if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
655: } else {
656: if (ip<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
657: if (ip%2) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
658: if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
659: }
660: if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
661: ctx->L = 16;
662: } else {
663: if (bs<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
664: ctx->L = bs;
665: }
666: if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
667: ctx->M = ctx->N/4;
668: } else {
669: if (ms<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
670: 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");
671: ctx->M = ms;
672: }
673: onpart = ctx->npart;
674: if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
675: ctx->npart = 1;
676: } else {
677: if (npart<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The npart argument must be > 0");
678: ctx->npart = npart;
679: if (npart>1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Current implementation does not allow partitions");
680: }
681: if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
682: ctx->L_max = 64;
683: } else {
684: if (bsmax<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
685: ctx->L_max = PetscMax(bsmax,ctx->L);
686: }
687: if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) { NEPCISSResetSubcomm(nep); }
688: ctx->isreal = realmats;
689: nep->state = NEP_STATE_INITIAL;
690: return(0);
691: }
693: /*@
694: NEPCISSSetSizes - Sets the values of various size parameters in the CISS solver.
696: Logically Collective on NEP
698: Input Parameters:
699: + nep - the eigenproblem solver context
700: . ip - number of integration points
701: . bs - block size
702: . ms - moment size
703: . npart - number of partitions when splitting the communicator
704: . bsmax - max block size
705: - realmats - T(z) is real for real z
707: Options Database Keys:
708: + -nep_ciss_integration_points - Sets the number of integration points
709: . -nep_ciss_blocksize - Sets the block size
710: . -nep_ciss_moments - Sets the moment size
711: . -nep_ciss_partitions - Sets the number of partitions
712: . -nep_ciss_maxblocksize - Sets the maximum block size
713: - -nep_ciss_realmats - T(z) is real for real z
715: Note:
716: The default number of partitions is 1. This means the internal KSP object is shared
717: among all processes of the NEP communicator. Otherwise, the communicator is split
718: into npart communicators, so that npart KSP solves proceed simultaneously.
720: The realmats flag can be set to true when T(.) is guaranteed to be real
721: when the argument is a real value, for example, when all matrices in
722: the split form are real. When set to true, the solver avoids some computations.
724: Level: advanced
726: .seealso: NEPCISSGetSizes()
727: @*/
728: PetscErrorCode NEPCISSSetSizes(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
729: {
740: PetscTryMethod(nep,"NEPCISSSetSizes_C",(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(nep,ip,bs,ms,npart,bsmax,realmats));
741: return(0);
742: }
744: static PetscErrorCode NEPCISSGetSizes_CISS(NEP nep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
745: {
746: NEP_CISS *ctx = (NEP_CISS*)nep->data;
749: if (ip) *ip = ctx->N;
750: if (bs) *bs = ctx->L;
751: if (ms) *ms = ctx->M;
752: if (npart) *npart = ctx->npart;
753: if (bsmax) *bsmax = ctx->L_max;
754: if (realmats) *realmats = ctx->isreal;
755: return(0);
756: }
758: /*@
759: NEPCISSGetSizes - Gets the values of various size parameters in the CISS solver.
761: Not Collective
763: Input Parameter:
764: . nep - the eigenproblem solver context
766: Output Parameters:
767: + ip - number of integration points
768: . bs - block size
769: . ms - moment size
770: . npart - number of partitions when splitting the communicator
771: . bsmax - max block size
772: - realmats - T(z) is real for real z
774: Level: advanced
776: .seealso: NEPCISSSetSizes()
777: @*/
778: PetscErrorCode NEPCISSGetSizes(NEP nep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
779: {
784: PetscUseMethod(nep,"NEPCISSGetSizes_C",(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(nep,ip,bs,ms,npart,bsmax,realmats));
785: return(0);
786: }
788: static PetscErrorCode NEPCISSSetThreshold_CISS(NEP nep,PetscReal delta,PetscReal spur)
789: {
790: NEP_CISS *ctx = (NEP_CISS*)nep->data;
793: if (delta == PETSC_DEFAULT) {
794: ctx->delta = 1e-12;
795: } else {
796: if (delta<=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
797: ctx->delta = delta;
798: }
799: if (spur == PETSC_DEFAULT) {
800: ctx->spurious_threshold = 1e-4;
801: } else {
802: if (spur<=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
803: ctx->spurious_threshold = spur;
804: }
805: return(0);
806: }
808: /*@
809: NEPCISSSetThreshold - Sets the values of various threshold parameters in
810: the CISS solver.
812: Logically Collective on NEP
814: Input Parameters:
815: + nep - the eigenproblem solver context
816: . delta - threshold for numerical rank
817: - spur - spurious threshold (to discard spurious eigenpairs)
819: Options Database Keys:
820: + -nep_ciss_delta - Sets the delta
821: - -nep_ciss_spurious_threshold - Sets the spurious threshold
823: Level: advanced
825: .seealso: NEPCISSGetThreshold()
826: @*/
827: PetscErrorCode NEPCISSSetThreshold(NEP nep,PetscReal delta,PetscReal spur)
828: {
835: PetscTryMethod(nep,"NEPCISSSetThreshold_C",(NEP,PetscReal,PetscReal),(nep,delta,spur));
836: return(0);
837: }
839: static PetscErrorCode NEPCISSGetThreshold_CISS(NEP nep,PetscReal *delta,PetscReal *spur)
840: {
841: NEP_CISS *ctx = (NEP_CISS*)nep->data;
844: if (delta) *delta = ctx->delta;
845: if (spur) *spur = ctx->spurious_threshold;
846: return(0);
847: }
849: /*@
850: NEPCISSGetThreshold - Gets the values of various threshold parameters in
851: the CISS solver.
853: Not Collective
855: Input Parameter:
856: . nep - the eigenproblem solver context
858: Output Parameters:
859: + delta - threshold for numerical rank
860: - spur - spurious threshold (to discard spurious eigenpairs)
862: Level: advanced
864: .seealso: NEPCISSSetThreshold()
865: @*/
866: PetscErrorCode NEPCISSGetThreshold(NEP nep,PetscReal *delta,PetscReal *spur)
867: {
872: PetscUseMethod(nep,"NEPCISSGetThreshold_C",(NEP,PetscReal*,PetscReal*),(nep,delta,spur));
873: return(0);
874: }
876: static PetscErrorCode NEPCISSSetRefinement_CISS(NEP nep,PetscInt inner,PetscInt blsize)
877: {
878: NEP_CISS *ctx = (NEP_CISS*)nep->data;
881: if (inner == PETSC_DEFAULT) {
882: ctx->refine_inner = 0;
883: } else {
884: if (inner<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
885: ctx->refine_inner = inner;
886: }
887: if (blsize == PETSC_DEFAULT) {
888: ctx->refine_blocksize = 0;
889: } else {
890: if (blsize<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
891: ctx->refine_blocksize = blsize;
892: }
893: return(0);
894: }
896: /*@
897: NEPCISSSetRefinement - Sets the values of various refinement parameters
898: in the CISS solver.
900: Logically Collective on NEP
902: Input Parameters:
903: + nep - the eigenproblem solver context
904: . inner - number of iterative refinement iterations (inner loop)
905: - blsize - number of iterative refinement iterations (blocksize loop)
907: Options Database Keys:
908: + -nep_ciss_refine_inner - Sets number of inner iterations
909: - -nep_ciss_refine_blocksize - Sets number of blocksize iterations
911: Level: advanced
913: .seealso: NEPCISSGetRefinement()
914: @*/
915: PetscErrorCode NEPCISSSetRefinement(NEP nep,PetscInt inner,PetscInt blsize)
916: {
923: PetscTryMethod(nep,"NEPCISSSetRefinement_C",(NEP,PetscInt,PetscInt),(nep,inner,blsize));
924: return(0);
925: }
927: static PetscErrorCode NEPCISSGetRefinement_CISS(NEP nep,PetscInt *inner,PetscInt *blsize)
928: {
929: NEP_CISS *ctx = (NEP_CISS*)nep->data;
932: if (inner) *inner = ctx->refine_inner;
933: if (blsize) *blsize = ctx->refine_blocksize;
934: return(0);
935: }
937: /*@
938: NEPCISSGetRefinement - Gets the values of various refinement parameters
939: in the CISS solver.
941: Not Collective
943: Input Parameter:
944: . nep - the eigenproblem solver context
946: Output Parameters:
947: + inner - number of iterative refinement iterations (inner loop)
948: - blsize - number of iterative refinement iterations (blocksize loop)
950: Level: advanced
952: .seealso: NEPCISSSetRefinement()
953: @*/
954: PetscErrorCode NEPCISSGetRefinement(NEP nep, PetscInt *inner, PetscInt *blsize)
955: {
960: PetscUseMethod(nep,"NEPCISSGetRefinement_C",(NEP,PetscInt*,PetscInt*),(nep,inner,blsize));
961: return(0);
962: }
964: static PetscErrorCode NEPCISSGetKSPs_CISS(NEP nep,PetscInt *nsolve,KSP **ksp)
965: {
967: NEP_CISS *ctx = (NEP_CISS*)nep->data;
968: PetscInt i;
969: PC pc;
972: if (!ctx->ksp) {
973: if (!ctx->subcomm) { /* initialize subcomm first */
974: NEPCISSSetUseConj(nep,&ctx->useconj);
975: NEPCISSSetUpSubComm(nep,&ctx->num_solve_point);
976: }
977: PetscMalloc1(ctx->num_solve_point,&ctx->ksp);
978: for (i=0;i<ctx->num_solve_point;i++) {
979: KSPCreate(PetscSubcommChild(ctx->subcomm),&ctx->ksp[i]);
980: PetscObjectIncrementTabLevel((PetscObject)ctx->ksp[i],(PetscObject)nep,1);
981: KSPSetOptionsPrefix(ctx->ksp[i],((PetscObject)nep)->prefix);
982: KSPAppendOptionsPrefix(ctx->ksp[i],"nep_ciss_");
983: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp[i]);
984: KSPSetErrorIfNotConverged(ctx->ksp[i],PETSC_TRUE);
985: KSPSetTolerances(ctx->ksp[i],SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
986: KSPGetPC(ctx->ksp[i],&pc);
987: KSPSetType(ctx->ksp[i],KSPPREONLY);
988: PCSetType(pc,PCLU);
989: }
990: }
991: if (nsolve) *nsolve = ctx->num_solve_point;
992: if (ksp) *ksp = ctx->ksp;
993: return(0);
994: }
996: /*@C
997: NEPCISSGetKSPs - Retrieve the array of linear solver objects associated with
998: the CISS solver.
1000: Not Collective
1002: Input Parameter:
1003: . nep - nonlinear eigenvalue solver
1005: Output Parameters:
1006: + nsolve - number of solver objects
1007: - ksp - array of linear solver object
1009: Notes:
1010: The number of KSP solvers is equal to the number of integration points divided by
1011: the number of partitions. This value is halved in the case of real matrices with
1012: a region centered at the real axis.
1014: Level: advanced
1016: .seealso: NEPCISSSetSizes()
1017: @*/
1018: PetscErrorCode NEPCISSGetKSPs(NEP nep,PetscInt *nsolve,KSP **ksp)
1019: {
1024: PetscUseMethod(nep,"NEPCISSGetKSPs_C",(NEP,PetscInt*,KSP**),(nep,nsolve,ksp));
1025: return(0);
1026: }
1028: PetscErrorCode NEPReset_CISS(NEP nep)
1029: {
1031: PetscInt i;
1032: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1035: BVDestroy(&ctx->S);
1036: BVDestroy(&ctx->V);
1037: BVDestroy(&ctx->Y);
1038: for (i=0;i<ctx->num_solve_point;i++) {
1039: KSPReset(ctx->ksp[i]);
1040: }
1041: return(0);
1042: }
1044: PetscErrorCode NEPSetFromOptions_CISS(PetscOptionItems *PetscOptionsObject,NEP nep)
1045: {
1047: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1048: PetscReal r1,r2;
1049: PetscInt i,i1,i2,i3,i4,i5,i6,i7;
1050: PetscBool b1;
1053: PetscOptionsHead(PetscOptionsObject,"NEP CISS Options");
1055: NEPCISSGetSizes(nep,&i1,&i2,&i3,&i4,&i5,&b1);
1056: PetscOptionsInt("-nep_ciss_integration_points","Number of integration points","NEPCISSSetSizes",i1,&i1,NULL);
1057: PetscOptionsInt("-nep_ciss_blocksize","Block size","NEPCISSSetSizes",i2,&i2,NULL);
1058: PetscOptionsInt("-nep_ciss_moments","Moment size","NEPCISSSetSizes",i3,&i3,NULL);
1059: PetscOptionsInt("-nep_ciss_partitions","Number of partitions","NEPCISSSetSizes",i4,&i4,NULL);
1060: PetscOptionsInt("-nep_ciss_maxblocksize","Maximum block size","NEPCISSSetSizes",i5,&i5,NULL);
1061: PetscOptionsBool("-nep_ciss_realmats","True if T(z) is real for real z","NEPCISSSetSizes",b1,&b1,NULL);
1062: NEPCISSSetSizes(nep,i1,i2,i3,i4,i5,b1);
1064: NEPCISSGetThreshold(nep,&r1,&r2);
1065: PetscOptionsReal("-nep_ciss_delta","Threshold for numerical rank","NEPCISSSetThreshold",r1,&r1,NULL);
1066: PetscOptionsReal("-nep_ciss_spurious_threshold","Threshold for the spurious eigenpairs","NEPCISSSetThreshold",r2,&r2,NULL);
1067: NEPCISSSetThreshold(nep,r1,r2);
1069: NEPCISSGetRefinement(nep,&i6,&i7);
1070: PetscOptionsInt("-nep_ciss_refine_inner","Number of inner iterative refinement iterations","NEPCISSSetRefinement",i6,&i6,NULL);
1071: PetscOptionsInt("-nep_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","NEPCISSSetRefinement",i7,&i7,NULL);
1072: NEPCISSSetRefinement(nep,i6,i7);
1074: PetscOptionsTail();
1076: if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
1077: RGSetFromOptions(nep->rg); /* this is necessary here to set useconj */
1078: if (!ctx->ksp) { NEPCISSGetKSPs(nep,&ctx->num_solve_point,&ctx->ksp); }
1079: for (i=0;i<ctx->num_solve_point;i++) {
1080: KSPSetFromOptions(ctx->ksp[i]);
1081: }
1082: PetscSubcommSetFromOptions(ctx->subcomm);
1083: return(0);
1084: }
1086: PetscErrorCode NEPDestroy_CISS(NEP nep)
1087: {
1089: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1092: NEPCISSResetSubcomm(nep);
1093: PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
1094: PetscFree(nep->data);
1095: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NULL);
1096: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NULL);
1097: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NULL);
1098: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NULL);
1099: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NULL);
1100: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NULL);
1101: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetKSPs_C",NULL);
1102: return(0);
1103: }
1105: PetscErrorCode NEPView_CISS(NEP nep,PetscViewer viewer)
1106: {
1108: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1109: PetscBool isascii;
1112: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1113: if (isascii) {
1114: 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);
1115: if (ctx->isreal) {
1116: PetscViewerASCIIPrintf(viewer," exploiting symmetry of integration points\n");
1117: }
1118: PetscViewerASCIIPrintf(viewer," threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold);
1119: PetscViewerASCIIPrintf(viewer," iterative refinement { inner: %D, blocksize: %D }\n",ctx->refine_inner, ctx->refine_blocksize);
1120: if (!ctx->ksp) { NEPCISSGetKSPs(nep,&ctx->num_solve_point,&ctx->ksp); }
1121: KSPView(ctx->ksp[0],viewer);
1122: }
1123: return(0);
1124: }
1126: PETSC_EXTERN PetscErrorCode NEPCreate_CISS(NEP nep)
1127: {
1129: NEP_CISS *ctx = (NEP_CISS*)nep->data;
1132: PetscNewLog(nep,&ctx);
1133: nep->data = ctx;
1134: /* set default values of parameters */
1135: ctx->N = 32;
1136: ctx->L = 16;
1137: ctx->M = ctx->N/4;
1138: ctx->delta = 1e-12;
1139: ctx->L_max = 64;
1140: ctx->spurious_threshold = 1e-4;
1141: ctx->isreal = PETSC_FALSE;
1142: ctx->npart = 1;
1144: nep->ops->solve = NEPSolve_CISS;
1145: nep->ops->setup = NEPSetUp_CISS;
1146: nep->ops->setfromoptions = NEPSetFromOptions_CISS;
1147: nep->ops->reset = NEPReset_CISS;
1148: nep->ops->destroy = NEPDestroy_CISS;
1149: nep->ops->view = NEPView_CISS;
1151: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NEPCISSSetSizes_CISS);
1152: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NEPCISSGetSizes_CISS);
1153: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NEPCISSSetThreshold_CISS);
1154: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NEPCISSGetThreshold_CISS);
1155: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NEPCISSSetRefinement_CISS);
1156: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NEPCISSGetRefinement_CISS);
1157: PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetKSPs_C",NEPCISSGetKSPs_CISS);
1158: return(0);
1159: }