Actual source code: nciss.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*

  3:    SLEPc eigensolver: "ciss"

  5:    Method: Contour Integral Spectral Slicing

  7:    Algorithm:

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

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

 15:    References:

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

 20:        [2] T. Sakurai and H. Tadano, "CIRR: a Rayleigh-Ritz type method with
 21:            contour integral for generalized eigenvalue problems", Hokkaido
 22:            Math. J. 36:745-757, 2007.

 24:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 25:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 26:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

 28:    This file is part of SLEPc.

 30:    SLEPc is free software: you can redistribute it and/or modify it under  the
 31:    terms of version 3 of the GNU Lesser General Public License as published by
 32:    the Free Software Foundation.

 34:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 35:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 36:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 37:    more details.

 39:    You  should have received a copy of the GNU Lesser General  Public  License
 40:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 41:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 42: */

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

 47: typedef struct {
 48:   /* parameters */
 49:   PetscInt     N;          /* number of integration points (32) */
 50:   PetscInt     L;          /* block size (16) */
 51:   PetscInt     M;          /* moment degree (N/4 = 4) */
 52:   PetscReal    delta;      /* threshold of singular value (1e-12) */
 53:   PetscInt     L_max;      /* maximum number of columns of the source matrix V */
 54:   PetscReal    spurious_threshold; /* discard spurious eigenpairs */
 55:   PetscBool    isreal;     /* T(z) is real for real z */
 56:   PetscInt     refine_inner;
 57:   PetscInt     refine_blocksize;
 58:   /* private data */
 59:   PetscReal    *sigma;     /* threshold for numerical rank */
 60:   PetscInt     num_subcomm;
 61:   PetscInt     subcomm_id;
 62:   PetscInt     num_solve_point;
 63:   PetscScalar  *weight;
 64:   PetscScalar  *omega;
 65:   PetscScalar  *pp;
 66:   BV           V;
 67:   BV           S;
 68:   BV           Y;
 69:   KSP          *ksp;
 70:   Mat          *kspMat;
 71:   PetscBool    useconj;
 72:   PetscReal    est_eig;
 73:   PetscSubcomm subcomm;
 74:   PetscBool    usest;
 75: } NEP_CISS;

 79: static PetscErrorCode SetSolverComm(NEP nep)
 80: {
 82:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
 83:   PetscInt       N = ctx->N;

 86:   if (ctx->useconj) N = N/2;
 87:   if (!ctx->subcomm) {
 88:     PetscSubcommCreate(PetscObjectComm((PetscObject)nep),&ctx->subcomm);
 89:     PetscSubcommSetNumber(ctx->subcomm,ctx->num_subcomm);
 90:     PetscSubcommSetType(ctx->subcomm,PETSC_SUBCOMM_INTERLACED);
 91:     PetscLogObjectMemory((PetscObject)nep,sizeof(PetscSubcomm));
 92:     PetscSubcommSetFromOptions(ctx->subcomm);
 93:   }
 94:   ctx->subcomm_id = ctx->subcomm->color;
 95:   ctx->num_solve_point = N / ctx->num_subcomm;
 96:   if ((N%ctx->num_subcomm) > ctx->subcomm_id) ctx->num_solve_point+=1;
 97:   return(0);
 98: }

102: static PetscErrorCode SetPathParameter(NEP nep)
103: {
105:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
106:   PetscInt       i;
107:   PetscScalar    center;
108:   PetscReal      theta,radius,vscale,rgscale;
109:   PetscBool      isellipse=PETSC_FALSE;

112:   PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse);
113:   RGGetScale(nep->rg,&rgscale);
114:   if (isellipse) {
115:     RGEllipseGetParameters(nep->rg,&center,&radius,&vscale);
116:   } else SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Region must be Ellipse");
117:   for (i=0;i<ctx->N;i++) {
118:     theta = ((2*PETSC_PI)/ctx->N)*(i+0.5);
119:     ctx->pp[i] = PetscCosReal(theta) + PETSC_i*vscale*PetscSinReal(theta);
120:     ctx->weight[i] = radius*(vscale*PetscCosReal(theta) + PETSC_i*PetscSinReal(theta))/(PetscReal)ctx->N;
121:     ctx->omega[i] = rgscale*(center + radius*ctx->pp[i]);
122:   }
123:   return(0);
124: }

128: static PetscErrorCode CISSVecSetRandom(BV V,PetscInt i0,PetscInt i1)
129: {
131:   PetscInt       i,j,nlocal;
132:   PetscScalar    *vdata;
133:   Vec            x;

136:   BVGetSizes(V,&nlocal,NULL,NULL);
137:   for (i=i0;i<i1;i++) {
138:     BVSetRandomColumn(V,i);
139:     BVGetColumn(V,i,&x);
140:     VecGetArray(x,&vdata);
141:     for (j=0;j<nlocal;j++) {
142:       vdata[j] = PetscRealPart(vdata[j]);
143:       if (PetscRealPart(vdata[j]) < 0.5) vdata[j] = -1.0;
144:       else vdata[j] = 1.0;
145:     }
146:     VecRestoreArray(x,&vdata);
147:     BVRestoreColumn(V,i,&x);
148:   }
149:   return(0);
150: }

