Actual source code: ciss.c

slepc-3.14.1 2020-12-08
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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/epsimpl.h>
 34: #include <slepcblaslapack.h>

 36: PetscLogEvent EPS_CISS_SVD;

 38: typedef struct {
 39:   /* parameters */
 40:   PetscInt          N;          /* number of integration points (32) */
 41:   PetscInt          L;          /* block size (16) */
 42:   PetscInt          M;          /* moment degree (N/4 = 4) */
 43:   PetscReal         delta;      /* threshold of singular value (1e-12) */
 44:   PetscInt          L_max;      /* maximum number of columns of the source matrix V */
 45:   PetscReal         spurious_threshold; /* discard spurious eigenpairs */
 46:   PetscBool         isreal;     /* A and B are real */
 47:   PetscInt          npart;      /* number of partitions */
 48:   PetscInt          refine_inner;
 49:   PetscInt          refine_blocksize;
 50:   /* private data */
 51:   PetscReal         *sigma;     /* threshold for numerical rank */
 52:   PetscInt          subcomm_id;
 53:   PetscInt          num_solve_point;
 54:   PetscScalar       *weight;
 55:   PetscScalar       *omega;
 56:   PetscScalar       *pp;
 57:   BV                V;
 58:   BV                S;
 59:   BV                pV;
 60:   BV                Y;
 61:   Vec               xsub;
 62:   Vec               xdup;
 63:   KSP               *ksp;       /* ksp array for storing factorizations at integration points */
 64:   PetscBool         useconj;
 65:   PetscReal         est_eig;
 66:   VecScatter        scatterin;
 67:   Mat               pA,pB;
 68:   PetscSubcomm      subcomm;
 69:   PetscBool         usest;
 70:   PetscBool         usest_set;  /* whether the user set the usest flag or not */
 71:   EPSCISSQuadRule   quad;
 72:   EPSCISSExtraction extraction;
 73: } EPS_CISS;

 75: /* destroy KSP objects when the number of solve points changes */
 76: PETSC_STATIC_INLINE PetscErrorCode EPSCISSResetSolvers(EPS eps)
 77: {
 79:   PetscInt       i;
 80:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

 83:   if (ctx->ksp) {
 84:     for (i=0;i<ctx->num_solve_point;i++) {
 85:       KSPDestroy(&ctx->ksp[i]);
 86:     }
 87:     PetscFree(ctx->ksp);
 88:   }
 89:   return(0);
 90: }

 92: /* clean PetscSubcomm object when the number of partitions changes */
 93: PETSC_STATIC_INLINE PetscErrorCode EPSCISSResetSubcomm(EPS eps)
 94: {
 96:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

 99:   EPSCISSResetSolvers(eps);
100:   PetscSubcommDestroy(&ctx->subcomm);
101:   return(0);
102: }

104: /* determine whether half of integration points can be avoided (use its conjugates);
105:    depends on isreal and the center of the region */
106: PETSC_STATIC_INLINE PetscErrorCode EPSCISSSetUseConj(EPS eps,PetscBool *useconj)
107: {
109:   PetscScalar    center;
110:   PetscReal      c,d;
111:   PetscBool      isellipse,isinterval;
112: #if defined(PETSC_USE_COMPLEX)
113:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
114: #endif

117:   *useconj = PETSC_FALSE;
118:   PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
119:   if (isellipse) {
120:     RGEllipseGetParameters(eps->rg,&center,NULL,NULL);
121: #if defined(PETSC_USE_COMPLEX)
122:     *useconj = (ctx->isreal && PetscImaginaryPart(center) == 0.0)? PETSC_TRUE: PETSC_FALSE;
123: #endif
124:   } else {
125:     PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval);
126:     if (isinterval) {
127:       RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d);
128: #if defined(PETSC_USE_COMPLEX)
129:       *useconj = (ctx->isreal && c==d)? PETSC_TRUE: PETSC_FALSE;
130: #endif
131:     }
132:   }
133:   return(0);
134: }

136: /* create PetscSubcomm object and determine num_solve_point (depends on npart, N, useconj) */
137: PETSC_STATIC_INLINE PetscErrorCode EPSCISSSetUpSubComm(EPS eps,PetscInt *num_solve_point)
138: {
140:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
141:   PetscInt       N = ctx->N;

144:   PetscSubcommCreate(PetscObjectComm((PetscObject)eps),&ctx->subcomm);
145:   PetscSubcommSetNumber(ctx->subcomm,ctx->npart);
146:   PetscSubcommSetType(ctx->subcomm,PETSC_SUBCOMM_INTERLACED);
147:   PetscLogObjectMemory((PetscObject)eps,sizeof(PetscSubcomm));
148:   ctx->subcomm_id = ctx->subcomm->color;
149:   EPSCISSSetUseConj(eps,&ctx->useconj);
150:   if (ctx->useconj) N = N/2;
151:   *num_solve_point = N / ctx->npart;
152:   if (N%ctx->npart > ctx->subcomm_id) (*num_solve_point)++;
153:   return(0);
154: }

156: static PetscErrorCode CISSRedundantMat(EPS eps)
157: {
159:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
160:   Mat            A,B;
161:   PetscInt       nmat;

164:   STGetNumMatrices(eps->st,&nmat);
165:   if (ctx->subcomm->n != 1) {
166:     STGetMatrix(eps->st,0,&A);
167:     MatDestroy(&ctx->pA);
168:     MatCreateRedundantMatrix(A,ctx->subcomm->n,PetscSubcommChild(ctx->subcomm),MAT_INITIAL_MATRIX,&ctx->pA);
169:     if (nmat>1) {
170:       STGetMatrix(eps->st,1,&B);
171:       MatDestroy(&ctx->pB);
172:       MatCreateRedundantMatrix(B,ctx->subcomm->n,PetscSubcommChild(ctx->subcomm),MAT_INITIAL_MATRIX,&ctx->pB);
173:     } else ctx->pB = NULL;
174:   } else {
175:     ctx->pA = NULL;
176:     ctx->pB = NULL;
177:   }
178:   return(0);
179: }

181: static PetscErrorCode CISSScatterVec(EPS eps)
182: {
184:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
185:   IS             is1,is2;
186:   Vec            v0;
187:   PetscInt       i,j,k,mstart,mend,mlocal;
188:   PetscInt       *idx1,*idx2,mloc_sub;

191:   VecDestroy(&ctx->xsub);
192:   MatCreateVecs(ctx->pA,&ctx->xsub,NULL);

194:   VecDestroy(&ctx->xdup);
195:   MatGetLocalSize(ctx->pA,&mloc_sub,NULL);
196:   VecCreateMPI(PetscSubcommContiguousParent(ctx->subcomm),mloc_sub,PETSC_DECIDE,&ctx->xdup);

198:   VecScatterDestroy(&ctx->scatterin);
199:   BVGetColumn(ctx->V,0,&v0);
200:   VecGetOwnershipRange(v0,&mstart,&mend);
201:   mlocal = mend - mstart;
202:   PetscMalloc2(ctx->subcomm->n*mlocal,&idx1,ctx->subcomm->n*mlocal,&idx2);
203:   j = 0;
204:   for (k=0;k<ctx->subcomm->n;k++) {
205:     for (i=mstart;i<mend;i++) {
206:       idx1[j]   = i;
207:       idx2[j++] = i + eps->n*k;
208:     }
209:   }
210:   ISCreateGeneral(PetscObjectComm((PetscObject)eps),ctx->subcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);
211:   ISCreateGeneral(PetscObjectComm((PetscObject)eps),ctx->subcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);
212:   VecScatterCreate(v0,is1,ctx->xdup,is2,&ctx->scatterin);
213:   ISDestroy(&is1);
214:   ISDestroy(&is2);
215:   PetscFree2(idx1,idx2);
216:   BVRestoreColumn(ctx->V,0,&v0);
217:   return(0);
218: }

220: static PetscErrorCode SetPathParameter(EPS eps)
221: {
223:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
224:   PetscInt       i,j;
225:   PetscScalar    center=0.0,tmp,tmp2,*omegai;
226:   PetscReal      theta,radius=1.0,vscale,a,b,c,d,max_w=0.0,rgscale;
227: #if defined(PETSC_USE_COMPLEX)
228:   PetscReal      start_ang,end_ang;
229: #endif
230:   PetscBool      isring=PETSC_FALSE,isellipse=PETSC_FALSE,isinterval=PETSC_FALSE;

233:   PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
234:   PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring);
235:   PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval);
236:   RGGetScale(eps->rg,&rgscale);
237:   PetscMalloc1(ctx->N+1l,&omegai);
238:   RGComputeContour(eps->rg,ctx->N,ctx->omega,omegai);
239:   if (isellipse) {
240:     RGEllipseGetParameters(eps->rg,&center,&radius,&vscale);
241:     for (i=0;i<ctx->N;i++) {
242: #if defined(PETSC_USE_COMPLEX)
243:       theta = 2.0*PETSC_PI*(i+0.5)/ctx->N;
244:       ctx->pp[i] = PetscCMPLX(PetscCosReal(theta),vscale*PetscSinReal(theta));
245:       ctx->weight[i] = rgscale*radius*(PetscCMPLX(vscale*PetscCosReal(theta),PetscSinReal(theta)))/(PetscReal)ctx->N;
246: #else
247:       theta = (PETSC_PI/ctx->N)*(i+0.5);
248:       ctx->pp[i] = PetscCosReal(theta);
249:       ctx->weight[i] = PetscCosReal((ctx->N-1)*theta)/ctx->N;
250:       ctx->omega[i] = rgscale*(center + radius*ctx->pp[i]);
251: #endif
252:     }
253:   } else if (ctx->quad == EPS_CISS_QUADRULE_CHEBYSHEV) {
254:     for (i=0;i<ctx->N;i++) {
255:       theta = (PETSC_PI/ctx->N)*(i+0.5);
256:       ctx->pp[i] = PetscCosReal(theta);
257:       ctx->weight[i] = PetscCosReal((ctx->N-1)*theta)/ctx->N;
258:     }
259:     if (isinterval) {
260:       RGIntervalGetEndpoints(eps->rg,&a,&b,&c,&d);
261:       if ((c!=d || c!=0.0) && (a!=b || a!=0.0)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Endpoints of the imaginary axis or the real axis must be both zero");
262:       for (i=0;i<ctx->N;i++) {
263:         if (c==d) ctx->omega[i] = ((b-a)*(ctx->pp[i]+1.0)/2.0+a)*rgscale;
264:         if (a==b) {
265: #if defined(PETSC_USE_COMPLEX)
266:           ctx->omega[i] = ((d-c)*(ctx->pp[i]+1.0)/2.0+c)*rgscale*PETSC_i;
267: #else
268:           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Integration points on a vertical line require complex arithmetic");
269: #endif
270:         }
271:       }
272:     }
273:     if (isring) {  /* only supported in complex scalars */
274: #if defined(PETSC_USE_COMPLEX)
275:       RGRingGetParameters(eps->rg,&center,&radius,&vscale,&start_ang,&end_ang,NULL);
276:       for (i=0;i<ctx->N;i++) {
277:         theta = (start_ang*2.0+(end_ang-start_ang)*(PetscRealPart(ctx->pp[i])+1.0))*PETSC_PI;
278:         ctx->omega[i] = rgscale*(center + radius*PetscCMPLX(PetscCosReal(theta),vscale*PetscSinReal(theta)));
279:       }
280: #endif
281:     }
282:   } else {
283:     if (isinterval) {
284:       RGIntervalGetEndpoints(eps->rg,&a,&b,&c,&d);
285:       center = rgscale*((b+a)/2.0+(d+c)/2.0*PETSC_PI);
286:       radius = PetscSqrtReal(PetscPowRealInt(rgscale*(b-a)/2.0,2)+PetscPowRealInt(rgscale*(d-c)/2.0,2));
287:     } else if (isring) {
288:       RGRingGetParameters(eps->rg,&center,&radius,NULL,NULL,NULL,NULL);
289:       center *= rgscale;
290:       radius *= rgscale;
291:     }
292:     for (i=0;i<ctx->N;i++) {
293:       ctx->pp[i] = (ctx->omega[i]-center)/radius;
294:       tmp = 1; tmp2 = 1;
295:       for (j=0;j<ctx->N;j++) {
296:         tmp *= ctx->omega[j];
297:         if (i != j) tmp2 *= ctx->omega[j]-ctx->omega[i];
298:       }
299:       ctx->weight[i] = tmp/tmp2;
300:       max_w = PetscMax(PetscAbsScalar(ctx->weight[i]),max_w);
301:     }
302:     for (i=0;i<ctx->N;i++) ctx->weight[i] /= (PetscScalar)max_w;
303:   }
304:   PetscFree(omegai);
305:   return(0);
306: }