154: static PetscErrorCode SolveLinearSystem(NEP nep,Mat T,Mat dT,BV V,PetscInt L_start,PetscInt L_end,PetscBool initksp)
155: {
157:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
158:   PetscInt       i,j,p_id;
159:   Mat            Fz;
160:   PC             pc;
161:   Vec            Bvj,vj,yj;
162:   KSP            ksp;

165:   if (ctx->usest) {
166:     NEPComputeFunction(nep,0,T,T);
167:     MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&Fz);
168:     KSPCreate(PetscObjectComm((PetscObject)nep),&ksp);
169:   }
170:   BVCreateVec(V,&Bvj);
171:   for (i=0;i<ctx->num_solve_point;i++) {
172:     p_id = i*ctx->subcomm->n + ctx->subcomm_id;
173:     NEPComputeFunction(nep,ctx->omega[p_id],T,T);
174:     NEPComputeJacobian(nep,ctx->omega[p_id],dT);
175:     if (!ctx->usest && initksp == PETSC_TRUE) {
176:       MatDuplicate(T,MAT_COPY_VALUES,&ctx->kspMat[i]);
177:       KSPSetOperators(ctx->ksp[i],ctx->kspMat[i],ctx->kspMat[i]);
178:       KSPSetType(ctx->ksp[i],KSPPREONLY);
179:       KSPGetPC(ctx->ksp[i],&pc);
180:       PCSetType(pc,PCLU);
181:       KSPSetFromOptions(ctx->ksp[i]);
182:     } else if (ctx->usest) {
183:       MatCopy(T,Fz,DIFFERENT_NONZERO_PATTERN);
184:       KSPSetOperators(ksp,Fz,Fz);
185:       KSPSetType(ksp,KSPPREONLY);
186:       KSPGetPC(ksp,&pc);
187:       PCSetType(pc,PCLU);
188:       KSPSetFromOptions(ksp);
189:     }
190:     for (j=L_start;j<L_end;j++) {
191:       BVGetColumn(V,j,&vj);
192:       BVGetColumn(ctx->Y,i*ctx->L_max+j,&yj);
193:       MatMult(dT,vj,Bvj);
194:       if (ctx->usest) {
195:         KSPSolve(ksp,Bvj,yj);
196:       } else {
197:         KSPSolve(ctx->ksp[i],Bvj,yj);
198:       }
199:       BVRestoreColumn(V,j,&vj);
200:       BVRestoreColumn(ctx->Y,i*ctx->L_max+j,&yj);
201:     }
202:     if (ctx->usest && i<ctx->num_solve_point-1) {  KSPReset(ksp); }
203:   }
204:   if (ctx->usest) {
205:     MatDestroy(&Fz);
206:     KSPDestroy(&ksp);
207:   }
208:   VecDestroy(&Bvj);
209:   return(0);
210: }

214: static PetscErrorCode EstimateNumberEigs(NEP nep,PetscInt *L_add)
215: {
217:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
218:   PetscInt       i,j,p_id;
219:   PetscScalar    tmp,m = 1,sum = 0.0;
220:   PetscReal      eta;
221:   Vec            v,vtemp,vj,yj;

224:   BVGetColumn(ctx->Y,0,&yj);
225:   VecDuplicate(yj,&v);
226:   BVRestoreColumn(ctx->Y,0,&yj);
227:   BVCreateVec(ctx->V,&vtemp);
228:   for (j=0;j<ctx->L;j++) {
229:     VecSet(v,0);
230:     for (i=0;i<ctx->num_solve_point; i++) {
231:       p_id = i*ctx->subcomm->n + ctx->subcomm_id;
232:       BVSetActiveColumns(ctx->Y,i*ctx->L_max+j,i*ctx->L_max+j+1);
233:       BVMultVec(ctx->Y,ctx->weight[p_id],1,v,&m);
234:     }
235:     BVGetColumn(ctx->V,j,&vj);
236:     VecDot(vj,v,&tmp);
237:     BVRestoreColumn(ctx->V,j,&vj);
238:     if (ctx->useconj) sum += PetscRealPart(tmp)*2;
239:     else sum += tmp;
240:   }
241:   ctx->est_eig = PetscAbsScalar(sum/(PetscReal)ctx->L);
242:   eta = PetscPowReal(10,-PetscLog10Real(nep->tol)/ctx->N);
243:   PetscInfo1(nep,"Estimation_#Eig %f\n",(double)ctx->est_eig);
244:   *L_add = (PetscInt)PetscCeilReal((ctx->est_eig*eta)/ctx->M) - ctx->L;
245:   if (*L_add < 0) *L_add = 0;
246:   if (*L_add>ctx->L_max-ctx->L) {
247:     PetscInfo(nep,"Number of eigenvalues around the contour path may be too large\n");
248:     *L_add = ctx->L_max-ctx->L;
249:   }
250:   VecDestroy(&v);
251:   VecDestroy(&vtemp);
252:   return(0);
253: }

257: static PetscErrorCode CalcMu(NEP nep, PetscScalar *Mu)
258: {
260:   PetscMPIInt    sub_size,len;
261:   PetscInt       i,j,k,s;
262:   PetscScalar    *m,*temp,*temp2,*ppk,alp;
263:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
264:   Mat            M;

267:   MPI_Comm_size(PetscSubcommChild(ctx->subcomm),&sub_size);
268:   PetscMalloc3(ctx->num_solve_point*ctx->L*(ctx->L+1),&temp,2*ctx->M*ctx->L*ctx->L,&temp2,ctx->num_solve_point,&ppk);
269:   MatCreateSeqDense(PETSC_COMM_SELF,ctx->L,ctx->L_max*ctx->num_solve_point,NULL,&M);
270:   for (i=0;i<2*ctx->M*ctx->L*ctx->L;i++) temp2[i] = 0;
271:   BVSetActiveColumns(ctx->Y,0,ctx->L_max*ctx->num_solve_point);
272:   BVSetActiveColumns(ctx->V,0,ctx->L);
273:   BVDot(ctx->Y,ctx->V,M);
274:   MatDenseGetArray(M,&m);
275:   for (i=0;i<ctx->num_solve_point;i++) {
276:     for (j=0;j<ctx->L;j++) {
277:       for (k=0;k<ctx->L;k++) {
278:         temp[k+j*ctx->L+i*ctx->L*ctx->L]=m[k+j*ctx->L+i*ctx->L*ctx->L_max];
279:       }
280:     }
281:   }
282:   MatDenseRestoreArray(M,&m);
283:   for (i=0;i<ctx->num_solve_point;i++) ppk[i] = 1;
284:   for (k=0;k<2*ctx->M;k++) {
285:     for (j=0;j<ctx->L;j++) {
286:       for (i=0;i<ctx->num_solve_point;i++) {
287:         alp = ppk[i]*ctx->weight[i*ctx->subcomm->n + ctx->subcomm_id];
288:         for (s=0;s<ctx->L;s++) {
289:           if (ctx->useconj) temp2[s+(j+k*ctx->L)*ctx->L] += PetscRealPart(alp*temp[s+(j+i*ctx->L)*ctx->L])*2;
290:           else temp2[s+(j+k*ctx->L)*ctx->L] += alp*temp[s+(j+i*ctx->L)*ctx->L];
291:         }
292:       }
293:     }
294:     for (i=0;i<ctx->num_solve_point;i++) 
295:       ppk[i] *= ctx->pp[i*ctx->subcomm->n + ctx->subcomm_id];
296:   }
297:   for (i=0;i<2*ctx->M*ctx->L*ctx->L;i++) temp2[i] /= sub_size;
298:   PetscMPIIntCast(2*ctx->M*ctx->L*ctx->L,&len);
299:   MPI_Allreduce(temp2,Mu,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)nep));
300:   PetscFree3(temp,temp2,ppk);
301:   MatDestroy(&M);
302:   return(0);
303: }

307: static PetscErrorCode BlockHankel(NEP nep,PetscScalar *Mu,PetscInt s,PetscScalar *H)
308: {
309:   NEP_CISS *ctx = (NEP_CISS*)nep->data;
310:   PetscInt  i,j,k,L=ctx->L,M=ctx->M;

313:   for (k=0;k<L*M;k++)
314:     for (j=0;j<M;j++) 
315:       for (i=0;i<L;i++) 
316:         H[j*L+i+k*L*M] = Mu[i+k*L+(j+s)*L*L];
317:   return(0);
318: }

322: static PetscErrorCode SVD_H0(NEP nep,PetscScalar *S,PetscInt *K)
323: {
324: #if defined(PETSC_MISSING_LAPACK_GESVD)
326:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable");
327: #else
329:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
330:   PetscInt       i,ml=ctx->L*ctx->M;
331:   PetscBLASInt   m,n,lda,ldu,ldvt,lwork,info;
332:   PetscScalar    *work;
333: #if defined(PETSC_USE_COMPLEX)
334:   PetscReal      *rwork;
335: #endif

338:   PetscMalloc1(5*ml,&work);
339: #if defined(PETSC_USE_COMPLEX)
340:   PetscMalloc1(5*ml,&rwork);
341: #endif
342:   PetscBLASIntCast(ml,&m);
343:   n = m; lda = m; ldu = m; ldvt = m; lwork = 5*m;
344:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
345: #if defined(PETSC_USE_COMPLEX)
346:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,S,&lda,ctx->sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,rwork,&info));
347: #else
348:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,S,&lda,ctx->sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,&info));
349: #endif
350:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESVD %d",info);
351:   PetscFPTrapPop();
352:   (*K) = 0;
353:   for (i=0;i<ml;i++) {
354:     if (ctx->sigma[i]/PetscMax(ctx->sigma[0],1)>ctx->delta) (*K)++;
355:   }
356:   PetscFree(work);
357: #if defined(PETSC_USE_COMPLEX)
358:   PetscFree(rwork);
359: #endif
360:   return(0);
361: #endif
362: }

366: static PetscErrorCode ConstructS(NEP nep)
367: {
369:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
370:   PetscInt       i,j,k,vec_local_size,p_id;
371:   Vec            v,sj,yj;
372:   PetscScalar    *ppk, *v_data, m = 1;

375:   BVGetSizes(ctx->Y,&vec_local_size,NULL,NULL);
376:   PetscMalloc1(ctx->num_solve_point,&ppk);
377:   for (i=0;i<ctx->num_solve_point;i++) ppk[i] = 1;
378:   BVGetColumn(ctx->Y,0,&yj);
379:   VecDuplicate(yj,&v);
380:   BVRestoreColumn(ctx->Y,0,&yj);
381:   for (k=0;k<ctx->M;k++) {
382:     for (j=0;j<ctx->L;j++) {
383:       VecSet(v,0);
384:       for (i=0;i<ctx->num_solve_point;i++) {
385:         p_id = i*ctx->subcomm->n + ctx->subcomm_id;
386:         BVSetActiveColumns(ctx->Y,i*ctx->L_max+j,i*ctx->L_max+j+1);
387:         BVMultVec(ctx->Y,ppk[i]*ctx->weight[p_id],1,v,&m);
388:       }
389:       if (ctx->useconj) {
390:         VecGetArray(v,&v_data);
391:         for (i=0;i<vec_local_size;i++) v_data[i] = PetscRealPart(v_data[i])*2;
392:         VecRestoreArray(v,&v_data);
393:       }
394:       BVGetColumn(ctx->S,k*ctx->L+j,&sj);
395:       VecCopy(v,sj);
396:       BVRestoreColumn(ctx->S,k*ctx->L+j,&sj);
397:     }
398:     for (i=0;i<ctx->num_solve_point;i++) {
399:       p_id = i*ctx->subcomm->n + ctx->subcomm_id;
400:       ppk[i] *= ctx->pp[p_id];
401:     }
402:   }
403:   PetscFree(ppk);
404:   VecDestroy(&v);
405:   return(0);
406: }