308: static PetscErrorCode CISSVecSetRandom(BV V,PetscInt i0,PetscInt i1)
309: {
311:   PetscInt       i,j,nlocal;
312:   PetscScalar    *vdata;
313:   Vec            x;

316:   BVGetSizes(V,&nlocal,NULL,NULL);
317:   for (i=i0;i<i1;i++) {
318:     BVSetRandomColumn(V,i);
319:     BVGetColumn(V,i,&x);
320:     VecGetArray(x,&vdata);
321:     for (j=0;j<nlocal;j++) {
322:       vdata[j] = PetscRealPart(vdata[j]);
323:       if (PetscRealPart(vdata[j]) < 0.5) vdata[j] = -1.0;
324:       else vdata[j] = 1.0;
325:     }
326:     VecRestoreArray(x,&vdata);
327:     BVRestoreColumn(V,i,&x);
328:   }
329:   return(0);
330: }

332: static PetscErrorCode VecScatterVecs(EPS eps,BV Vin,PetscInt n)
333: {
334:   PetscErrorCode    ierr;
335:   EPS_CISS          *ctx = (EPS_CISS*)eps->data;
336:   PetscInt          i;
337:   Vec               vi,pvi;
338:   const PetscScalar *array;

341:   for (i=0;i<n;i++) {
342:     BVGetColumn(Vin,i,&vi);
343:     VecScatterBegin(ctx->scatterin,vi,ctx->xdup,INSERT_VALUES,SCATTER_FORWARD);
344:     VecScatterEnd(ctx->scatterin,vi,ctx->xdup,INSERT_VALUES,SCATTER_FORWARD);
345:     BVRestoreColumn(Vin,i,&vi);
346:     VecGetArrayRead(ctx->xdup,&array);
347:     VecPlaceArray(ctx->xsub,array);
348:     BVGetColumn(ctx->pV,i,&pvi);
349:     VecCopy(ctx->xsub,pvi);
350:     BVRestoreColumn(ctx->pV,i,&pvi);
351:     VecResetArray(ctx->xsub);
352:     VecRestoreArrayRead(ctx->xdup,&array);
353:   }
354:   return(0);
355: }

357: static PetscErrorCode SolveLinearSystem(EPS eps,Mat A,Mat B,BV V,PetscInt L_start,PetscInt L_end,PetscBool initksp)
358: {
360:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
361:   PetscInt       i,p_id;
362:   Mat            Fz,kspMat,MV,BMV=NULL,MC;
363:   KSP            ksp;

366:   if (!ctx->ksp) { EPSCISSGetKSPs(eps,&ctx->num_solve_point,&ctx->ksp); }
367:   if (ctx->usest) {
368:     MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Fz);
369:   }
370:   BVSetActiveColumns(V,L_start,L_end);
371:   BVGetMat(V,&MV);
372:   if (B) {
373:     MatProductCreate(B,MV,NULL,&BMV);
374:     MatProductSetType(BMV,MATPRODUCT_AB);
375:     MatProductSetFromOptions(BMV);
376:     MatProductSymbolic(BMV);
377:   }
378:   for (i=0;i<ctx->num_solve_point;i++) {
379:     p_id = i*ctx->subcomm->n + ctx->subcomm_id;
380:     if (!ctx->usest && initksp) {
381:       MatDuplicate(A,MAT_COPY_VALUES,&kspMat);
382:       if (B) {
383:         MatAXPY(kspMat,-ctx->omega[p_id],B,DIFFERENT_NONZERO_PATTERN);
384:       } else {
385:         MatShift(kspMat,-ctx->omega[p_id]);
386:       }
387:       KSPSetOperators(ctx->ksp[i],kspMat,kspMat);
388:       MatDestroy(&kspMat);
389:     } else if (ctx->usest) {
390:       STSetShift(eps->st,ctx->omega[p_id]);
391:       STGetKSP(eps->st,&ksp);
392:     }
393:     BVSetActiveColumns(ctx->Y,i*ctx->L_max+L_start,i*ctx->L_max+L_end);
394:     BVGetMat(ctx->Y,&MC);
395:     if (B) {
396:       MatProductNumeric(BMV);
397:       if (ctx->usest) {
398:         KSPMatSolve(ksp,BMV,MC);
399:       } else {
400:         KSPMatSolve(ctx->ksp[i],BMV,MC);
401:       }
402:     } else {
403:       if (ctx->usest) {
404:         KSPMatSolve(ksp,MV,MC);
405:       } else {
406:         KSPMatSolve(ctx->ksp[i],MV,MC);
407:       }
408:     }
409:     if (ctx->usest && i<ctx->num_solve_point-1) { KSPReset(ksp); }
410:     BVRestoreMat(ctx->Y,&MC);
411:   }
412:   MatDestroy(&BMV);
413:   BVRestoreMat(V,&MV);
414:   if (ctx->usest) { MatDestroy(&Fz); }
415:   return(0);
416: }

418: #if defined(PETSC_USE_COMPLEX)
419: static PetscErrorCode EstimateNumberEigs(EPS eps,PetscInt *L_add)
420: {
422:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
423:   PetscInt       i,j,p_id;
424:   PetscScalar    tmp,m = 1,sum = 0.0;
425:   PetscReal      eta;
426:   Vec            v,vtemp,vj;

429:   BVCreateVec(ctx->Y,&v);
430:   BVCreateVec(ctx->V,&vtemp);
431:   for (j=0;j<ctx->L;j++) {
432:     VecSet(v,0);
433:     for (i=0;i<ctx->num_solve_point; i++) {
434:       p_id = i*ctx->subcomm->n + ctx->subcomm_id;
435:       BVSetActiveColumns(ctx->Y,i*ctx->L_max+j,i*ctx->L_max+j+1);
436:       BVMultVec(ctx->Y,ctx->weight[p_id],1,v,&m);
437:     }
438:     BVGetColumn(ctx->V,j,&vj);
439:     if (ctx->pA) {
440:       VecSet(vtemp,0);
441:       VecScatterBegin(ctx->scatterin,v,vtemp,ADD_VALUES,SCATTER_REVERSE);
442:       VecScatterEnd(ctx->scatterin,v,vtemp,ADD_VALUES,SCATTER_REVERSE);
443:       VecDot(vj,vtemp,&tmp);
444:     } else {
445:       VecDot(vj,v,&tmp);
446:     }
447:     BVRestoreColumn(ctx->V,j,&vj);
448:     if (ctx->useconj) sum += PetscRealPart(tmp)*2;
449:     else sum += tmp;
450:   }
451:   ctx->est_eig = PetscAbsScalar(sum/(PetscReal)ctx->L);
452:   eta = PetscPowReal(10.0,-PetscLog10Real(eps->tol)/ctx->N);
453:   PetscInfo1(eps,"Estimation_#Eig %f\n",(double)ctx->est_eig);
454:   *L_add = (PetscInt)PetscCeilReal((ctx->est_eig*eta)/ctx->M) - ctx->L;
455:   if (*L_add < 0) *L_add = 0;
456:   if (*L_add>ctx->L_max-ctx->L) {
457:     PetscInfo(eps,"Number of eigenvalues around the contour path may be too large\n");
458:     *L_add = ctx->L_max-ctx->L;
459:   }
460:   VecDestroy(&v);
461:   VecDestroy(&vtemp);
462:   return(0);
463: }
464: #endif