410: static PetscErrorCode isGhost(NEP nep,PetscInt ld,PetscInt nv,PetscBool *fl)
411: {
413:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
414:   PetscInt       i,j;
415:   PetscScalar    *pX;
416:   PetscReal      *tau,s1,s2,tau_max=0.0;

419:   PetscMalloc1(nv,&tau);
420:   DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
421:   DSGetArray(nep->ds,DS_MAT_X,&pX);

423:   for (i=0;i<nv;i++) {
424:     s1 = 0;
425:     s2 = 0;
426:     for (j=0;j<nv;j++) {
427:       s1 += PetscAbsScalar(PetscPowScalar(pX[i*ld+j],2));
428:       s2 += PetscPowReal(PetscAbsScalar(pX[i*ld+j]),2)/ctx->sigma[j];
429:     }
430:     tau[i] = s1/s2;
431:     tau_max = PetscMax(tau_max,tau[i]);
432:   }
433:   DSRestoreArray(nep->ds,DS_MAT_X,&pX);
434:   for (i=0;i<nv;i++) {
435:     tau[i] /= tau_max;
436:   }
437:   for (i=0;i<nv;i++) {
438:     if (tau[i]>=ctx->spurious_threshold) fl[i] = PETSC_TRUE;
439:     else fl[i] = PETSC_FALSE;
440:   }
441:   PetscFree(tau);
442:   return(0);
443: }

447: PetscErrorCode NEPSetUp_CISS(NEP nep)
448: {
450:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
451:   PetscInt       i,nwork;
452:   PetscBool      istrivial,isellipse,flg;
453:   PetscScalar    center;

456:   if (!nep->ncv) nep->ncv = ctx->L_max*ctx->M;
457:   else {
458:     ctx->L_max = nep->ncv/ctx->M;
459:     if (ctx->L_max == 0) {
460:       ctx->L_max = 1;
461:       nep->ncv = ctx->L_max*ctx->M;
462:     }
463:     if (ctx->L > ctx->L_max) ctx->L = ctx->L_max;
464:   }
465:   if (!nep->max_it) nep->max_it = 1;
466:   if (!nep->mpd) nep->mpd = nep->ncv;
467:   if (!nep->which) nep->which = NEP_ALL;
468:   if (nep->stopping!=NEPStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support user-defined stopping test");

470:   /* check region */
471:   RGIsTrivial(nep->rg,&istrivial);
472:   if (istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
473:   RGGetComplement(nep->rg,&flg);
474:   if (flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
475:   PetscObjectTypeCompare((PetscObject)nep->rg,RGELLIPSE,&isellipse);
476:   if (!isellipse) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Currently only implemented for elliptic or arc regions");
477:   RGEllipseGetParameters(nep->rg,&center,NULL,NULL);
478:   if (ctx->isreal && PetscImaginaryPart(center) == 0.0) ctx->useconj = PETSC_TRUE;
479:   else ctx->useconj = PETSC_FALSE;

481:   /* create split comm */
482:   ctx->num_subcomm = 1;
483:   SetSolverComm(nep);

485:   NEPAllocateSolution(nep,0);
486:   PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma);
487:   PetscLogObjectMemory((PetscObject)nep,3*ctx->N*sizeof(PetscScalar)+ctx->L_max*ctx->N*sizeof(PetscReal));

489:   /* allocate basis vectors */
490:   BVDuplicateResize(nep->V,ctx->L_max*ctx->M,&ctx->S);
491:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->S);
492:   BVDuplicateResize(nep->V,ctx->L_max,&ctx->V);
493:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->V);

495:   if (!ctx->usest) {
496:     PetscMalloc2(ctx->num_solve_point,&ctx->ksp,ctx->num_solve_point,&ctx->kspMat);
497:     PetscLogObjectMemory((PetscObject)nep,ctx->num_solve_point*sizeof(KSP)+ctx->num_solve_point*sizeof(Mat));
498:     for (i=0;i<ctx->num_solve_point;i++) {
499:       KSPCreate(PetscSubcommChild(ctx->subcomm),&ctx->ksp[i]);
500:       PetscObjectIncrementTabLevel((PetscObject)ctx->ksp[i],(PetscObject)nep,1);
501:       PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp[i]);
502:       KSPSetOptionsPrefix(ctx->ksp[i],((PetscObject)nep)->prefix);
503:       KSPAppendOptionsPrefix(ctx->ksp[i],"nep_ciss_");
504:       KSPSetErrorIfNotConverged(ctx->ksp[i],PETSC_TRUE);
505:     }
506:   }

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

510:   DSSetType(nep->ds,DSGNHEP);
511:   DSAllocate(nep->ds,nep->ncv);
512:   nwork = (nep->fui==NEP_USER_INTERFACE_SPLIT)? 2: 1;
513:   NEPSetWorkVecs(nep,nwork);
514:   return(0);
515: }

519: PetscErrorCode NEPSolve_CISS(NEP nep)
520: {
522:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
523:   Mat            X,M;
524:   PetscInt       i,j,ld,L_add=0,nv=0,L_base=ctx->L,inner,*inside;
525:   PetscScalar    *Mu,*H0,*H1,*rr,*temp,center;
526:   PetscReal      error,max_error,radius,rgscale;
527:   PetscBool      *fl1;
528:   Vec            si;
529:   SlepcSC        sc;
530:   PetscRandom    rand;

533:   DSGetSlepcSC(nep->ds,&sc);
534:   sc->comparison    = SlepcCompareLargestMagnitude;
535:   sc->comparisonctx = NULL;
536:   sc->map           = NULL;
537:   sc->mapobj        = NULL;
538:   DSGetLeadingDimension(nep->ds,&ld);
539:   SetPathParameter(nep);
540:   CISSVecSetRandom(ctx->V,0,ctx->L);
541:   BVGetRandomContext(ctx->V,&rand);

543:   SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,0,ctx->L,PETSC_TRUE);
544:   EstimateNumberEigs(nep,&L_add);
545:   if (L_add>0) {
546:     PetscInfo2(nep,"Changing L %D -> %D by Estimate #Eig\n",ctx->L,ctx->L+L_add);
547:     CISSVecSetRandom(ctx->V,ctx->L,ctx->L+L_add);
548:     SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
549:     ctx->L += L_add;
550:   }
551:   PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
552:   for (i=0;i<ctx->refine_blocksize;i++) {
553:     CalcMu(nep,Mu);
554:     BlockHankel(nep,Mu,0,H0);
555:     SVD_H0(nep,H0,&nv);
556:     if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
557:     L_add = L_base;
558:     if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
559:     PetscInfo2(nep,"Changing L %D -> %D by SVD(H0)\n",ctx->L,ctx->L+L_add);
560:     CISSVecSetRandom(ctx->V,ctx->L,ctx->L+L_add);
561:     SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
562:     ctx->L += L_add;
563:   }
564:   PetscFree2(Mu,H0);

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

569:   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);
570:   while (nep->reason == NEP_CONVERGED_ITERATING) {
571:     nep->its++;
572:     for (inner=0;inner<=ctx->refine_inner;inner++) {
573:       CalcMu(nep,Mu);
574:       BlockHankel(nep,Mu,0,H0);
575:       SVD_H0(nep,H0,&nv);
576:       if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
577:         ConstructS(nep);
578:         BVSetActiveColumns(ctx->S,0,ctx->L);
579:         BVCopy(ctx->S,ctx->V);
580:         SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,0,ctx->L,PETSC_FALSE);
581:       } else break;
582:     }

584:     nep->nconv = 0;
585:     if (nv == 0) break;
586:     BlockHankel(nep,Mu,0,H0);
587:     BlockHankel(nep,Mu,1,H1);
588:     DSSetDimensions(nep->ds,nv,0,0,0);
589:     DSSetState(nep->ds,DS_STATE_RAW);
590:     DSGetArray(nep->ds,DS_MAT_A,&temp);
591:     for (j=0;j<nv;j++)
592:       for (i=0;i<nv;i++)
593:         temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
594:     DSRestoreArray(nep->ds,DS_MAT_A,&temp);
595:     DSGetArray(nep->ds,DS_MAT_B,&temp);
596:     for (j=0;j<nv;j++) 
597:       for (i=0;i<nv;i++)
598:         temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
599:     DSRestoreArray(nep->ds,DS_MAT_B,&temp);
600:     DSSolve(nep->ds,nep->eigr,nep->eigi);
601:     DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
602:     for (i=0;i<nv;i++){
603:       nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
604: #if !defined(PETSC_USE_COMPLEX)
605:       nep->eigi[i] = nep->eigi[i]*radius*rgscale;
606: #endif
607:     }
608:     PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr);
609:     isGhost(nep,ld,nv,fl1);
610:     RGCheckInside(nep->rg,nv,nep->eigr,nep->eigi,inside);
611:     for (i=0;i<nv;i++) {
612:       if (fl1[i] && inside[i]>0) {
613:         rr[i] = 1.0;
614:         nep->nconv++;
615:       } else rr[i] = 0.0;
616:     }
617:     DSSort(nep->ds,nep->eigr,nep->eigi,rr,NULL,&nep->nconv);
618:     for (i=0;i<nv;i++){
619:       nep->eigr[i] = (nep->eigr[i]*radius+center)*rgscale;
620: #if !defined(PETSC_USE_COMPLEX)
621:       nep->eigi[i] = nep->eigi[i]*radius*rgscale;
622: #endif
623:     }
624:     PetscFree3(fl1,inside,rr);
625:     BVSetActiveColumns(nep->V,0,nv);
626:     ConstructS(nep);
627:     BVSetActiveColumns(ctx->S,0,nv);
628:     BVCopy(ctx->S,nep->V);