466: static PetscErrorCode CalcMu(EPS eps,PetscScalar *Mu)
467: {
469:   PetscMPIInt    sub_size,len;
470:   PetscInt       i,j,k,s;
471:   PetscScalar    *m,*temp,*temp2,*ppk,alp;
472:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
473:   Mat            M;

476:   MPI_Comm_size(PetscSubcommChild(ctx->subcomm),&sub_size);
477:   PetscMalloc3(ctx->num_solve_point*ctx->L*(ctx->L+1),&temp,2*ctx->M*ctx->L*ctx->L,&temp2,ctx->num_solve_point,&ppk);
478:   MatCreateSeqDense(PETSC_COMM_SELF,ctx->L,ctx->L_max*ctx->num_solve_point,NULL,&M);
479:   for (i=0;i<2*ctx->M*ctx->L*ctx->L;i++) temp2[i] = 0;
480:   BVSetActiveColumns(ctx->Y,0,ctx->L_max*ctx->num_solve_point);
481:   if (ctx->pA) {
482:     BVSetActiveColumns(ctx->pV,0,ctx->L);
483:     BVDot(ctx->Y,ctx->pV,M);
484:   } else {
485:     BVSetActiveColumns(ctx->V,0,ctx->L);
486:     BVDot(ctx->Y,ctx->V,M);
487:   }
488:   MatDenseGetArray(M,&m);
489:   for (i=0;i<ctx->num_solve_point;i++) {
490:     for (j=0;j<ctx->L;j++) {
491:       for (k=0;k<ctx->L;k++) {
492:         temp[k+j*ctx->L+i*ctx->L*ctx->L]=m[k+j*ctx->L+i*ctx->L*ctx->L_max];
493:       }
494:     }
495:   }
496:   MatDenseRestoreArray(M,&m);
497:   for (i=0;i<ctx->num_solve_point;i++) ppk[i] = 1;
498:   for (k=0;k<2*ctx->M;k++) {
499:     for (j=0;j<ctx->L;j++) {
500:       for (i=0;i<ctx->num_solve_point;i++) {
501:         alp = ppk[i]*ctx->weight[i*ctx->subcomm->n + ctx->subcomm_id];
502:         for (s=0;s<ctx->L;s++) {
503:           if (ctx->useconj) temp2[s+(j+k*ctx->L)*ctx->L] += PetscRealPart(alp*temp[s+(j+i*ctx->L)*ctx->L])*2;
504:           else temp2[s+(j+k*ctx->L)*ctx->L] += alp*temp[s+(j+i*ctx->L)*ctx->L];
505:         }
506:       }
507:     }
508:     for (i=0;i<ctx->num_solve_point;i++)
509:       ppk[i] *= ctx->pp[i*ctx->subcomm->n + ctx->subcomm_id];
510:   }
511:   for (i=0;i<2*ctx->M*ctx->L*ctx->L;i++) temp2[i] /= sub_size;
512:   PetscMPIIntCast(2*ctx->M*ctx->L*ctx->L,&len);
513:   MPIU_Allreduce(temp2,Mu,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)eps));
514:   PetscFree3(temp,temp2,ppk);
515:   MatDestroy(&M);
516:   return(0);
517: }

519: static PetscErrorCode BlockHankel(EPS eps,PetscScalar *Mu,PetscInt s,PetscScalar *H)
520: {
521:   EPS_CISS *ctx = (EPS_CISS*)eps->data;
522:   PetscInt i,j,k,L=ctx->L,M=ctx->M;

525:   for (k=0;k<L*M;k++)
526:     for (j=0;j<M;j++)
527:       for (i=0;i<L;i++)
528:         H[j*L+i+k*L*M] = Mu[i+k*L+(j+s)*L*L];
529:   return(0);
530: }

532: static PetscErrorCode SVD_H0(EPS eps,PetscScalar *S,PetscInt *K)
533: {
535:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
536:   PetscInt       i,ml=ctx->L*ctx->M;
537:   PetscBLASInt   m,n,lda,ldu,ldvt,lwork,info;
538:   PetscScalar    *work;
539: #if defined(PETSC_USE_COMPLEX)
540:   PetscReal      *rwork;
541: #endif

544:   PetscLogEventBegin(EPS_CISS_SVD,eps,0,0,0);
545:   PetscMalloc1(5*ml,&work);
546: #if defined(PETSC_USE_COMPLEX)
547:   PetscMalloc1(5*ml,&rwork);
548: #endif
549:   PetscBLASIntCast(ml,&m);
550:   n = m; lda = m; ldu = m; ldvt = m; lwork = 5*m;
551:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
552: #if defined(PETSC_USE_COMPLEX)
553:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,S,&lda,ctx->sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,rwork,&info));
554: #else
555:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,S,&lda,ctx->sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,&info));
556: #endif
557:   SlepcCheckLapackInfo("gesvd",info);
558:   PetscFPTrapPop();
559:   (*K) = 0;
560:   for (i=0;i<ml;i++) {
561:     if (ctx->sigma[i]/PetscMax(ctx->sigma[0],1)>ctx->delta) (*K)++;
562:   }
563:   PetscFree(work);
564: #if defined(PETSC_USE_COMPLEX)
565:   PetscFree(rwork);
566: #endif
567:   PetscLogEventEnd(EPS_CISS_SVD,eps,0,0,0);
568:   return(0);
569: }

571: static PetscErrorCode ConstructS(EPS eps)
572: {
574:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
575:   PetscInt       i,j,k,vec_local_size,p_id;
576:   Vec            v,sj;
577:   PetscScalar    *ppk, *v_data, m = 1;

580:   BVGetSizes(ctx->Y,&vec_local_size,NULL,NULL);
581:   PetscMalloc1(ctx->num_solve_point,&ppk);
582:   for (i=0;i<ctx->num_solve_point;i++) ppk[i] = 1;
583:   BVCreateVec(ctx->Y,&v);
584:   for (k=0;k<ctx->M;k++) {
585:     for (j=0;j<ctx->L;j++) {
586:       VecSet(v,0);
587:       for (i=0;i<ctx->num_solve_point;i++) {
588:         p_id = i*ctx->subcomm->n + ctx->subcomm_id;
589:         BVSetActiveColumns(ctx->Y,i*ctx->L_max+j,i*ctx->L_max+j+1);
590:         BVMultVec(ctx->Y,ppk[i]*ctx->weight[p_id],1.0,v,&m);
591:       }
592:       if (ctx->useconj) {
593:         VecGetArray(v,&v_data);
594:         for (i=0;i<vec_local_size;i++) v_data[i] = PetscRealPart(v_data[i])*2;
595:         VecRestoreArray(v,&v_data);
596:       }
597:       BVGetColumn(ctx->S,k*ctx->L+j,&sj);
598:       if (ctx->pA) {
599:         VecSet(sj,0);
600:         VecScatterBegin(ctx->scatterin,v,sj,ADD_VALUES,SCATTER_REVERSE);
601:         VecScatterEnd(ctx->scatterin,v,sj,ADD_VALUES,SCATTER_REVERSE);
602:       } else {
603:         VecCopy(v,sj);
604:       }
605:       BVRestoreColumn(ctx->S,k*ctx->L+j,&sj);
606:     }
607:     for (i=0;i<ctx->num_solve_point;i++) {
608:       p_id = i*ctx->subcomm->n + ctx->subcomm_id;
609:       ppk[i] *= ctx->pp[p_id];
610:     }
611:   }
612:   PetscFree(ppk);
613:   VecDestroy(&v);
614:   return(0);
615: }

617: static PetscErrorCode SVD_S(BV S,PetscInt ml,PetscReal delta,PetscReal *sigma,PetscInt *K)
618: {
620:   PetscInt       i,j,k,local_size;
621:   PetscMPIInt    len;
622:   PetscScalar    *work,*temp,*B,*tempB,*s_data,*Q1,*Q2,*temp2,alpha=1,beta=0;
623:   PetscBLASInt   l,m,n,lda,ldu,ldvt,lwork,info,ldb,ldc;
624: #if defined(PETSC_USE_COMPLEX)
625:   PetscReal      *rwork;
626: #endif

629:   BVGetSizes(S,&local_size,NULL,NULL);
630:   BVGetArray(S,&s_data);
631:   PetscMalloc7(ml*ml,&temp,ml*ml,&temp2,local_size*ml,&Q1,local_size*ml,&Q2,ml*ml,&B,ml*ml,&tempB,5*ml,&work);
632:   PetscArrayzero(B,ml*ml);
633: #if defined(PETSC_USE_COMPLEX)
634:   PetscMalloc1(5*ml,&rwork);
635: #endif
636:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);

638:   for (i=0;i<ml;i++) B[i*ml+i]=1;

640:   for (k=0;k<2;k++) {
641:     PetscBLASIntCast(local_size,&m);
642:     PetscBLASIntCast(ml,&l);
643:     n = l; lda = m; ldb = m; ldc = l;
644:     if (k == 0) {
645:       PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&l,&n,&m,&alpha,s_data,&lda,s_data,&ldb,&beta,temp,&ldc));
646:     } else if ((k%2)==1) {
647:       PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&l,&n,&m,&alpha,Q1,&lda,Q1,&ldb,&beta,temp,&ldc));
648:     } else {
649:       PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&l,&n,&m,&alpha,Q2,&lda,Q2,&ldb,&beta,temp,&ldc));
650:     }
651:     PetscArrayzero(temp2,ml*ml);
652:     PetscMPIIntCast(ml*ml,&len);
653:     MPIU_Allreduce(temp,temp2,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)S));

655:     PetscBLASIntCast(ml,&m);
656:     n = m; lda = m; lwork = 5*m, ldu = 1; ldvt = 1;
657: #if defined(PETSC_USE_COMPLEX)
658:     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&m,&n,temp2,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,rwork,&info));
659: #else
660:     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&m,&n,temp2,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,&info));
661: #endif
662:     SlepcCheckLapackInfo("gesvd",info);

664:     PetscBLASIntCast(local_size,&l);
665:     PetscBLASIntCast(ml,&n);
666:     m = n; lda = l; ldb = m; ldc = l;
667:     if (k==0) {
668:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&l,&n,&m,&alpha,s_data,&lda,temp2,&ldb,&beta,Q1,&ldc));
669:     } else if ((k%2)==1) {
670:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&l,&n,&m,&alpha,Q1,&lda,temp2,&ldb,&beta,Q2,&ldc));
671:     } else {
672:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&l,&n,&m,&alpha,Q2,&lda,temp2,&ldb,&beta,Q1,&ldc));
673:     }

675:     PetscBLASIntCast(ml,&l);
676:     m = l; n = l; lda = l; ldb = m; ldc = l;
677:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&l,&n,&m,&alpha,B,&lda,temp2,&ldb,&beta,tempB,&ldc));
678:     for (i=0;i<ml;i++) {
679:       sigma[i] = sqrt(sigma[i]);
680:       for (j=0;j<local_size;j++) {
681:         if ((k%2)==1) Q2[j+i*local_size]/=sigma[i];
682:         else Q1[j+i*local_size]/=sigma[i];
683:       }
684:       for (j=0;j<ml;j++) {
685:         B[j+i*ml]=tempB[j+i*ml]*sigma[i];
686:       }
687:     }
688:   }

690:   PetscBLASIntCast(ml,&m);
691:   n = m; lda = m; ldu=1; ldvt=1;
692: #if defined(PETSC_USE_COMPLEX)
693:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","O",&m,&n,B,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,rwork,&info));
694: #else
695:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","O",&m,&n,B,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,&info));
696: #endif
697:   SlepcCheckLapackInfo("gesvd",info);