630:     DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
631:     DSGetMat(nep->ds,DS_MAT_X,&X);
632:     BVMultInPlace(ctx->S,X,0,nep->nconv);
633:     BVMultInPlace(nep->V,X,0,nep->nconv);
634:     MatDestroy(&X);
635:     max_error = 0.0;
636:     for (i=0;i<nep->nconv;i++) {
637:       BVGetColumn(nep->V,i,&si);
638:       VecNormalize(si,NULL);
639:       NEPComputeResidualNorm_Private(nep,nep->eigr[i],si,nep->work,&error);
640:       (*nep->converged)(nep,nep->eigr[i],0,error,&error,nep->convergedctx);
641:       BVRestoreColumn(nep->V,i,&si);
642:       max_error = PetscMax(max_error,error);
643:     }
644:     if (max_error <= nep->tol) nep->reason = NEP_CONVERGED_TOL;
645:     else if (nep->its > nep->max_it) nep->reason = NEP_DIVERGED_ITS;
646:     else {
647:       if (nep->nconv > ctx->L) nv = nep->nconv;
648:       else if (ctx->L > nv) nv = ctx->L;
649:       MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M);
650:       MatDenseGetArray(M,&temp);
651:       for (i=0;i<ctx->L*nv;i++) {
652:         PetscRandomGetValue(rand,&temp[i]);
653:         temp[i] = PetscRealPart(temp[i]);
654:       }
655:       MatDenseRestoreArray(M,&temp);
656:       BVSetActiveColumns(ctx->S,0,nv);
657:       BVMultInPlace(ctx->S,M,0,ctx->L);
658:       MatDestroy(&M);
659:       BVSetActiveColumns(ctx->S,0,ctx->L);
660:       BVCopy(ctx->S,ctx->V);
661:       SolveLinearSystem(nep,nep->function,nep->jacobian,ctx->V,0,ctx->L,PETSC_FALSE);
662:     }
663:   }
664:   PetscFree3(Mu,H0,H1);  
665:   return(0);
666: }

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

675:   if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
676:     if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
677:   } else {
678:     if (ip<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
679:     if (ip%2) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
680:     if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
681:   }
682:   if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
683:     ctx->L = 16;
684:   } else {
685:     if (bs<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
686:     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");
687:     ctx->L = bs;
688:   }
689:   if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
690:     ctx->M = ctx->N/4;
691:   } else {
692:     if (ms<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
693:     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");
694:     ctx->M = ms;
695:   }
696:   if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
697:     ctx->num_subcomm = 1;
698:   } else {
699:     if (npart<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The npart argument must be > 0");
700:     ctx->num_subcomm = npart;
701:   }
702:   if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
703:     ctx->L = 256;
704:   } else {
705:     if (bsmax<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
706:     if (bsmax<ctx->L) ctx->L_max = ctx->L;
707:     else ctx->L_max = bsmax;
708:   }
709:   ctx->isreal = realmats;
710:   return(0);
711: }

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

718:    Logically Collective on NEP

720:    Input Parameters:
721: +  nep   - the eigenproblem solver context
722: .  ip    - number of integration points
723: .  bs    - block size
724: .  ms    - moment size
725: .  npart - number of partitions when splitting the communicator
726: .  bsmax - max block size
727: -  realmats - T(z) is real for real z

729:    Options Database Keys:
730: +  -nep_ciss_integration_points - Sets the number of integration points
731: .  -nep_ciss_blocksize - Sets the block size
732: .  -nep_ciss_moments - Sets the moment size
733: .  -nep_ciss_partitions - Sets the number of partitions
734: .  -nep_ciss_maxblocksize - Sets the maximum block size
735: -  -nep_ciss_realmats - T(z) is real for real z

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

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

746:    Level: advanced

748: .seealso: NEPCISSGetSizes()
749: @*/
750: PetscErrorCode NEPCISSSetSizes(NEP nep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
751: {

762:   PetscTryMethod(nep,"NEPCISSSetSizes_C",(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(nep,ip,bs,ms,npart,bsmax,realmats));
763:   return(0);
764: }

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

773:   if (ip) *ip = ctx->N;
774:   if (bs) *bs = ctx->L;
775:   if (ms) *ms = ctx->M;
776:   if (npart) *npart = ctx->num_subcomm;
777:   if (bsmax) *bsmax = ctx->L_max;
778:   if (realmats) *realmats = ctx->isreal;
779:   return(0);
780: }

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

787:    Not Collective

789:    Input Parameter:
790: .  nep - the eigenproblem solver context

792:    Output Parameters:
793: +  ip    - number of integration points
794: .  bs    - block size
795: .  ms    - moment size
796: .  npart - number of partitions when splitting the communicator
797: .  bsmax - max block size
798: -  realmats - T(z) is real for real z

800:    Level: advanced

802: .seealso: NEPCISSSetSizes()
803: @*/
804: PetscErrorCode NEPCISSGetSizes(NEP nep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
805: {

810:   PetscUseMethod(nep,"NEPCISSGetSizes_C",(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(nep,ip,bs,ms,npart,bsmax,realmats));
811:   return(0);
812: }

816: static PetscErrorCode NEPCISSSetThreshold_CISS(NEP nep,PetscReal delta,PetscReal spur)
817: {
818:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

821:   if (delta == PETSC_DEFAULT) {
822:     ctx->delta = 1e-12;
823:   } else {
824:     if (delta<=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
825:     ctx->delta = delta;
826:   }
827:   if (spur == PETSC_DEFAULT) {
828:     ctx->spurious_threshold = 1e-4;
829:   } else {
830:     if (spur<=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
831:     ctx->spurious_threshold = spur;
832:   }
833:   return(0);
834: }

838: /*@
839:    NEPCISSSetThreshold - Sets the values of various threshold parameters in
840:    the CISS solver.

842:    Logically Collective on NEP

844:    Input Parameters:
845: +  nep   - the eigenproblem solver context
846: .  delta - threshold for numerical rank
847: -  spur  - spurious threshold (to discard spurious eigenpairs)

849:    Options Database Keys:
850: +  -nep_ciss_delta - Sets the delta
851: -  -nep_ciss_spurious_threshold - Sets the spurious threshold

853:    Level: advanced

855: .seealso: NEPCISSGetThreshold()
856: @*/
857: PetscErrorCode NEPCISSSetThreshold(NEP nep,PetscReal delta,PetscReal spur)
858: {

865:   PetscTryMethod(nep,"NEPCISSSetThreshold_C",(NEP,PetscReal,PetscReal),(nep,delta,spur));
866:   return(0);
867: }

871: static PetscErrorCode NEPCISSGetThreshold_CISS(NEP nep,PetscReal *delta,PetscReal *spur)
872: {
873:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

876:   if (delta) *delta = ctx->delta;
877:   if (spur)  *spur = ctx->spurious_threshold;
878:   return(0);
879: }

883: /*@
884:    NEPCISSGetThreshold - Gets the values of various threshold parameters in
885:    the CISS solver.

887:    Not Collective

889:    Input Parameter:
890: .  nep - the eigenproblem solver context

892:    Output Parameters:
893: +  delta - threshold for numerical rank
894: -  spur  - spurious threshold (to discard spurious eigenpairs)

896:    Level: advanced

898: .seealso: NEPCISSSetThreshold()
899: @*/
900: PetscErrorCode NEPCISSGetThreshold(NEP nep,PetscReal *delta,PetscReal *spur)
901: {

906:   PetscUseMethod(nep,"NEPCISSGetThreshold_C",(NEP,PetscReal*,PetscReal*),(nep,delta,spur));
907:   return(0);
908: }

912: static PetscErrorCode NEPCISSSetRefinement_CISS(NEP nep,PetscInt inner,PetscInt blsize)
913: {
914:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

917:   if (inner == PETSC_DEFAULT) {
918:     ctx->refine_inner = 0;
919:   } else {
920:     if (inner<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
921:     ctx->refine_inner = inner;
922:   }
923:   if (blsize == PETSC_DEFAULT) {
924:     ctx->refine_blocksize = 0;
925:   } else {
926:     if (blsize<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
927:     ctx->refine_blocksize = blsize;
928:   }
929:   return(0);
930: }

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

938:    Logically Collective on NEP

940:    Input Parameters:
941: +  nep    - the eigenproblem solver context
942: .  inner  - number of iterative refinement iterations (inner loop)
943: -  blsize - number of iterative refinement iterations (blocksize loop)

945:    Options Database Keys:
946: +  -nep_ciss_refine_inner - Sets number of inner iterations
947: -  -nep_ciss_refine_blocksize - Sets number of blocksize iterations

949:    Level: advanced

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

961:   PetscTryMethod(nep,"NEPCISSSetRefinement_C",(NEP,PetscInt,PetscInt),(nep,inner,blsize));
962:   return(0);
963: }

967: static PetscErrorCode NEPCISSGetRefinement_CISS(NEP nep,PetscInt *inner,PetscInt *blsize)
968: {
969:   NEP_CISS *ctx = (NEP_CISS*)nep->data;

972:   if (inner)  *inner = ctx->refine_inner;
973:   if (blsize) *blsize = ctx->refine_blocksize;
974:   return(0);
975: }

979: /*@
980:    NEPCISSGetRefinement - Gets the values of various refinement parameters
981:    in the CISS solver.

983:    Not Collective

985:    Input Parameter:
986: .  nep - the eigenproblem solver context

988:    Output Parameters:
989: +  inner  - number of iterative refinement iterations (inner loop)
990: -  blsize - number of iterative refinement iterations (blocksize loop)

992:    Level: advanced

994: .seealso: NEPCISSSetRefinement()
995: @*/
996: PetscErrorCode NEPCISSGetRefinement(NEP nep, PetscInt *inner, PetscInt *blsize)
997: {

1002:   PetscUseMethod(nep,"NEPCISSGetRefinement_C",(NEP,PetscInt*,PetscInt*),(nep,inner,blsize));
1003:   return(0);
1004: }

1008: PetscErrorCode NEPReset_CISS(NEP nep)
1009: {
1011:   PetscInt       i;
1012:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;

1015:   PetscSubcommDestroy(&ctx->subcomm);
1016:   PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
1017:   BVDestroy(&ctx->S);
1018:   BVDestroy(&ctx->V);
1019:   BVDestroy(&ctx->Y);
1020:   if (!ctx->usest) {
1021:     for (i=0;i<ctx->num_solve_point;i++) {
1022:       KSPDestroy(&ctx->ksp[i]);
1023:     }
1024:     for (i=0;i<ctx->num_solve_point;i++) {
1025:       MatDestroy(&ctx->kspMat[i]);
1026:     }
1027:     PetscFree2(ctx->ksp,ctx->kspMat);
1028:   }
1029:   return(0);
1030: }

1034: PetscErrorCode NEPSetFromOptions_CISS(PetscOptionItems *PetscOptionsObject,NEP nep)
1035: {
1037:   PetscReal      r1,r2;
1038:   PetscInt       i1,i2,i3,i4,i5,i6,i7;
1039:   PetscBool      b1;

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

1044:   NEPCISSGetSizes(nep,&i1,&i2,&i3,&i4,&i5,&b1);
1045:   PetscOptionsInt("-nep_ciss_integration_points","CISS number of integration points","NEPCISSSetSizes",i1,&i1,NULL);
1046:   PetscOptionsInt("-nep_ciss_blocksize","CISS block size","NEPCISSSetSizes",i2,&i2,NULL);
1047:   PetscOptionsInt("-nep_ciss_moments","CISS moment size","NEPCISSSetSizes",i3,&i3,NULL);
1048:   PetscOptionsInt("-nep_ciss_partitions","CISS number of partitions","NEPCISSSetSizes",i4,&i4,NULL);
1049:   PetscOptionsInt("-nep_ciss_maxblocksize","CISS maximum block size","NEPCISSSetSizes",i5,&i5,NULL);
1050:   PetscOptionsBool("-nep_ciss_realmats","CISS flag indicating that T(z) is real for real z","NEPCISSSetSizes",b1,&b1,NULL);
1051:   NEPCISSSetSizes(nep,i1,i2,i3,i4,i5,b1);

1053:   NEPCISSGetThreshold(nep,&r1,&r2);
1054:   PetscOptionsReal("-nep_ciss_delta","CISS threshold for numerical rank","NEPCISSSetThreshold",r1,&r1,NULL);
1055:   PetscOptionsReal("-nep_ciss_spurious_threshold","CISS threshold for the spurious eigenpairs","NEPCISSSetThreshold",r2,&r2,NULL);
1056:   NEPCISSSetThreshold(nep,r1,r2);

1058:   NEPCISSGetRefinement(nep,&i6,&i7);
1059:   PetscOptionsInt("-nep_ciss_refine_inner","CISS number of inner iterative refinement iterations","NEPCISSSetRefinement",i6,&i6,NULL);
1060:   PetscOptionsInt("-nep_ciss_refine_blocksize","CISS number of blocksize iterative refinement iterations","NEPCISSSetRefinement",i7,&i7,NULL);
1061:   NEPCISSSetRefinement(nep,i6,i7);

1063:   PetscOptionsTail();
1064:   return(0);
1065: }

1069: PetscErrorCode NEPDestroy_CISS(NEP nep)
1070: {

1074:   PetscFree(nep->data);
1075:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NULL);
1076:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NULL);
1077:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NULL);
1078:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NULL);
1079:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NULL);
1080:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NULL);
1081:   return(0);
1082: }