699:   PetscBLASIntCast(local_size,&l);
700:   PetscBLASIntCast(ml,&n);
701:   m = n; lda = l; ldb = m; ldc = l;
702:   if ((k%2)==1) {
703:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&l,&n,&m,&alpha,Q1,&lda,B,&ldb,&beta,s_data,&ldc));
704:   } else {
705:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&l,&n,&m,&alpha,Q2,&lda,B,&ldb,&beta,s_data,&ldc));
706:   }

708:   PetscFPTrapPop();
709:   BVRestoreArray(S,&s_data);

711:   (*K) = 0;
712:   for (i=0;i<ml;i++) {
713:     if (sigma[i]/PetscMax(sigma[0],1)>delta) (*K)++;
714:   }
715:   PetscFree7(temp,temp2,Q1,Q2,B,tempB,work);
716: #if defined(PETSC_USE_COMPLEX)
717:   PetscFree(rwork);
718: #endif
719:   return(0);
720: }

722: static PetscErrorCode isGhost(EPS eps,PetscInt ld,PetscInt nv,PetscBool *fl)
723: {
725:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
726:   PetscInt       i,j;
727:   PetscScalar    *pX;
728:   PetscReal      *tau,s1,s2,tau_max=0.0;

731:   PetscMalloc1(nv,&tau);
732:   DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
733:   DSGetArray(eps->ds,DS_MAT_X,&pX);

735:   for (i=0;i<nv;i++) {
736:     s1 = 0;
737:     s2 = 0;
738:     for (j=0;j<nv;j++) {
739:       s1 += PetscAbsScalar(PetscPowScalarInt(pX[i*ld+j],2));
740:       s2 += PetscPowRealInt(PetscAbsScalar(pX[i*ld+j]),2)/ctx->sigma[j];
741:     }
742:     tau[i] = s1/s2;
743:     tau_max = PetscMax(tau_max,tau[i]);
744:   }
745:   DSRestoreArray(eps->ds,DS_MAT_X,&pX);
746:   for (i=0;i<nv;i++) {
747:     tau[i] /= tau_max;
748:   }
749:   for (i=0;i<nv;i++) {
750:     if (tau[i]>=ctx->spurious_threshold) fl[i] = PETSC_TRUE;
751:     else fl[i] = PETSC_FALSE;
752:   }
753:   PetscFree(tau);
754:   return(0);
755: }

757: static PetscErrorCode rescale_eig(EPS eps,PetscInt nv)
758: {
760:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
761:   PetscInt       i;
762:   PetscScalar    center;
763:   PetscReal      radius,a,b,c,d,rgscale;
764: #if defined(PETSC_USE_COMPLEX)
765:   PetscReal      start_ang,end_ang,vscale,theta;
766: #endif
767:   PetscBool      isring,isellipse,isinterval;

770:   PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
771:   PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring);
772:   PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval);
773:   RGGetScale(eps->rg,&rgscale);
774:   if (isinterval) {
775:     RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d);
776:     if (c==d) {
777:       for (i=0;i<nv;i++) {
778: #if defined(PETSC_USE_COMPLEX)
779:         eps->eigr[i] = PetscRealPart(eps->eigr[i]);
780: #else
781:         eps->eigi[i] = 0;
782: #endif
783:       }
784:     }
785:   }
786:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
787:     if (isellipse) {
788:       RGEllipseGetParameters(eps->rg,&center,&radius,NULL);
789:       for (i=0;i<nv;i++) eps->eigr[i] = rgscale*(center + radius*eps->eigr[i]);
790:     } else if (isinterval) {
791:       RGIntervalGetEndpoints(eps->rg,&a,&b,&c,&d);
792:       if (ctx->quad == EPS_CISS_QUADRULE_CHEBYSHEV) {
793:         for (i=0;i<nv;i++) {
794:           if (c==d) eps->eigr[i] = ((b-a)*(eps->eigr[i]+1.0)/2.0+a)*rgscale;
795:           if (a==b) {
796: #if defined(PETSC_USE_COMPLEX)
797:             eps->eigr[i] = ((d-c)*(eps->eigr[i]+1.0)/2.0+c)*rgscale*PETSC_i;
798: #else
799:             SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Integration points on a vertical line require complex arithmetic");
800: #endif
801:           }
802:         }
803:       } else {
804:         center = (b+a)/2.0+(d+c)/2.0*PETSC_PI;
805:         radius = PetscSqrtReal(PetscPowRealInt((b-a)/2.0,2)+PetscPowRealInt((d-c)/2.0,2));
806:         for (i=0;i<nv;i++) eps->eigr[i] = center + radius*eps->eigr[i];
807:       }
808:     } else if (isring) {  /* only supported in complex scalars */
809: #if defined(PETSC_USE_COMPLEX)
810:       RGRingGetParameters(eps->rg,&center,&radius,&vscale,&start_ang,&end_ang,NULL);
811:       if (ctx->quad == EPS_CISS_QUADRULE_CHEBYSHEV) {
812:         for (i=0;i<nv;i++) {
813:           theta = (start_ang*2.0+(end_ang-start_ang)*(PetscRealPart(eps->eigr[i])+1.0))*PETSC_PI;
814:           eps->eigr[i] = rgscale*center + (rgscale*radius+PetscImaginaryPart(eps->eigr[i]))*PetscCMPLX(PetscCosReal(theta),vscale*PetscSinReal(theta));
815:         }
816:       } else {
817:         for (i=0;i<nv;i++) eps->eigr[i] = rgscale*(center + radius*eps->eigr[i]);
818:       }
819: #endif
820:     }
821:   }
822:   return(0);
823: }

825: PetscErrorCode EPSSetUp_CISS(EPS eps)
826: {
828:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
829:   PetscBool      istrivial,isring,isellipse,isinterval,flg,useconj;
830:   PetscReal      c,d;
831:   Mat            A;

834:   if (eps->ncv==PETSC_DEFAULT) {
835:     eps->ncv = ctx->L_max*ctx->M;
836:     if (eps->ncv>eps->n) {
837:       eps->ncv = eps->n;
838:       ctx->L_max = eps->ncv/ctx->M;
839:       if (!ctx->L_max) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cannot adjust solver parameters, try setting a smaller value of M (moment size)");
840:     }
841:   } else {
842:     EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
843:     ctx->L_max = eps->ncv/ctx->M;
844:     if (!ctx->L_max) {
845:       ctx->L_max = 1;
846:       eps->ncv = ctx->L_max*ctx->M;
847:     }
848:   }
849:   ctx->L = PetscMin(ctx->L,ctx->L_max);
850:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
851:   if (eps->mpd==PETSC_DEFAULT) eps->mpd = eps->ncv;
852:   if (!eps->which) eps->which = EPS_ALL;
853:   if (eps->which!=EPS_ALL) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only computing all eigenvalues");
854:   EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_EXTRACTION | EPS_FEATURE_STOPPING | EPS_FEATURE_TWOSIDED);

856:   /* check region */
857:   RGIsTrivial(eps->rg,&istrivial);
858:   if (istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
859:   RGGetComplement(eps->rg,&flg);
860:   if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
861:   PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
862:   PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring);
863:   PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval);
864:   if (!isellipse && !isring && !isinterval) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Currently only implemented for interval, elliptic or ring regions");
865:   /* if useconj has changed, then reset subcomm data */
866:   EPSCISSSetUseConj(eps,&useconj);
867:   if (useconj!=ctx->useconj) { EPSCISSResetSubcomm(eps); }

869: #if !defined(PETSC_USE_COMPLEX)
870:   if (isring) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Ring region only supported for complex scalars");
871: #endif
872:   if (isinterval) {
873:     RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d);
874: #if !defined(PETSC_USE_COMPLEX)
875:     if (c!=d || c!=0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"In real scalars, endpoints of the imaginary axis must be both zero");
876: #endif
877:     if (!ctx->quad && c==d) ctx->quad = EPS_CISS_QUADRULE_CHEBYSHEV;
878:   }
879:   if (!ctx->quad) ctx->quad = EPS_CISS_QUADRULE_TRAPEZOIDAL;

881:   /* create split comm */
882:   if (!ctx->subcomm) { EPSCISSSetUpSubComm(eps,&ctx->num_solve_point); }

884:   EPSAllocateSolution(eps,0);
885:   if (ctx->weight) { PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma); }
886:   PetscMalloc4(ctx->N,&ctx->weight,ctx->N+1,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma);
887:   PetscLogObjectMemory((PetscObject)eps,3*ctx->N*sizeof(PetscScalar)+ctx->L_max*ctx->N*sizeof(PetscReal));

889:   /* allocate basis vectors */
890:   BVDestroy(&ctx->S);
891:   BVDuplicateResize(eps->V,ctx->L_max*ctx->M,&ctx->S);
892:   PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->S);
893:   BVDestroy(&ctx->V);
894:   BVDuplicateResize(eps->V,ctx->L_max,&ctx->V);
895:   PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->V);

897:   STGetMatrix(eps->st,0,&A);
898:   PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);
899:   if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Matrix type shell is not supported in this solver");

901:   if (!ctx->usest_set) ctx->usest = (ctx->npart>1)? PETSC_FALSE: PETSC_TRUE;
902:   if (ctx->usest && ctx->npart>1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The usest flag is not supported when partitions > 1");

904:   CISSRedundantMat(eps);
905:   if (ctx->pA) {
906:     CISSScatterVec(eps);
907:     BVDestroy(&ctx->pV);
908:     BVCreate(PetscObjectComm((PetscObject)ctx->xsub),&ctx->pV);
909:     BVSetSizesFromVec(ctx->pV,ctx->xsub,eps->n);
910:     BVSetFromOptions(ctx->pV);
911:     BVResize(ctx->pV,ctx->L_max,PETSC_FALSE);
912:     PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->pV);
913:   }

915:   EPSCheckSinvertCondition(eps,ctx->usest," (with the usest flag set)");

917:   BVDestroy(&ctx->Y);
918:   if (ctx->pA) {
919:     BVCreate(PetscObjectComm((PetscObject)ctx->xsub),&ctx->Y);
920:     BVSetSizesFromVec(ctx->Y,ctx->xsub,eps->n);
921:     BVSetFromOptions(ctx->Y);
922:     BVResize(ctx->Y,ctx->num_solve_point*ctx->L_max,PETSC_FALSE);
923:   } else {
924:     BVDuplicateResize(eps->V,ctx->num_solve_point*ctx->L_max,&ctx->Y);
925:   }
926:   PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->Y);

928:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
929:     DSSetType(eps->ds,DSGNHEP);
930:   } else if (eps->isgeneralized) {
931:     if (eps->ishermitian && eps->ispositive) {
932:       DSSetType(eps->ds,DSGHEP);
933:     } else {
934:       DSSetType(eps->ds,DSGNHEP);
935:     }
936:   } else {
937:     if (eps->ishermitian) {
938:       DSSetType(eps->ds,DSHEP);
939:     } else {
940:       DSSetType(eps->ds,DSNHEP);
941:     }
942:   }
943:   DSAllocate(eps->ds,eps->ncv);
944:   EPSSetWorkVecs(eps,2);

946: #if !defined(PETSC_USE_COMPLEX)
947:   if (!eps->ishermitian) { PetscInfo(eps,"Warning: complex eigenvalues are not calculated exactly without --with-scalar-type=complex in PETSc\n"); }
948: #endif
949:   return(0);
950: }

952: PetscErrorCode EPSSetUpSort_CISS(EPS eps)
953: {
955:   SlepcSC        sc;

958:   /* fill sorting criterion context */
959:   eps->sc->comparison    = SlepcCompareSmallestReal;
960:   eps->sc->comparisonctx = NULL;
961:   eps->sc->map           = NULL;
962:   eps->sc->mapobj        = NULL;

964:   /* fill sorting criterion for DS */
965:   DSGetSlepcSC(eps->ds,&sc);
966:   sc->comparison    = SlepcCompareLargestMagnitude;
967:   sc->comparisonctx = NULL;
968:   sc->map           = NULL;
969:   sc->mapobj        = NULL;
970:   return(0);
971: }

973: PetscErrorCode EPSSolve_CISS(EPS eps)
974: {
976:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
977:   Mat            A,B,X,M,pA,pB;
978:   PetscInt       i,j,ld,nmat,L_add=0,nv=0,L_base=ctx->L,inner,nlocal,*inside;
979:   PetscScalar    *Mu,*H0,*H1=NULL,*rr,*temp;
980:   PetscReal      error,max_error,norm;
981:   PetscBool      *fl1;
982:   Vec            si,w[3];
983:   PetscRandom    rand;
984: #if defined(PETSC_USE_COMPLEX)
985:   PetscBool      isellipse;
986: #endif

989:   w[0] = eps->work[0];
990:   w[1] = NULL;
991:   w[2] = eps->work[1];
992:   VecGetLocalSize(w[0],&nlocal);
993:   DSGetLeadingDimension(eps->ds,&ld);
994:   STGetNumMatrices(eps->st,&nmat);
995:   STGetMatrix(eps->st,0,&A);
996:   if (nmat>1) { STGetMatrix(eps->st,1,&B); }
997:   else B = NULL;
998:   SetPathParameter(eps);
999:   CISSVecSetRandom(ctx->V,0,ctx->L);
1000:   BVGetRandomContext(ctx->V,&rand);

1002:   if (ctx->pA) {
1003:     VecScatterVecs(eps,ctx->V,ctx->L);
1004:     SolveLinearSystem(eps,ctx->pA,ctx->pB,ctx->pV,0,ctx->L,PETSC_TRUE);
1005:   } else {
1006:     SolveLinearSystem(eps,A,B,ctx->V,0,ctx->L,PETSC_TRUE);
1007:   }
1008: #if defined(PETSC_USE_COMPLEX)
1009:   PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
1010:   if (isellipse) {
1011:     EstimateNumberEigs(eps,&L_add);
1012:   } else {
1013:     L_add = 0;
1014:   }
1015: #else
1016:   L_add = 0;
1017: #endif
1018:   if (L_add>0) {
1019:     PetscInfo2(eps,"Changing L %D -> %D by Estimate #Eig\n",ctx->L,ctx->L+L_add);
1020:     CISSVecSetRandom(ctx->V,ctx->L,ctx->L+L_add);
1021:     if (ctx->pA) {
1022:       VecScatterVecs(eps,ctx->V,ctx->L+L_add);
1023:       SolveLinearSystem(eps,ctx->pA,ctx->pB,ctx->pV,ctx->L,ctx->L+L_add,PETSC_FALSE);
1024:     } else {
1025:       SolveLinearSystem(eps,A,B,ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
1026:     }
1027:     ctx->L += L_add;
1028:   }
1029:   PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
1030:   for (i=0;i<ctx->refine_blocksize;i++) {
1031:     CalcMu(eps,Mu);
1032:     BlockHankel(eps,Mu,0,H0);
1033:     SVD_H0(eps,H0,&nv);
1034:     if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
1035:     L_add = L_base;
1036:     if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
1037:     PetscInfo2(eps,"Changing L %D -> %D by SVD(H0)\n",ctx->L,ctx->L+L_add);
1038:     CISSVecSetRandom(ctx->V,ctx->L,ctx->L+L_add);
1039:     if (ctx->pA) {
1040:       VecScatterVecs(eps,ctx->V,ctx->L+L_add);
1041:       SolveLinearSystem(eps,ctx->pA,ctx->pB,ctx->pV,ctx->L,ctx->L+L_add,PETSC_FALSE);
1042:     } else {
1043:       SolveLinearSystem(eps,A,B,ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
1044:     }
1045:     ctx->L += L_add;
1046:   }
1047:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
1048:     PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1);
1049:   }

1051:   while (eps->reason == EPS_CONVERGED_ITERATING) {
1052:     eps->its++;
1053:     for (inner=0;inner<=ctx->refine_inner;inner++) {
1054:       if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
1055:         CalcMu(eps,Mu);
1056:         BlockHankel(eps,Mu,0,H0);
1057:         SVD_H0(eps,H0,&nv);
1058:         break;
1059:       } else {
1060:         ConstructS(eps);
1061:         BVSetActiveColumns(ctx->S,0,ctx->L);
1062:         BVCopy(ctx->S,ctx->V);
1063:         PetscLogEventBegin(EPS_CISS_SVD,eps,0,0,0);
1064:         SVD_S(ctx->S,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
1065:         PetscLogEventEnd(EPS_CISS_SVD,eps,0,0,0);
1066:         if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
1067:           if (ctx->pA) {
1068:             VecScatterVecs(eps,ctx->V,ctx->L);
1069:             SolveLinearSystem(eps,ctx->pA,ctx->pB,ctx->pV,0,ctx->L,PETSC_FALSE);
1070:           } else {
1071:             SolveLinearSystem(eps,A,B,ctx->V,0,ctx->L,PETSC_FALSE);
1072:           }
1073:         } else break;
1074:       }
1075:     }
1076:     eps->nconv = 0;
1077:     if (nv == 0) eps->reason = EPS_CONVERGED_TOL;
1078:     else {
1079:       DSSetDimensions(eps->ds,nv,0,0,0);
1080:       DSSetState(eps->ds,DS_STATE_RAW);

1082:       if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
1083:         BlockHankel(eps,Mu,0,H0);
1084:         BlockHankel(eps,Mu,1,H1);
1085:         DSGetArray(eps->ds,DS_MAT_A,&temp);
1086:         for (j=0;j<nv;j++) {
1087:           for (i=0;i<nv;i++) {
1088:             temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
1089:           }
1090:         }
1091:         DSRestoreArray(eps->ds,DS_MAT_A,&temp);
1092:         DSGetArray(eps->ds,DS_MAT_B,&temp);
1093:         for (j=0;j<nv;j++) {
1094:           for (i=0;i<nv;i++) {
1095:             temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
1096:           }
1097:         }
1098:         DSRestoreArray(eps->ds,DS_MAT_B,&temp);
1099:       } else {
1100:         BVSetActiveColumns(ctx->S,0,nv);
1101:         DSGetMat(eps->ds,DS_MAT_A,&pA);
1102:         MatZeroEntries(pA);
1103:         BVMatProject(ctx->S,A,ctx->S,pA);
1104:         DSRestoreMat(eps->ds,DS_MAT_A,&pA);
1105:         if (B) {
1106:           DSGetMat(eps->ds,DS_MAT_B,&pB);
1107:           MatZeroEntries(pB);
1108:           BVMatProject(ctx->S,B,ctx->S,pB);
1109:           DSRestoreMat(eps->ds,DS_MAT_B,&pB);
1110:         }
1111:       }

1113:       DSSolve(eps->ds,eps->eigr,eps->eigi);
1114:       DSSynchronize(eps->ds,eps->eigr,eps->eigi);

1116:       PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr);
1117:       rescale_eig(eps,nv);
1118:       isGhost(eps,ld,nv,fl1);
1119:       RGCheckInside(eps->rg,nv,eps->eigr,eps->eigi,inside);
1120:       for (i=0;i<nv;i++) {
1121:         if (fl1[i] && inside[i]>=0) {
1122:           rr[i] = 1.0;
1123:           eps->nconv++;
1124:         } else rr[i] = 0.0;
1125:       }
1126:       DSSort(eps->ds,eps->eigr,eps->eigi,rr,NULL,&eps->nconv);
1127:       DSSynchronize(eps->ds,eps->eigr,eps->eigi);
1128:       rescale_eig(eps,nv);
1129:       PetscFree3(fl1,inside,rr);
1130:       BVSetActiveColumns(eps->V,0,nv);
1131:       if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
1132:         ConstructS(eps);
1133:         BVSetActiveColumns(ctx->S,0,ctx->L);
1134:         BVCopy(ctx->S,ctx->V);
1135:         BVSetActiveColumns(ctx->S,0,nv);
1136:       }
1137:       BVCopy(ctx->S,eps->V);

1139:       DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
1140:       DSGetMat(eps->ds,DS_MAT_X,&X);
1141:       BVMultInPlace(ctx->S,X,0,eps->nconv);
1142:       if (eps->ishermitian) {
1143:         BVMultInPlace(eps->V,X,0,eps->nconv);
1144:       }
1145:       MatDestroy(&X);
1146:       max_error = 0.0;
1147:       for (i=0;i<eps->nconv;i++) {
1148:         BVGetColumn(ctx->S,i,&si);
1149:         EPSComputeResidualNorm_Private(eps,PETSC_FALSE,eps->eigr[i],eps->eigi[i],si,NULL,w,&error);
1150:         if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {  /* vector is not normalized */
1151:           VecNorm(si,NORM_2,&norm);
1152:           error /= norm;
1153:         }
1154:         (*eps->converged)(eps,eps->eigr[i],eps->eigi[i],error,&error,eps->convergedctx);
1155:         BVRestoreColumn(ctx->S,i,&si);
1156:         max_error = PetscMax(max_error,error);
1157:       }

1159:       if (max_error <= eps->tol) eps->reason = EPS_CONVERGED_TOL;
1160:       else if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
1161:       else {
1162:         if (eps->nconv > ctx->L) {
1163:           MatCreateSeqDense(PETSC_COMM_SELF,eps->nconv,ctx->L,NULL,&M);
1164:           MatDenseGetArray(M,&temp);
1165:           for (i=0;i<ctx->L*eps->nconv;i++) {
1166:             PetscRandomGetValue(rand,&temp[i]);
1167:             temp[i] = PetscRealPart(temp[i]);
1168:           }
1169:           MatDenseRestoreArray(M,&temp);
1170:           BVSetActiveColumns(ctx->S,0,eps->nconv);
1171:           BVMultInPlace(ctx->S,M,0,ctx->L);
1172:           MatDestroy(&M);
1173:           BVSetActiveColumns(ctx->S,0,ctx->L);
1174:           BVCopy(ctx->S,ctx->V);
1175:         }
1176:         if (ctx->pA) {
1177:           VecScatterVecs(eps,ctx->V,ctx->L);
1178:           SolveLinearSystem(eps,ctx->pA,ctx->pB,ctx->pV,0,ctx->L,PETSC_FALSE);
1179:         } else {
1180:           SolveLinearSystem(eps,A,B,ctx->V,0,ctx->L,PETSC_FALSE);
1181:         }
1182:       }
1183:     }
1184:   }
1185:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
1186:     PetscFree(H1);
1187:   }
1188:   PetscFree2(Mu,H0);
1189:   return(0);
1190: }

1192: PetscErrorCode EPSComputeVectors_CISS(EPS eps)
1193: {
1195:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1196:   PetscInt       n;
1197:   Mat            Z,B=NULL;

1200:   if (eps->ishermitian) {
1201:     if (eps->isgeneralized && !eps->ispositive) {
1202:       EPSComputeVectors_Indefinite(eps);
1203:     } else {
1204:       EPSComputeVectors_Hermitian(eps);
1205:     }
1206:     return(0);
1207:   }
1208:   DSGetDimensions(eps->ds,&n,NULL,NULL,NULL,NULL);
1209:   BVSetActiveColumns(eps->V,0,n);

1211:   /* right eigenvectors */
1212:   DSVectors(eps->ds,DS_MAT_X,NULL,NULL);

1214:   /* V = V * Z */
1215:   DSGetMat(eps->ds,DS_MAT_X,&Z);
1216:   BVMultInPlace(eps->V,Z,0,n);
1217:   MatDestroy(&Z);
1218:   BVSetActiveColumns(eps->V,0,eps->nconv);

1220:   /* normalize */
1221:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
1222:     if (eps->isgeneralized && eps->ishermitian && eps->ispositive) {
1223:       STGetMatrix(eps->st,1,&B);
1224:       BVSetMatrix(eps->V,B,PETSC_FALSE);
1225:     }
1226:     BVNormalize(eps->V,NULL);
1227:     if (B) { BVSetMatrix(eps->V,NULL,PETSC_FALSE); }
1228:   }
1229:   return(0);
1230: }