1086: PetscErrorCode NEPView_CISS(NEP nep,PetscViewer viewer)
1087: {
1089:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;
1090:   PetscBool      isascii;

1093:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1094:   if (isascii) {
1095:     PetscViewerASCIIPrintf(viewer,"  CISS: 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);
1096:     if (ctx->isreal) {
1097:       PetscViewerASCIIPrintf(viewer,"  CISS: exploiting symmetry of integration points\n");
1098:     }
1099:     PetscViewerASCIIPrintf(viewer,"  CISS: threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold);
1100:     PetscViewerASCIIPrintf(viewer,"  CISS: iterative refinement  { inner: %D, blocksize: %D }\n",ctx->refine_inner, ctx->refine_blocksize);
1101:     PetscViewerASCIIPushTab(viewer);
1102:     if (!ctx->usest && ctx->ksp[0]) { KSPView(ctx->ksp[0],viewer); }
1103:     PetscViewerASCIIPopTab(viewer);
1104:   }
1105:   return(0);
1106: }

1110: PETSC_EXTERN PetscErrorCode NEPCreate_CISS(NEP nep)
1111: {
1113:   NEP_CISS       *ctx = (NEP_CISS*)nep->data;

1116:   PetscNewLog(nep,&ctx);
1117:   nep->data = ctx;
1118:   nep->ops->solve          = NEPSolve_CISS;
1119:   nep->ops->setup          = NEPSetUp_CISS;
1120:   nep->ops->setfromoptions = NEPSetFromOptions_CISS;
1121:   nep->ops->reset          = NEPReset_CISS;
1122:   nep->ops->destroy        = NEPDestroy_CISS;
1123:   nep->ops->view           = NEPView_CISS;
1124:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetSizes_C",NEPCISSSetSizes_CISS);
1125:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetSizes_C",NEPCISSGetSizes_CISS);
1126:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetThreshold_C",NEPCISSSetThreshold_CISS);
1127:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetThreshold_C",NEPCISSGetThreshold_CISS);
1128:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSSetRefinement_C",NEPCISSSetRefinement_CISS);
1129:   PetscObjectComposeFunction((PetscObject)nep,"NEPCISSGetRefinement_C",NEPCISSGetRefinement_CISS);
1130:   /* set default values of parameters */
1131:   ctx->N       = 32;
1132:   ctx->L       = 16;
1133:   ctx->M       = ctx->N/4;
1134:   ctx->delta   = 1e-12;
1135:   ctx->L_max   = 64;
1136:   ctx->spurious_threshold = 1e-4;
1137:   ctx->usest   = PETSC_FALSE;
1138:   ctx->isreal  = PETSC_FALSE;
1139:   ctx->refine_inner = 0;
1140:   ctx->refine_blocksize = 0;
1141:   ctx->num_subcomm = 1;
1142:   return(0);
1143: }