1232: static PetscErrorCode EPSCISSSetSizes_CISS(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
1233: {
1235:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1236:   PetscInt       oN,onpart;

1239:   oN = ctx->N;
1240:   if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
1241:     if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
1242:   } else {
1243:     if (ip<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
1244:     if (ip%2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
1245:     if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
1246:   }
1247:   if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
1248:     ctx->L = 16;
1249:   } else {
1250:     if (bs<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
1251:     ctx->L = bs;
1252:   }
1253:   if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
1254:     ctx->M = ctx->N/4;
1255:   } else {
1256:     if (ms<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
1257:     if (ms>ctx->N) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be less than or equal to the number of integration points");
1258:     ctx->M = ms;
1259:   }
1260:   onpart = ctx->npart;
1261:   if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
1262:     ctx->npart = 1;
1263:   } else {
1264:     if (npart<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The npart argument must be > 0");
1265:     ctx->npart = npart;
1266:   }
1267:   if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
1268:     ctx->L_max = 64;
1269:   } else {
1270:     if (bsmax<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
1271:     ctx->L_max = PetscMax(bsmax,ctx->L);
1272:   }
1273:   if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) { EPSCISSResetSubcomm(eps); }
1274:   ctx->isreal = realmats;
1275:   eps->state = EPS_STATE_INITIAL;
1276:   return(0);
1277: }

1279: /*@
1280:    EPSCISSSetSizes - Sets the values of various size parameters in the CISS solver.

1282:    Logically Collective on eps

1284:    Input Parameters:
1285: +  eps   - the eigenproblem solver context
1286: .  ip    - number of integration points
1287: .  bs    - block size
1288: .  ms    - moment size
1289: .  npart - number of partitions when splitting the communicator
1290: .  bsmax - max block size
1291: -  realmats - A and B are real

1293:    Options Database Keys:
1294: +  -eps_ciss_integration_points - Sets the number of integration points
1295: .  -eps_ciss_blocksize - Sets the block size
1296: .  -eps_ciss_moments - Sets the moment size
1297: .  -eps_ciss_partitions - Sets the number of partitions
1298: .  -eps_ciss_maxblocksize - Sets the maximum block size
1299: -  -eps_ciss_realmats - A and B are real

1301:    Note:
1302:    The default number of partitions is 1. This means the internal KSP object is shared
1303:    among all processes of the EPS communicator. Otherwise, the communicator is split
1304:    into npart communicators, so that npart KSP solves proceed simultaneously.

1306:    Level: advanced

1308: .seealso: EPSCISSGetSizes()
1309: @*/
1310: PetscErrorCode EPSCISSSetSizes(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
1311: {

1322:   PetscTryMethod(eps,"EPSCISSSetSizes_C",(EPS,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(eps,ip,bs,ms,npart,bsmax,realmats));
1323:   return(0);
1324: }

1326: static PetscErrorCode EPSCISSGetSizes_CISS(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
1327: {
1328:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1331:   if (ip) *ip = ctx->N;
1332:   if (bs) *bs = ctx->L;
1333:   if (ms) *ms = ctx->M;
1334:   if (npart) *npart = ctx->npart;
1335:   if (bsmax) *bsmax = ctx->L_max;
1336:   if (realmats) *realmats = ctx->isreal;
1337:   return(0);
1338: }

1340: /*@
1341:    EPSCISSGetSizes - Gets the values of various size parameters in the CISS solver.

1343:    Not Collective

1345:    Input Parameter:
1346: .  eps - the eigenproblem solver context

1348:    Output Parameters:
1349: +  ip    - number of integration points
1350: .  bs    - block size
1351: .  ms    - moment size
1352: .  npart - number of partitions when splitting the communicator
1353: .  bsmax - max block size
1354: -  realmats - A and B are real

1356:    Level: advanced

1358: .seealso: EPSCISSSetSizes()
1359: @*/
1360: PetscErrorCode EPSCISSGetSizes(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
1361: {

1366:   PetscUseMethod(eps,"EPSCISSGetSizes_C",(EPS,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(eps,ip,bs,ms,npart,bsmax,realmats));
1367:   return(0);
1368: }

1370: static PetscErrorCode EPSCISSSetThreshold_CISS(EPS eps,PetscReal delta,PetscReal spur)
1371: {
1372:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1375:   if (delta == PETSC_DEFAULT) {
1376:     ctx->delta = 1e-12;
1377:   } else {
1378:     if (delta<=0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
1379:     ctx->delta = delta;
1380:   }
1381:   if (spur == PETSC_DEFAULT) {
1382:     ctx->spurious_threshold = 1e-4;
1383:   } else {
1384:     if (spur<=0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
1385:     ctx->spurious_threshold = spur;
1386:   }
1387:   return(0);
1388: }

1390: /*@
1391:    EPSCISSSetThreshold - Sets the values of various threshold parameters in
1392:    the CISS solver.

1394:    Logically Collective on eps

1396:    Input Parameters:
1397: +  eps   - the eigenproblem solver context
1398: .  delta - threshold for numerical rank
1399: -  spur  - spurious threshold (to discard spurious eigenpairs)

1401:    Options Database Keys:
1402: +  -eps_ciss_delta - Sets the delta
1403: -  -eps_ciss_spurious_threshold - Sets the spurious threshold

1405:    Level: advanced

1407: .seealso: EPSCISSGetThreshold()
1408: @*/
1409: PetscErrorCode EPSCISSSetThreshold(EPS eps,PetscReal delta,PetscReal spur)
1410: {

1417:   PetscTryMethod(eps,"EPSCISSSetThreshold_C",(EPS,PetscReal,PetscReal),(eps,delta,spur));
1418:   return(0);
1419: }

1421: static PetscErrorCode EPSCISSGetThreshold_CISS(EPS eps,PetscReal *delta,PetscReal *spur)
1422: {
1423:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1426:   if (delta) *delta = ctx->delta;
1427:   if (spur)  *spur = ctx->spurious_threshold;
1428:   return(0);
1429: }

1431: /*@
1432:    EPSCISSGetThreshold - Gets the values of various threshold parameters
1433:    in the CISS solver.

1435:    Not Collective

1437:    Input Parameter:
1438: .  eps - the eigenproblem solver context

1440:    Output Parameters:
1441: +  delta - threshold for numerical rank
1442: -  spur  - spurious threshold (to discard spurious eigenpairs)

1444:    Level: advanced

1446: .seealso: EPSCISSSetThreshold()
1447: @*/
1448: PetscErrorCode EPSCISSGetThreshold(EPS eps,PetscReal *delta,PetscReal *spur)
1449: {

1454:   PetscUseMethod(eps,"EPSCISSGetThreshold_C",(EPS,PetscReal*,PetscReal*),(eps,delta,spur));
1455:   return(0);
1456: }

1458: static PetscErrorCode EPSCISSSetRefinement_CISS(EPS eps,PetscInt inner,PetscInt blsize)
1459: {
1460:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1463:   if (inner == PETSC_DEFAULT) {
1464:     ctx->refine_inner = 0;
1465:   } else {
1466:     if (inner<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
1467:     ctx->refine_inner = inner;
1468:   }
1469:   if (blsize == PETSC_DEFAULT) {
1470:     ctx->refine_blocksize = 0;
1471:   } else {
1472:     if (blsize<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
1473:     ctx->refine_blocksize = blsize;
1474:   }
1475:   return(0);
1476: }

1478: /*@
1479:    EPSCISSSetRefinement - Sets the values of various refinement parameters
1480:    in the CISS solver.

1482:    Logically Collective on eps

1484:    Input Parameters:
1485: +  eps    - the eigenproblem solver context
1486: .  inner  - number of iterative refinement iterations (inner loop)
1487: -  blsize - number of iterative refinement iterations (blocksize loop)

1489:    Options Database Keys:
1490: +  -eps_ciss_refine_inner - Sets number of inner iterations
1491: -  -eps_ciss_refine_blocksize - Sets number of blocksize iterations

1493:    Level: advanced

1495: .seealso: EPSCISSGetRefinement()
1496: @*/
1497: PetscErrorCode EPSCISSSetRefinement(EPS eps,PetscInt inner,PetscInt blsize)
1498: {

1505:   PetscTryMethod(eps,"EPSCISSSetRefinement_C",(EPS,PetscInt,PetscInt),(eps,inner,blsize));
1506:   return(0);
1507: }

1509: static PetscErrorCode EPSCISSGetRefinement_CISS(EPS eps,PetscInt *inner,PetscInt *blsize)
1510: {
1511:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1514:   if (inner)  *inner = ctx->refine_inner;
1515:   if (blsize) *blsize = ctx->refine_blocksize;
1516:   return(0);
1517: }

1519: /*@
1520:    EPSCISSGetRefinement - Gets the values of various refinement parameters
1521:    in the CISS solver.

1523:    Not Collective

1525:    Input Parameter:
1526: .  eps - the eigenproblem solver context

1528:    Output Parameters:
1529: +  inner  - number of iterative refinement iterations (inner loop)
1530: -  blsize - number of iterative refinement iterations (blocksize loop)

1532:    Level: advanced

1534: .seealso: EPSCISSSetRefinement()
1535: @*/
1536: PetscErrorCode EPSCISSGetRefinement(EPS eps, PetscInt *inner, PetscInt *blsize)
1537: {

1542:   PetscUseMethod(eps,"EPSCISSGetRefinement_C",(EPS,PetscInt*,PetscInt*),(eps,inner,blsize));
1543:   return(0);
1544: }

1546: static PetscErrorCode EPSCISSSetUseST_CISS(EPS eps,PetscBool usest)
1547: {
1548:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1551:   ctx->usest     = usest;
1552:   ctx->usest_set = PETSC_TRUE;
1553:   eps->state     = EPS_STATE_INITIAL;
1554:   return(0);
1555: }

1557: /*@
1558:    EPSCISSSetUseST - Sets a flag indicating that the CISS solver will
1559:    use the ST object for the linear solves.

1561:    Logically Collective on eps

1563:    Input Parameters:
1564: +  eps    - the eigenproblem solver context
1565: -  usest  - boolean flag to use the ST object or not

1567:    Options Database Keys:
1568: .  -eps_ciss_usest <bool> - whether the ST object will be used or not

1570:    Level: advanced

1572: .seealso: EPSCISSGetUseST()
1573: @*/
1574: PetscErrorCode EPSCISSSetUseST(EPS eps,PetscBool usest)
1575: {

1581:   PetscTryMethod(eps,"EPSCISSSetUseST_C",(EPS,PetscBool),(eps,usest));
1582:   return(0);
1583: }

1585: static PetscErrorCode EPSCISSGetUseST_CISS(EPS eps,PetscBool *usest)
1586: {
1587:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1590:   *usest = ctx->usest;
1591:   return(0);
1592: }

1594: /*@
1595:    EPSCISSGetUseST - Gets the flag for using the ST object
1596:    in the CISS solver.

1598:    Not Collective

1600:    Input Parameter:
1601: .  eps - the eigenproblem solver context

1603:    Output Parameters:
1604: .  usest - boolean flag indicating if the ST object is being used

1606:    Level: advanced

1608: .seealso: EPSCISSSetUseST()
1609: @*/
1610: PetscErrorCode EPSCISSGetUseST(EPS eps,PetscBool *usest)
1611: {

1617:   PetscUseMethod(eps,"EPSCISSGetUseST_C",(EPS,PetscBool*),(eps,usest));
1618:   return(0);
1619: }

1621: static PetscErrorCode EPSCISSSetQuadRule_CISS(EPS eps,EPSCISSQuadRule quad)
1622: {
1623:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1626:   ctx->quad = quad;
1627:   return(0);
1628: }

1630: /*@
1631:    EPSCISSSetQuadRule - Sets the quadrature rule used in the CISS solver.

1633:    Logically Collective on eps

1635:    Input Parameters:
1636: +  eps  - the eigenproblem solver context
1637: -  quad - the quadrature rule

1639:    Options Database Key:
1640: .  -eps_ciss_quadrule - Sets the quadrature rule (either 'trapezoidal' or
1641:                            'chebyshev')

1643:    Notes:
1644:    By default, the trapezoidal rule is used (EPS_CISS_QUADRULE_TRAPEZOIDAL).

1646:    If the 'chebyshev' option is specified (EPS_CISS_QUADRULE_CHEBYSHEV), then
1647:    Chebyshev points are used as quadrature points.

1649:    Level: advanced

1651: .seealso: EPSCISSGetQuadRule(), EPSCISSQuadRule
1652: @*/
1653: PetscErrorCode EPSCISSSetQuadRule(EPS eps,EPSCISSQuadRule quad)
1654: {

1660:   PetscTryMethod(eps,"EPSCISSSetQuadRule_C",(EPS,EPSCISSQuadRule),(eps,quad));
1661:   return(0);
1662: }

1664: static PetscErrorCode EPSCISSGetQuadRule_CISS(EPS eps,EPSCISSQuadRule *quad)
1665: {
1666:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1669:   *quad = ctx->quad;
1670:   return(0);
1671: }

1673: /*@
1674:    EPSCISSGetQuadRule - Gets the quadrature rule used in the CISS solver.

1676:    Not Collective

1678:    Input Parameter:
1679: .  eps - the eigenproblem solver context

1681:    Output Parameters:
1682: .  quad - quadrature rule

1684:    Level: advanced

1686: .seealso: EPSCISSSetQuadRule() EPSCISSQuadRule
1687: @*/
1688: PetscErrorCode EPSCISSGetQuadRule(EPS eps,EPSCISSQuadRule *quad)
1689: {

1695:   PetscUseMethod(eps,"EPSCISSGetQuadRule_C",(EPS,EPSCISSQuadRule*),(eps,quad));
1696:   return(0);
1697: }

1699: static PetscErrorCode EPSCISSSetExtraction_CISS(EPS eps,EPSCISSExtraction extraction)
1700: {
1701:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1704:   ctx->extraction = extraction;
1705:   return(0);
1706: }

1708: /*@
1709:    EPSCISSSetExtraction - Sets the extraction technique used in the CISS solver.

1711:    Logically Collective on eps

1713:    Input Parameters:
1714: +  eps        - the eigenproblem solver context
1715: -  extraction - the extraction technique

1717:    Options Database Key:
1718: .  -eps_ciss_extraction - Sets the extraction technique (either 'ritz' or
1719:                            'hankel')

1721:    Notes:
1722:    By default, the Rayleigh-Ritz extraction is used (EPS_CISS_EXTRACTION_RITZ).

1724:    If the 'hankel' option is specified (EPS_CISS_EXTRACTION_HANKEL), then
1725:    the Block Hankel method is used for extracting eigenpairs.

1727:    Level: advanced

1729: .seealso: EPSCISSGetExtraction(), EPSCISSExtraction
1730: @*/
1731: PetscErrorCode EPSCISSSetExtraction(EPS eps,EPSCISSExtraction extraction)
1732: {

1738:   PetscTryMethod(eps,"EPSCISSSetExtraction_C",(EPS,EPSCISSExtraction),(eps,extraction));
1739:   return(0);
1740: }

1742: static PetscErrorCode EPSCISSGetExtraction_CISS(EPS eps,EPSCISSExtraction *extraction)
1743: {
1744:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1747:   *extraction = ctx->extraction;
1748:   return(0);
1749: }

1751: /*@
1752:    EPSCISSGetExtraction - Gets the extraction technique used in the CISS solver.

1754:    Not Collective

1756:    Input Parameter:
1757: .  eps - the eigenproblem solver context

1759:    Output Parameters:
1760: .  extraction - extraction technique

1762:    Level: advanced

1764: .seealso: EPSCISSSetExtraction() EPSCISSExtraction
1765: @*/
1766: PetscErrorCode EPSCISSGetExtraction(EPS eps,EPSCISSExtraction *extraction)
1767: {

1773:   PetscUseMethod(eps,"EPSCISSGetExtraction_C",(EPS,EPSCISSExtraction*),(eps,extraction));
1774:   return(0);
1775: }

1777: static PetscErrorCode EPSCISSGetKSPs_CISS(EPS eps,PetscInt *nsolve,KSP **ksp)
1778: {
1780:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1781:   PetscInt       i;
1782:   PC             pc;

1785:   if (!ctx->ksp) {
1786:     if (!ctx->subcomm) {  /* initialize subcomm first */
1787:       EPSCISSSetUseConj(eps,&ctx->useconj);
1788:       EPSCISSSetUpSubComm(eps,&ctx->num_solve_point);
1789:     }
1790:     PetscMalloc1(ctx->num_solve_point,&ctx->ksp);
1791:     for (i=0;i<ctx->num_solve_point;i++) {
1792:       KSPCreate(PetscSubcommChild(ctx->subcomm),&ctx->ksp[i]);
1793:       PetscObjectIncrementTabLevel((PetscObject)ctx->ksp[i],(PetscObject)eps,1);
1794:       KSPSetOptionsPrefix(ctx->ksp[i],((PetscObject)eps)->prefix);
1795:       KSPAppendOptionsPrefix(ctx->ksp[i],"eps_ciss_");
1796:       PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->ksp[i]);
1797:       PetscObjectSetOptions((PetscObject)ctx->ksp[i],((PetscObject)eps)->options);
1798:       KSPSetErrorIfNotConverged(ctx->ksp[i],PETSC_TRUE);
1799:       KSPSetTolerances(ctx->ksp[i],SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
1800:       KSPGetPC(ctx->ksp[i],&pc);
1801:       KSPSetType(ctx->ksp[i],KSPPREONLY);
1802:       PCSetType(pc,PCLU);
1803:     }
1804:   }
1805:   if (nsolve) *nsolve = ctx->num_solve_point;
1806:   if (ksp)    *ksp    = ctx->ksp;
1807:   return(0);
1808: }

1810: /*@C
1811:    EPSCISSGetKSPs - Retrieve the array of linear solver objects associated with
1812:    the CISS solver.

1814:    Not Collective

1816:    Input Parameter:
1817: .  eps - the eigenproblem solver solver

1819:    Output Parameters:
1820: +  nsolve - number of solver objects
1821: -  ksp - array of linear solver object

1823:    Notes:
1824:    The number of KSP solvers is equal to the number of integration points divided by
1825:    the number of partitions. This value is halved in the case of real matrices with
1826:    a region centered at the real axis.

1828:    Level: advanced

1830: .seealso: EPSCISSSetSizes()
1831: @*/
1832: PetscErrorCode EPSCISSGetKSPs(EPS eps,PetscInt *nsolve,KSP **ksp)
1833: {

1838:   PetscUseMethod(eps,"EPSCISSGetKSPs_C",(EPS,PetscInt*,KSP**),(eps,nsolve,ksp));
1839:   return(0);
1840: }

1842: PetscErrorCode EPSReset_CISS(EPS eps)
1843: {
1845:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1846:   PetscInt       i;

1849:   BVDestroy(&ctx->S);
1850:   BVDestroy(&ctx->V);
1851:   BVDestroy(&ctx->Y);
1852:   if (!ctx->usest) {
1853:     for (i=0;i<ctx->num_solve_point;i++) {
1854:       KSPReset(ctx->ksp[i]);
1855:     }
1856:   }
1857:   VecScatterDestroy(&ctx->scatterin);
1858:   VecDestroy(&ctx->xsub);
1859:   VecDestroy(&ctx->xdup);
1860:   if (ctx->pA) {
1861:     MatDestroy(&ctx->pA);
1862:     MatDestroy(&ctx->pB);
1863:     BVDestroy(&ctx->pV);
1864:   }
1865:   return(0);
1866: }

1868: PetscErrorCode EPSSetFromOptions_CISS(PetscOptionItems *PetscOptionsObject,EPS eps)
1869: {
1870:   PetscErrorCode    ierr;
1871:   PetscReal         r3,r4;
1872:   PetscInt          i,i1,i2,i3,i4,i5,i6,i7;
1873:   PetscBool         b1,b2,flg;
1874:   EPS_CISS          *ctx = (EPS_CISS*)eps->data;
1875:   EPSCISSQuadRule   quad;
1876:   EPSCISSExtraction extraction;

1879:   PetscOptionsHead(PetscOptionsObject,"EPS CISS Options");

1881:     EPSCISSGetSizes(eps,&i1,&i2,&i3,&i4,&i5,&b1);
1882:     PetscOptionsInt("-eps_ciss_integration_points","Number of integration points","EPSCISSSetSizes",i1,&i1,NULL);
1883:     PetscOptionsInt("-eps_ciss_blocksize","Block size","EPSCISSSetSizes",i2,&i2,NULL);
1884:     PetscOptionsInt("-eps_ciss_moments","Moment size","EPSCISSSetSizes",i3,&i3,NULL);
1885:     PetscOptionsInt("-eps_ciss_partitions","Number of partitions","EPSCISSSetSizes",i4,&i4,NULL);
1886:     PetscOptionsInt("-eps_ciss_maxblocksize","Maximum block size","EPSCISSSetSizes",i5,&i5,NULL);
1887:     PetscOptionsBool("-eps_ciss_realmats","True if A and B are real","EPSCISSSetSizes",b1,&b1,NULL);
1888:     EPSCISSSetSizes(eps,i1,i2,i3,i4,i5,b1);

1890:     EPSCISSGetThreshold(eps,&r3,&r4);
1891:     PetscOptionsReal("-eps_ciss_delta","Threshold for numerical rank","EPSCISSSetThreshold",r3,&r3,NULL);
1892:     PetscOptionsReal("-eps_ciss_spurious_threshold","Threshold for the spurious eigenpairs","EPSCISSSetThreshold",r4,&r4,NULL);
1893:     EPSCISSSetThreshold(eps,r3,r4);

1895:     EPSCISSGetRefinement(eps,&i6,&i7);
1896:     PetscOptionsInt("-eps_ciss_refine_inner","Number of inner iterative refinement iterations","EPSCISSSetRefinement",i6,&i6,NULL);
1897:     PetscOptionsInt("-eps_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","EPSCISSSetRefinement",i7,&i7,NULL);
1898:     EPSCISSSetRefinement(eps,i6,i7);

1900:     EPSCISSGetUseST(eps,&b2);
1901:     PetscOptionsBool("-eps_ciss_usest","Use ST for linear solves","EPSCISSSetUseST",b2,&b2,&flg);
1902:     if (flg) { EPSCISSSetUseST(eps,b2); }

1904:     PetscOptionsEnum("-eps_ciss_quadrule","Quadrature rule","EPSCISSSetQuadRule",EPSCISSQuadRules,(PetscEnum)ctx->quad,(PetscEnum*)&quad,&flg);
1905:     if (flg) { EPSCISSSetQuadRule(eps,quad); }

1907:     PetscOptionsEnum("-eps_ciss_extraction","Extraction technique","EPSCISSSetExtraction",EPSCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg);
1908:     if (flg) { EPSCISSSetExtraction(eps,extraction); }

1910:   PetscOptionsTail();

1912:   if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
1913:   RGSetFromOptions(eps->rg); /* this is necessary here to set useconj */
1914:   if (!ctx->ksp) { EPSCISSGetKSPs(eps,&ctx->num_solve_point,&ctx->ksp); }
1915:   for (i=0;i<ctx->num_solve_point;i++) {
1916:     KSPSetFromOptions(ctx->ksp[i]);
1917:   }
1918:   PetscSubcommSetFromOptions(ctx->subcomm);
1919:   return(0);
1920: }

1922: PetscErrorCode EPSDestroy_CISS(EPS eps)
1923: {
1925:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

1928:   EPSCISSResetSubcomm(eps);
1929:   PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
1930:   PetscFree(eps->data);
1931:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",NULL);
1932:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",NULL);
1933:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",NULL);
1934:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",NULL);
1935:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",NULL);
1936:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",NULL);
1937:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetUseST_C",NULL);
1938:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetUseST_C",NULL);
1939:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetQuadRule_C",NULL);
1940:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetQuadRule_C",NULL);
1941:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetExtraction_C",NULL);
1942:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetExtraction_C",NULL);
1943:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetKSPs_C",NULL);
1944:   return(0);
1945: }

1947: PetscErrorCode EPSView_CISS(EPS eps,PetscViewer viewer)
1948: {
1950:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1951:   PetscBool      isascii;

1954:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1955:   if (isascii) {
1956:     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);
1957:     if (ctx->isreal) {
1958:       PetscViewerASCIIPrintf(viewer,"  exploiting symmetry of integration points\n");
1959:     }
1960:     PetscViewerASCIIPrintf(viewer,"  threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold);
1961:     PetscViewerASCIIPrintf(viewer,"  iterative refinement { inner: %D, blocksize: %D }\n",ctx->refine_inner, ctx->refine_blocksize);
1962:     PetscViewerASCIIPrintf(viewer,"  extraction: %s\n",EPSCISSExtractions[ctx->extraction]);
1963:     PetscViewerASCIIPrintf(viewer,"  quadrature rule: %s\n",EPSCISSQuadRules[ctx->quad]);
1964:     if (ctx->usest) {
1965:       PetscViewerASCIIPrintf(viewer,"  using ST for linear solves\n");
1966:     } else {
1967:       if (!ctx->ksp) { EPSCISSGetKSPs(eps,&ctx->num_solve_point,&ctx->ksp); }
1968:       PetscViewerASCIIPushTab(viewer);
1969:       KSPView(ctx->ksp[0],viewer);
1970:       PetscViewerASCIIPopTab(viewer);
1971:     }
1972:   }
1973:   return(0);
1974: }

1976: PetscErrorCode EPSSetDefaultST_CISS(EPS eps)
1977: {
1979:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1980:   PetscBool      usest = ctx->usest;

1983:   if (!((PetscObject)eps->st)->type_name) {
1984:     if (!ctx->usest_set) usest = (ctx->npart>1)? PETSC_FALSE: PETSC_TRUE;
1985:     if (usest) {
1986:       STSetType(eps->st,STSINVERT);
1987:     } else {
1988:       /* we are not going to use ST, so avoid factorizing the matrix */
1989:       STSetType(eps->st,STSHIFT);
1990:     }
1991:   }
1992:   return(0);
1993: }

1995: SLEPC_EXTERN PetscErrorCode EPSCreate_CISS(EPS eps)
1996: {
1998:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

2001:   PetscNewLog(eps,&ctx);
2002:   eps->data = ctx;

2004:   eps->useds = PETSC_TRUE;
2005:   eps->categ = EPS_CATEGORY_CONTOUR;

2007:   eps->ops->solve          = EPSSolve_CISS;
2008:   eps->ops->setup          = EPSSetUp_CISS;
2009:   eps->ops->setupsort      = EPSSetUpSort_CISS;
2010:   eps->ops->setfromoptions = EPSSetFromOptions_CISS;
2011:   eps->ops->destroy        = EPSDestroy_CISS;
2012:   eps->ops->reset          = EPSReset_CISS;
2013:   eps->ops->view           = EPSView_CISS;
2014:   eps->ops->computevectors = EPSComputeVectors_CISS;
2015:   eps->ops->setdefaultst   = EPSSetDefaultST_CISS;

2017:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",EPSCISSSetSizes_CISS);
2018:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",EPSCISSGetSizes_CISS);
2019:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",EPSCISSSetThreshold_CISS);
2020:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",EPSCISSGetThreshold_CISS);
2021:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",EPSCISSSetRefinement_CISS);
2022:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",EPSCISSGetRefinement_CISS);
2023:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetUseST_C",EPSCISSSetUseST_CISS);
2024:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetUseST_C",EPSCISSGetUseST_CISS);
2025:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetQuadRule_C",EPSCISSSetQuadRule_CISS);
2026:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetQuadRule_C",EPSCISSGetQuadRule_CISS);
2027:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetExtraction_C",EPSCISSSetExtraction_CISS);
2028:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetExtraction_C",EPSCISSGetExtraction_CISS);
2029:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetKSPs_C",EPSCISSGetKSPs_CISS);

2031:   /* log events */
2032:   PetscLogEventRegister("EPSCISS_SVD",EPS_CLASSID,&EPS_CISS_SVD);

2034:   /* set default values of parameters */
2035:   ctx->N                  = 32;
2036:   ctx->L                  = 16;
2037:   ctx->M                  = ctx->N/4;
2038:   ctx->delta              = 1e-12;
2039:   ctx->L_max              = 64;
2040:   ctx->spurious_threshold = 1e-4;
2041:   ctx->usest              = PETSC_TRUE;
2042:   ctx->usest_set          = PETSC_FALSE;
2043:   ctx->isreal             = PETSC_FALSE;
2044:   ctx->refine_inner       = 0;
2045:   ctx->refine_blocksize   = 0;
2046:   ctx->npart              = 1;
2047:   ctx->quad               = (EPSCISSQuadRule)0;
2048:   ctx->extraction         = EPS_CISS_EXTRACTION_RITZ;
2049:   return(0);
2050: }