Actual source code: ciss.c

slepc-3.16.0 2021-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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 <slepc/private/slepccontour.h>
 35: #include <slepcblaslapack.h>

 37: typedef struct {
 38:   /* user parameters */
 39:   PetscInt          N;          /* number of integration points (32) */
 40:   PetscInt          L;          /* block size (16) */
 41:   PetscInt          M;          /* moment degree (N/4 = 4) */
 42:   PetscReal         delta;      /* threshold of singular value (1e-12) */
 43:   PetscInt          L_max;      /* maximum number of columns of the source matrix V */
 44:   PetscReal         spurious_threshold; /* discard spurious eigenpairs */
 45:   PetscBool         isreal;     /* A and B are real */
 46:   PetscInt          npart;      /* number of partitions */
 47:   PetscInt          refine_inner;
 48:   PetscInt          refine_blocksize;
 49:   EPSCISSQuadRule   quad;
 50:   EPSCISSExtraction extraction;
 51:   PetscBool         usest;
 52:   /* private data */
 53:   SlepcContourData  contour;
 54:   PetscReal         *sigma;     /* threshold for numerical rank */
 55:   PetscScalar       *weight;
 56:   PetscScalar       *omega;
 57:   PetscScalar       *pp;
 58:   BV                V;
 59:   BV                S;
 60:   BV                pV;
 61:   BV                Y;
 62:   PetscBool         useconj;
 63:   PetscBool         usest_set;  /* whether the user set the usest flag or not */
 64:   PetscObjectId     rgid;
 65:   PetscObjectState  rgstate;
 66: } EPS_CISS;

 68: static PetscErrorCode EPSCISSSolveSystem(EPS eps,Mat A,Mat B,BV V,PetscInt L_start,PetscInt L_end,PetscBool initksp)
 69: {
 70:   PetscErrorCode   ierr;
 71:   EPS_CISS         *ctx = (EPS_CISS*)eps->data;
 72:   SlepcContourData contour = ctx->contour;
 73:   PetscInt         i,p_id;
 74:   Mat              Fz,kspMat,MV,BMV=NULL,MC;
 75:   KSP              ksp;
 76:   const char       *prefix;

 79:   if (!contour || !contour->ksp) { EPSCISSGetKSPs(eps,NULL,NULL); }
 80:   if (ctx->usest) {
 81:     MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Fz);
 82:   }
 83:   BVSetActiveColumns(V,L_start,L_end);
 84:   BVGetMat(V,&MV);
 85:   if (B) {
 86:     MatProductCreate(B,MV,NULL,&BMV);
 87:     MatProductSetType(BMV,MATPRODUCT_AB);
 88:     MatProductSetFromOptions(BMV);
 89:     MatProductSymbolic(BMV);
 90:   }
 91:   for (i=0;i<contour->npoints;i++) {
 92:     p_id = i*contour->subcomm->n + contour->subcomm->color;
 93:     if (!ctx->usest && initksp) {
 94:       MatDuplicate(A,MAT_COPY_VALUES,&kspMat);
 95:       if (B) {
 96:         MatAXPY(kspMat,-ctx->omega[p_id],B,UNKNOWN_NONZERO_PATTERN);
 97:       } else {
 98:         MatShift(kspMat,-ctx->omega[p_id]);
 99:       }
100:       KSPSetOperators(contour->ksp[i],kspMat,kspMat);
101:       /* set Mat prefix to be the same as KSP to enable setting command-line options (e.g. MUMPS) */
102:       KSPGetOptionsPrefix(contour->ksp[i],&prefix);
103:       MatSetOptionsPrefix(kspMat,prefix);
104:       MatDestroy(&kspMat);
105:     } else if (ctx->usest) {
106:       STSetShift(eps->st,ctx->omega[p_id]);
107:       STGetKSP(eps->st,&ksp);
108:     }
109:     BVSetActiveColumns(ctx->Y,i*ctx->L_max+L_start,i*ctx->L_max+L_end);
110:     BVGetMat(ctx->Y,&MC);
111:     if (B) {
112:       MatProductNumeric(BMV);
113:       if (ctx->usest) {
114:         KSPMatSolve(ksp,BMV,MC);
115:       } else {
116:         KSPMatSolve(contour->ksp[i],BMV,MC);
117:       }
118:     } else {
119:       if (ctx->usest) {
120:         KSPMatSolve(ksp,MV,MC);
121:       } else {
122:         KSPMatSolve(contour->ksp[i],MV,MC);
123:       }
124:     }
125:     if (ctx->usest && i<contour->npoints-1) { KSPReset(ksp); }
126:     BVRestoreMat(ctx->Y,&MC);
127:   }
128:   MatDestroy(&BMV);
129:   BVRestoreMat(V,&MV);
130:   if (ctx->usest) { MatDestroy(&Fz); }
131:   return(0);
132: }

134: static PetscErrorCode rescale_eig(EPS eps,PetscInt nv)
135: {
137:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
138:   PetscInt       i;
139:   PetscScalar    center;
140:   PetscReal      radius,a,b,c,d,rgscale;
141: #if defined(PETSC_USE_COMPLEX)
142:   PetscReal      start_ang,end_ang,vscale,theta;
143: #endif
144:   PetscBool      isring,isellipse,isinterval;

147:   PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
148:   PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring);
149:   PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval);
150:   RGGetScale(eps->rg,&rgscale);
151:   if (isinterval) {
152:     RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d);
153:     if (c==d) {
154:       for (i=0;i<nv;i++) {
155: #if defined(PETSC_USE_COMPLEX)
156:         eps->eigr[i] = PetscRealPart(eps->eigr[i]);
157: #else
158:         eps->eigi[i] = 0;
159: #endif
160:       }
161:     }
162:   }
163:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
164:     if (isellipse) {
165:       RGEllipseGetParameters(eps->rg,&center,&radius,NULL);
166:       for (i=0;i<nv;i++) eps->eigr[i] = rgscale*(center + radius*eps->eigr[i]);
167:     } else if (isinterval) {
168:       RGIntervalGetEndpoints(eps->rg,&a,&b,&c,&d);
169:       if (ctx->quad == EPS_CISS_QUADRULE_CHEBYSHEV) {
170:         for (i=0;i<nv;i++) {
171:           if (c==d) eps->eigr[i] = ((b-a)*(eps->eigr[i]+1.0)/2.0+a)*rgscale;
172:           if (a==b) {
173: #if defined(PETSC_USE_COMPLEX)
174:             eps->eigr[i] = ((d-c)*(eps->eigr[i]+1.0)/2.0+c)*rgscale*PETSC_i;
175: #else
176:             SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Integration points on a vertical line require complex arithmetic");
177: #endif
178:           }
179:         }
180:       } else {
181:         center = (b+a)/2.0+(d+c)/2.0*PETSC_PI;
182:         radius = PetscSqrtReal(PetscPowRealInt((b-a)/2.0,2)+PetscPowRealInt((d-c)/2.0,2));
183:         for (i=0;i<nv;i++) eps->eigr[i] = center + radius*eps->eigr[i];
184:       }
185:     } else if (isring) {  /* only supported in complex scalars */
186: #if defined(PETSC_USE_COMPLEX)
187:       RGRingGetParameters(eps->rg,&center,&radius,&vscale,&start_ang,&end_ang,NULL);
188:       if (ctx->quad == EPS_CISS_QUADRULE_CHEBYSHEV) {
189:         for (i=0;i<nv;i++) {
190:           theta = (start_ang*2.0+(end_ang-start_ang)*(PetscRealPart(eps->eigr[i])+1.0))*PETSC_PI;
191:           eps->eigr[i] = rgscale*center + (rgscale*radius+PetscImaginaryPart(eps->eigr[i]))*PetscCMPLX(PetscCosReal(theta),vscale*PetscSinReal(theta));
192:         }
193:       } else {
194:         for (i=0;i<nv;i++) eps->eigr[i] = rgscale*(center + radius*eps->eigr[i]);
195:       }
196: #endif
197:     }
198:   }
199:   return(0);
200: }

202: PetscErrorCode EPSSetUp_CISS(EPS eps)
203: {
204:   PetscErrorCode   ierr;
205:   EPS_CISS         *ctx = (EPS_CISS*)eps->data;
206:   SlepcContourData contour;
207:   PetscBool        istrivial,isring,isellipse,isinterval,flg;
208:   PetscReal        c,d;
209:   PetscRandom      rand;
210:   PetscObjectId    id;
211:   PetscObjectState state;
212:   Mat              A[2];
213:   Vec              v0;

216:   if (eps->ncv==PETSC_DEFAULT) {
217:     eps->ncv = ctx->L_max*ctx->M;
218:     if (eps->ncv>eps->n) {
219:       eps->ncv = eps->n;
220:       ctx->L_max = eps->ncv/ctx->M;
221:       if (!ctx->L_max) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cannot adjust solver parameters, try setting a smaller value of M (moment size)");
222:     }
223:   } else {
224:     EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
225:     ctx->L_max = eps->ncv/ctx->M;
226:     if (!ctx->L_max) {
227:       ctx->L_max = 1;
228:       eps->ncv = ctx->L_max*ctx->M;
229:     }
230:   }
231:   ctx->L = PetscMin(ctx->L,ctx->L_max);
232:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = 5;
233:   if (eps->mpd==PETSC_DEFAULT) eps->mpd = eps->ncv;
234:   if (!eps->which) eps->which = EPS_ALL;
235:   if (eps->which!=EPS_ALL) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only computing all eigenvalues");
236:   EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_EXTRACTION | EPS_FEATURE_STOPPING | EPS_FEATURE_TWOSIDED);

238:   /* check region */
239:   RGIsTrivial(eps->rg,&istrivial);
240:   if (istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
241:   RGGetComplement(eps->rg,&flg);
242:   if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
243:   PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
244:   PetscObjectTypeCompare((PetscObject)eps->rg,RGRING,&isring);
245:   PetscObjectTypeCompare((PetscObject)eps->rg,RGINTERVAL,&isinterval);
246:   if (!isellipse && !isring && !isinterval) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Currently only implemented for interval, elliptic or ring regions");

248:   /* if the region has changed, then reset contour data */
249:   PetscObjectGetId((PetscObject)eps->rg,&id);
250:   PetscObjectStateGet((PetscObject)eps->rg,&state);
251:   if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
252:     SlepcContourDataDestroy(&ctx->contour);
253:     PetscInfo(eps,"Resetting the contour data structure due to a change of region\n");
254:     ctx->rgid = id; ctx->rgstate = state;
255:   }

257: #if !defined(PETSC_USE_COMPLEX)
258:   if (isring) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Ring region only supported for complex scalars");
259: #endif
260:   if (isinterval) {
261:     RGIntervalGetEndpoints(eps->rg,NULL,NULL,&c,&d);
262: #if !defined(PETSC_USE_COMPLEX)
263:     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");
264: #endif
265:     if (!ctx->quad && c==d) ctx->quad = EPS_CISS_QUADRULE_CHEBYSHEV;
266:   }
267:   if (!ctx->quad) ctx->quad = EPS_CISS_QUADRULE_TRAPEZOIDAL;

269:   /* create contour data structure */
270:   if (!ctx->contour) {
271:     RGCanUseConjugates(eps->rg,ctx->isreal,&ctx->useconj);
272:     SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)eps,&ctx->contour);
273:   }

275:   EPSAllocateSolution(eps,0);
276:   BVGetRandomContext(eps->V,&rand);  /* make sure the random context is available when duplicating */
277:   if (ctx->weight) { PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma); }
278:   PetscMalloc4(ctx->N,&ctx->weight,ctx->N+1,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma);
279:   PetscLogObjectMemory((PetscObject)eps,3*ctx->N*sizeof(PetscScalar)+ctx->L_max*ctx->N*sizeof(PetscReal));

281:   /* allocate basis vectors */
282:   BVDestroy(&ctx->S);
283:   BVDuplicateResize(eps->V,ctx->L_max*ctx->M,&ctx->S);
284:   PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->S);
285:   BVDestroy(&ctx->V);
286:   BVDuplicateResize(eps->V,ctx->L_max,&ctx->V);
287:   PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->V);

289:   STGetMatrix(eps->st,0,&A[0]);
290:   MatIsShell(A[0],&flg);
291:   if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Matrix type shell is not supported in this solver");
292:   if (eps->isgeneralized) { STGetMatrix(eps->st,0,&A[1]); }

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

297:   contour = ctx->contour;
298:   SlepcContourRedundantMat(contour,eps->isgeneralized?2:1,A);
299:   if (contour->pA) {
300:     BVGetColumn(ctx->V,0,&v0);
301:     SlepcContourScatterCreate(contour,v0);
302:     BVRestoreColumn(ctx->V,0,&v0);
303:     BVDestroy(&ctx->pV);
304:     BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV);
305:     BVSetSizesFromVec(ctx->pV,contour->xsub,eps->n);
306:     BVSetFromOptions(ctx->pV);
307:     BVResize(ctx->pV,ctx->L_max,PETSC_FALSE);
308:     PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->pV);
309:   }

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

314:   BVDestroy(&ctx->Y);
315:   if (contour->pA) {
316:     BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y);
317:     BVSetSizesFromVec(ctx->Y,contour->xsub,eps->n);
318:     BVSetFromOptions(ctx->Y);
319:     BVResize(ctx->Y,contour->npoints*ctx->L_max,PETSC_FALSE);
320:   } else {
321:     BVDuplicateResize(eps->V,contour->npoints*ctx->L_max,&ctx->Y);
322:   }
323:   PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->Y);

325:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
326:     DSSetType(eps->ds,DSGNHEP);
327:   } else if (eps->isgeneralized) {
328:     if (eps->ishermitian && eps->ispositive) {
329:       DSSetType(eps->ds,DSGHEP);
330:     } else {
331:       DSSetType(eps->ds,DSGNHEP);
332:     }
333:   } else {
334:     if (eps->ishermitian) {
335:       DSSetType(eps->ds,DSHEP);
336:     } else {
337:       DSSetType(eps->ds,DSNHEP);
338:     }
339:   }
340:   DSAllocate(eps->ds,eps->ncv);

342: #if !defined(PETSC_USE_COMPLEX)
343:   EPSSetWorkVecs(eps,3);
344:   if (!eps->ishermitian) { PetscInfo(eps,"Warning: complex eigenvalues are not calculated exactly without --with-scalar-type=complex in PETSc\n"); }
345: #else
346:   EPSSetWorkVecs(eps,2);
347: #endif
348:   return(0);
349: }

351: PetscErrorCode EPSSetUpSort_CISS(EPS eps)
352: {
354:   SlepcSC        sc;

357:   /* fill sorting criterion context */
358:   eps->sc->comparison    = SlepcCompareSmallestReal;
359:   eps->sc->comparisonctx = NULL;
360:   eps->sc->map           = NULL;
361:   eps->sc->mapobj        = NULL;

363:   /* fill sorting criterion for DS */
364:   DSGetSlepcSC(eps->ds,&sc);
365:   sc->comparison    = SlepcCompareLargestMagnitude;
366:   sc->comparisonctx = NULL;
367:   sc->map           = NULL;
368:   sc->mapobj        = NULL;
369:   return(0);
370: }

372: PetscErrorCode EPSSolve_CISS(EPS eps)
373: {
374:   PetscErrorCode   ierr;
375:   EPS_CISS         *ctx = (EPS_CISS*)eps->data;
376:   SlepcContourData contour = ctx->contour;
377:   Mat              A,B,X,M,pA,pB;
378:   PetscInt         i,j,ld,nmat,L_add=0,nv=0,L_base=ctx->L,inner,nlocal,*inside;
379:   PetscScalar      *Mu,*H0,*H1=NULL,*rr,*temp;
380:   PetscReal        error,max_error,norm;
381:   PetscBool        *fl1;
382:   Vec              si,si1=NULL,w[3];
383:   PetscRandom      rand;
384: #if defined(PETSC_USE_COMPLEX)
385:   PetscBool        isellipse;
386:   PetscReal        est_eig,eta;
387: #else
388:   PetscReal        normi;
389: #endif

392:   w[0] = eps->work[0];
393: #if defined(PETSC_USE_COMPLEX)
394:   w[1] = NULL;
395: #else
396:   w[1] = eps->work[2];
397: #endif
398:   w[2] = eps->work[1];
399:   VecGetLocalSize(w[0],&nlocal);
400:   DSGetLeadingDimension(eps->ds,&ld);
401:   STGetNumMatrices(eps->st,&nmat);
402:   STGetMatrix(eps->st,0,&A);
403:   if (nmat>1) { STGetMatrix(eps->st,1,&B); }
404:   else B = NULL;
405:   RGComputeQuadrature(eps->rg,ctx->quad==EPS_CISS_QUADRULE_CHEBYSHEV?RG_QUADRULE_CHEBYSHEV:RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight);
406:   BVSetActiveColumns(ctx->V,0,ctx->L);
407:   BVSetRandomSign(ctx->V);
408:   BVGetRandomContext(ctx->V,&rand);

410:   if (contour->pA) {
411:     BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
412:     EPSCISSSolveSystem(eps,contour->pA[0],contour->pA[1],ctx->pV,0,ctx->L,PETSC_TRUE);
413:   } else {
414:     EPSCISSSolveSystem(eps,A,B,ctx->V,0,ctx->L,PETSC_TRUE);
415:   }
416: #if defined(PETSC_USE_COMPLEX)
417:   PetscObjectTypeCompare((PetscObject)eps->rg,RGELLIPSE,&isellipse);
418:   if (isellipse) {
419:     BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L_max,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig);
420:     PetscInfo1(eps,"Estimated eigenvalue count: %f\n",(double)est_eig);
421:     eta = PetscPowReal(10.0,-PetscLog10Real(eps->tol)/ctx->N);
422:     L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
423:     if (L_add>ctx->L_max-ctx->L) {
424:       PetscInfo(eps,"Number of eigenvalues inside the contour path may be too large\n");
425:       L_add = ctx->L_max-ctx->L;
426:     }
427:   }
428: #endif
429:   if (L_add>0) {
430:     PetscInfo2(eps,"Changing L %D -> %D by Estimate #Eig\n",ctx->L,ctx->L+L_add);
431:     BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
432:     BVSetRandomSign(ctx->V);
433:     if (contour->pA) {
434:       BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
435:       EPSCISSSolveSystem(eps,contour->pA[0],contour->pA[1],ctx->pV,ctx->L,ctx->L+L_add,PETSC_FALSE);
436:     } else {
437:       EPSCISSSolveSystem(eps,A,B,ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
438:     }
439:     ctx->L += L_add;
440:   }
441:   PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
442:   for (i=0;i<ctx->refine_blocksize;i++) {
443:     BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj);
444:     CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
445:     PetscLogEventBegin(EPS_CISS_SVD,eps,0,0,0);
446:     SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
447:     PetscLogEventEnd(EPS_CISS_SVD,eps,0,0,0);
448:     if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
449:     L_add = L_base;
450:     if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
451:     PetscInfo2(eps,"Changing L %D -> %D by SVD(H0)\n",ctx->L,ctx->L+L_add);
452:     BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add);
453:     BVSetRandomSign(ctx->V);
454:     if (contour->pA) {
455:       BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
456:       EPSCISSSolveSystem(eps,contour->pA[0],contour->pA[1],ctx->pV,ctx->L,ctx->L+L_add,PETSC_FALSE);
457:     } else {
458:       EPSCISSSolveSystem(eps,A,B,ctx->V,ctx->L,ctx->L+L_add,PETSC_FALSE);
459:     }
460:     ctx->L += L_add;
461:     if (L_add) {
462:       PetscFree2(Mu,H0);
463:       PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0);
464:     }
465:   }
466:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
467:     PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1);
468:   }

470:   while (eps->reason == EPS_CONVERGED_ITERATING) {
471:     eps->its++;
472:     for (inner=0;inner<=ctx->refine_inner;inner++) {
473:       if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
474:         BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj);
475:         CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
476:         PetscLogEventBegin(EPS_CISS_SVD,eps,0,0,0);
477:         SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv);
478:         PetscLogEventEnd(EPS_CISS_SVD,eps,0,0,0);
479:         break;
480:       } else {
481:         BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
482:         BVSetActiveColumns(ctx->S,0,ctx->L);
483:         BVSetActiveColumns(ctx->V,0,ctx->L);
484:         BVCopy(ctx->S,ctx->V);
485:         BVSVDAndRank(ctx->S,ctx->M,ctx->L,ctx->delta,BV_SVD_METHOD_REFINE,H0,ctx->sigma,&nv);
486:         if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
487:           if (contour->pA) {
488:             BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
489:             EPSCISSSolveSystem(eps,contour->pA[0],contour->pA[1],ctx->pV,0,ctx->L,PETSC_FALSE);
490:           } else {
491:             EPSCISSSolveSystem(eps,A,B,ctx->V,0,ctx->L,PETSC_FALSE);
492:           }
493:         } else break;
494:       }
495:     }
496:     eps->nconv = 0;
497:     if (nv == 0) eps->reason = EPS_CONVERGED_TOL;
498:     else {
499:       DSSetDimensions(eps->ds,nv,0,0);
500:       DSSetState(eps->ds,DS_STATE_RAW);

502:       if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
503:         CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0);
504:         CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1);
505:         DSGetArray(eps->ds,DS_MAT_A,&temp);
506:         for (j=0;j<nv;j++) {
507:           for (i=0;i<nv;i++) {
508:             temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
509:           }
510:         }
511:         DSRestoreArray(eps->ds,DS_MAT_A,&temp);
512:         DSGetArray(eps->ds,DS_MAT_B,&temp);
513:         for (j=0;j<nv;j++) {
514:           for (i=0;i<nv;i++) {
515:             temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
516:           }
517:         }
518:         DSRestoreArray(eps->ds,DS_MAT_B,&temp);
519:       } else {
520:         BVSetActiveColumns(ctx->S,0,nv);
521:         DSGetMat(eps->ds,DS_MAT_A,&pA);
522:         MatZeroEntries(pA);
523:         BVMatProject(ctx->S,A,ctx->S,pA);
524:         DSRestoreMat(eps->ds,DS_MAT_A,&pA);
525:         if (B) {
526:           DSGetMat(eps->ds,DS_MAT_B,&pB);
527:           MatZeroEntries(pB);
528:           BVMatProject(ctx->S,B,ctx->S,pB);
529:           DSRestoreMat(eps->ds,DS_MAT_B,&pB);
530:         }
531:       }

533:       DSSolve(eps->ds,eps->eigr,eps->eigi);
534:       DSSynchronize(eps->ds,eps->eigr,eps->eigi);

536:       PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr);
537:       rescale_eig(eps,nv);
538:       DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
539:       DSGetMat(eps->ds,DS_MAT_X,&X);
540:       SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1);
541:       MatDestroy(&X);
542:       RGCheckInside(eps->rg,nv,eps->eigr,eps->eigi,inside);
543:       for (i=0;i<nv;i++) {
544:         if (fl1[i] && inside[i]>=0) {
545:           rr[i] = 1.0;
546:           eps->nconv++;
547:         } else rr[i] = 0.0;
548:       }
549:       DSSort(eps->ds,eps->eigr,eps->eigi,rr,NULL,&eps->nconv);
550:       DSSynchronize(eps->ds,eps->eigr,eps->eigi);
551:       rescale_eig(eps,nv);
552:       PetscFree3(fl1,inside,rr);
553:       BVSetActiveColumns(eps->V,0,nv);
554:       if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
555:         BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L_max,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj);
556:         BVSetActiveColumns(ctx->S,0,ctx->L);
557:         BVCopy(ctx->S,ctx->V);
558:         BVSetActiveColumns(ctx->S,0,nv);
559:       }
560:       BVCopy(ctx->S,eps->V);

562:       DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
563:       DSGetMat(eps->ds,DS_MAT_X,&X);
564:       BVMultInPlace(ctx->S,X,0,eps->nconv);
565:       if (eps->ishermitian) {
566:         BVMultInPlace(eps->V,X,0,eps->nconv);
567:       }
568:       MatDestroy(&X);
569:       max_error = 0.0;
570:       for (i=0;i<eps->nconv;i++) {
571:         BVGetColumn(ctx->S,i,&si);
572: #if !defined(PETSC_USE_COMPLEX)
573:         if (eps->eigi[i]!=0.0) { BVGetColumn(ctx->S,i+1,&si1); }
574: #endif
575:         EPSComputeResidualNorm_Private(eps,PETSC_FALSE,eps->eigr[i],eps->eigi[i],si,si1,w,&error);
576:         if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {  /* vector is not normalized */
577:           VecNorm(si,NORM_2,&norm);
578: #if !defined(PETSC_USE_COMPLEX)
579:           if (eps->eigi[i]!=0.0) {
580:             VecNorm(si1,NORM_2,&normi);
581:             norm = SlepcAbsEigenvalue(norm,normi);
582:           }
583: #endif
584:           error /= norm;
585:         }
586:         (*eps->converged)(eps,eps->eigr[i],eps->eigi[i],error,&error,eps->convergedctx);
587:         BVRestoreColumn(ctx->S,i,&si);
588: #if !defined(PETSC_USE_COMPLEX)
589:         if (eps->eigi[i]!=0.0) {
590:           BVRestoreColumn(ctx->S,i+1,&si1);
591:           i++;
592:         }
593: #endif
594:         max_error = PetscMax(max_error,error);
595:       }

597:       if (max_error <= eps->tol) eps->reason = EPS_CONVERGED_TOL;
598:       else if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
599:       else {
600:         if (eps->nconv > ctx->L) nv = eps->nconv;
601:         else if (ctx->L > nv) nv = ctx->L;
602:         MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M);
603:         MatSetRandom(M,rand);
604:         BVSetActiveColumns(ctx->S,0,nv);
605:         BVMultInPlace(ctx->S,M,0,ctx->L);
606:         MatDestroy(&M);
607:         BVSetActiveColumns(ctx->S,0,ctx->L);
608:         BVSetActiveColumns(ctx->V,0,ctx->L);
609:         BVCopy(ctx->S,ctx->V);
610:         if (contour->pA) {
611:           BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup);
612:           EPSCISSSolveSystem(eps,contour->pA[0],contour->pA[1],ctx->pV,0,ctx->L,PETSC_FALSE);
613:         } else {
614:           EPSCISSSolveSystem(eps,A,B,ctx->V,0,ctx->L,PETSC_FALSE);
615:         }
616:       }
617:     }
618:   }
619:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
620:     PetscFree(H1);
621:   }
622:   PetscFree2(Mu,H0);
623:   return(0);
624: }

626: PetscErrorCode EPSComputeVectors_CISS(EPS eps)
627: {
629:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
630:   PetscInt       n;
631:   Mat            Z,B=NULL;

634:   if (eps->ishermitian) {
635:     if (eps->isgeneralized && !eps->ispositive) {
636:       EPSComputeVectors_Indefinite(eps);
637:     } else {
638:       EPSComputeVectors_Hermitian(eps);
639:     }
640:     if (eps->isgeneralized && eps->ispositive && ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
641:       /* normalize to have unit B-norm */
642:       STGetMatrix(eps->st,1,&B);
643:       BVSetMatrix(eps->V,B,PETSC_FALSE);
644:       BVNormalize(eps->V,NULL);
645:       BVSetMatrix(eps->V,NULL,PETSC_FALSE);
646:     }
647:     return(0);
648:   }
649:   DSGetDimensions(eps->ds,&n,NULL,NULL,NULL);
650:   BVSetActiveColumns(eps->V,0,n);

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

655:   /* V = V * Z */
656:   DSGetMat(eps->ds,DS_MAT_X,&Z);
657:   BVMultInPlace(eps->V,Z,0,n);
658:   MatDestroy(&Z);
659:   BVSetActiveColumns(eps->V,0,eps->nconv);

661:   /* normalize */
662:   if (ctx->extraction == EPS_CISS_EXTRACTION_HANKEL) {
663:     BVNormalize(eps->V,NULL);
664:   }
665:   return(0);
666: }

668: static PetscErrorCode EPSCISSSetSizes_CISS(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
669: {
671:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
672:   PetscInt       oN,oL,oM,oLmax,onpart;

675:   oN = ctx->N;
676:   if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
677:     if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
678:   } else {
679:     if (ip<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
680:     if (ip%2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
681:     if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
682:   }
683:   oL = ctx->L;
684:   if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
685:     ctx->L = 16;
686:   } else {
687:     if (bs<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
688:     ctx->L = bs;
689:   }
690:   oM = ctx->M;
691:   if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
692:     ctx->M = ctx->N/4;
693:   } else {
694:     if (ms<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
695:     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");
696:     ctx->M = ms;
697:   }
698:   onpart = ctx->npart;
699:   if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
700:     ctx->npart = 1;
701:   } else {
702:     if (npart<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The npart argument must be > 0");
703:     ctx->npart = npart;
704:   }
705:   oLmax = ctx->L_max;
706:   if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
707:     ctx->L_max = 64;
708:   } else {
709:     if (bsmax<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
710:     ctx->L_max = PetscMax(bsmax,ctx->L);
711:   }
712:   if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
713:     SlepcContourDataDestroy(&ctx->contour);
714:     PetscInfo(eps,"Resetting the contour data structure due to a change of parameters\n");
715:     eps->state = EPS_STATE_INITIAL;
716:   }
717:   ctx->isreal = realmats;
718:   if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) eps->state = EPS_STATE_INITIAL;
719:   return(0);
720: }

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

725:    Logically Collective on eps

727:    Input Parameters:
728: +  eps   - the eigenproblem solver context
729: .  ip    - number of integration points
730: .  bs    - block size
731: .  ms    - moment size
732: .  npart - number of partitions when splitting the communicator
733: .  bsmax - max block size
734: -  realmats - A and B are real

736:    Options Database Keys:
737: +  -eps_ciss_integration_points - Sets the number of integration points
738: .  -eps_ciss_blocksize - Sets the block size
739: .  -eps_ciss_moments - Sets the moment size
740: .  -eps_ciss_partitions - Sets the number of partitions
741: .  -eps_ciss_maxblocksize - Sets the maximum block size
742: -  -eps_ciss_realmats - A and B are real

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

749:    Level: advanced

751: .seealso: EPSCISSGetSizes()
752: @*/
753: PetscErrorCode EPSCISSSetSizes(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
754: {

765:   PetscTryMethod(eps,"EPSCISSSetSizes_C",(EPS,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(eps,ip,bs,ms,npart,bsmax,realmats));
766:   return(0);
767: }

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

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

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

786:    Not Collective

788:    Input Parameter:
789: .  eps - the eigenproblem solver context

791:    Output Parameters:
792: +  ip    - number of integration points
793: .  bs    - block size
794: .  ms    - moment size
795: .  npart - number of partitions when splitting the communicator
796: .  bsmax - max block size
797: -  realmats - A and B are real

799:    Level: advanced

801: .seealso: EPSCISSSetSizes()
802: @*/
803: PetscErrorCode EPSCISSGetSizes(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
804: {

809:   PetscUseMethod(eps,"EPSCISSGetSizes_C",(EPS,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(eps,ip,bs,ms,npart,bsmax,realmats));
810:   return(0);
811: }

813: static PetscErrorCode EPSCISSSetThreshold_CISS(EPS eps,PetscReal delta,PetscReal spur)
814: {
815:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

818:   if (delta == PETSC_DEFAULT) {
819:     ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
820:   } else {
821:     if (delta<=0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
822:     ctx->delta = delta;
823:   }
824:   if (spur == PETSC_DEFAULT) {
825:     ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
826:   } else {
827:     if (spur<=0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
828:     ctx->spurious_threshold = spur;
829:   }
830:   return(0);
831: }

833: /*@
834:    EPSCISSSetThreshold - Sets the values of various threshold parameters in
835:    the CISS solver.

837:    Logically Collective on eps

839:    Input Parameters:
840: +  eps   - the eigenproblem solver context
841: .  delta - threshold for numerical rank
842: -  spur  - spurious threshold (to discard spurious eigenpairs)

844:    Options Database Keys:
845: +  -eps_ciss_delta - Sets the delta
846: -  -eps_ciss_spurious_threshold - Sets the spurious threshold

848:    Level: advanced

850: .seealso: EPSCISSGetThreshold()
851: @*/
852: PetscErrorCode EPSCISSSetThreshold(EPS eps,PetscReal delta,PetscReal spur)
853: {

860:   PetscTryMethod(eps,"EPSCISSSetThreshold_C",(EPS,PetscReal,PetscReal),(eps,delta,spur));
861:   return(0);
862: }

864: static PetscErrorCode EPSCISSGetThreshold_CISS(EPS eps,PetscReal *delta,PetscReal *spur)
865: {
866:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

869:   if (delta) *delta = ctx->delta;
870:   if (spur)  *spur = ctx->spurious_threshold;
871:   return(0);
872: }

874: /*@
875:    EPSCISSGetThreshold - Gets the values of various threshold parameters
876:    in the CISS solver.

878:    Not Collective

880:    Input Parameter:
881: .  eps - the eigenproblem solver context

883:    Output Parameters:
884: +  delta - threshold for numerical rank
885: -  spur  - spurious threshold (to discard spurious eigenpairs)

887:    Level: advanced

889: .seealso: EPSCISSSetThreshold()
890: @*/
891: PetscErrorCode EPSCISSGetThreshold(EPS eps,PetscReal *delta,PetscReal *spur)
892: {

897:   PetscUseMethod(eps,"EPSCISSGetThreshold_C",(EPS,PetscReal*,PetscReal*),(eps,delta,spur));
898:   return(0);
899: }

901: static PetscErrorCode EPSCISSSetRefinement_CISS(EPS eps,PetscInt inner,PetscInt blsize)
902: {
903:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

906:   if (inner == PETSC_DEFAULT) {
907:     ctx->refine_inner = 0;
908:   } else {
909:     if (inner<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
910:     ctx->refine_inner = inner;
911:   }
912:   if (blsize == PETSC_DEFAULT) {
913:     ctx->refine_blocksize = 0;
914:   } else {
915:     if (blsize<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
916:     ctx->refine_blocksize = blsize;
917:   }
918:   return(0);
919: }

921: /*@
922:    EPSCISSSetRefinement - Sets the values of various refinement parameters
923:    in the CISS solver.

925:    Logically Collective on eps

927:    Input Parameters:
928: +  eps    - the eigenproblem solver context
929: .  inner  - number of iterative refinement iterations (inner loop)
930: -  blsize - number of iterative refinement iterations (blocksize loop)

932:    Options Database Keys:
933: +  -eps_ciss_refine_inner - Sets number of inner iterations
934: -  -eps_ciss_refine_blocksize - Sets number of blocksize iterations

936:    Level: advanced

938: .seealso: EPSCISSGetRefinement()
939: @*/
940: PetscErrorCode EPSCISSSetRefinement(EPS eps,PetscInt inner,PetscInt blsize)
941: {

948:   PetscTryMethod(eps,"EPSCISSSetRefinement_C",(EPS,PetscInt,PetscInt),(eps,inner,blsize));
949:   return(0);
950: }

952: static PetscErrorCode EPSCISSGetRefinement_CISS(EPS eps,PetscInt *inner,PetscInt *blsize)
953: {
954:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

957:   if (inner)  *inner = ctx->refine_inner;
958:   if (blsize) *blsize = ctx->refine_blocksize;
959:   return(0);
960: }

962: /*@
963:    EPSCISSGetRefinement - Gets the values of various refinement parameters
964:    in the CISS solver.

966:    Not Collective

968:    Input Parameter:
969: .  eps - the eigenproblem solver context

971:    Output Parameters:
972: +  inner  - number of iterative refinement iterations (inner loop)
973: -  blsize - number of iterative refinement iterations (blocksize loop)

975:    Level: advanced

977: .seealso: EPSCISSSetRefinement()
978: @*/
979: PetscErrorCode EPSCISSGetRefinement(EPS eps, PetscInt *inner, PetscInt *blsize)
980: {

985:   PetscUseMethod(eps,"EPSCISSGetRefinement_C",(EPS,PetscInt*,PetscInt*),(eps,inner,blsize));
986:   return(0);
987: }

989: static PetscErrorCode EPSCISSSetUseST_CISS(EPS eps,PetscBool usest)
990: {
991:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

994:   ctx->usest     = usest;
995:   ctx->usest_set = PETSC_TRUE;
996:   eps->state     = EPS_STATE_INITIAL;
997:   return(0);
998: }

1000: /*@
1001:    EPSCISSSetUseST - Sets a flag indicating that the CISS solver will
1002:    use the ST object for the linear solves.

1004:    Logically Collective on eps

1006:    Input Parameters:
1007: +  eps    - the eigenproblem solver context
1008: -  usest  - boolean flag to use the ST object or not

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

1013:    Level: advanced

1015: .seealso: EPSCISSGetUseST()
1016: @*/
1017: PetscErrorCode EPSCISSSetUseST(EPS eps,PetscBool usest)
1018: {

1024:   PetscTryMethod(eps,"EPSCISSSetUseST_C",(EPS,PetscBool),(eps,usest));
1025:   return(0);
1026: }

1028: static PetscErrorCode EPSCISSGetUseST_CISS(EPS eps,PetscBool *usest)
1029: {
1030:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1033:   *usest = ctx->usest;
1034:   return(0);
1035: }

1037: /*@
1038:    EPSCISSGetUseST - Gets the flag for using the ST object
1039:    in the CISS solver.

1041:    Not Collective

1043:    Input Parameter:
1044: .  eps - the eigenproblem solver context

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

1049:    Level: advanced

1051: .seealso: EPSCISSSetUseST()
1052: @*/
1053: PetscErrorCode EPSCISSGetUseST(EPS eps,PetscBool *usest)
1054: {

1060:   PetscUseMethod(eps,"EPSCISSGetUseST_C",(EPS,PetscBool*),(eps,usest));
1061:   return(0);
1062: }

1064: static PetscErrorCode EPSCISSSetQuadRule_CISS(EPS eps,EPSCISSQuadRule quad)
1065: {
1066:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1069:   if (ctx->quad != quad) {
1070:     ctx->quad  = quad;
1071:     eps->state = EPS_STATE_INITIAL;
1072:   }
1073:   return(0);
1074: }

1076: /*@
1077:    EPSCISSSetQuadRule - Sets the quadrature rule used in the CISS solver.

1079:    Logically Collective on eps

1081:    Input Parameters:
1082: +  eps  - the eigenproblem solver context
1083: -  quad - the quadrature rule

1085:    Options Database Key:
1086: .  -eps_ciss_quadrule - Sets the quadrature rule (either 'trapezoidal' or
1087:                            'chebyshev')

1089:    Notes:
1090:    By default, the trapezoidal rule is used (EPS_CISS_QUADRULE_TRAPEZOIDAL).

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

1095:    Level: advanced

1097: .seealso: EPSCISSGetQuadRule(), EPSCISSQuadRule
1098: @*/
1099: PetscErrorCode EPSCISSSetQuadRule(EPS eps,EPSCISSQuadRule quad)
1100: {

1106:   PetscTryMethod(eps,"EPSCISSSetQuadRule_C",(EPS,EPSCISSQuadRule),(eps,quad));
1107:   return(0);
1108: }

1110: static PetscErrorCode EPSCISSGetQuadRule_CISS(EPS eps,EPSCISSQuadRule *quad)
1111: {
1112:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1115:   *quad = ctx->quad;
1116:   return(0);
1117: }

1119: /*@
1120:    EPSCISSGetQuadRule - Gets the quadrature rule used in the CISS solver.

1122:    Not Collective

1124:    Input Parameter:
1125: .  eps - the eigenproblem solver context

1127:    Output Parameters:
1128: .  quad - quadrature rule

1130:    Level: advanced

1132: .seealso: EPSCISSSetQuadRule() EPSCISSQuadRule
1133: @*/
1134: PetscErrorCode EPSCISSGetQuadRule(EPS eps,EPSCISSQuadRule *quad)
1135: {

1141:   PetscUseMethod(eps,"EPSCISSGetQuadRule_C",(EPS,EPSCISSQuadRule*),(eps,quad));
1142:   return(0);
1143: }

1145: static PetscErrorCode EPSCISSSetExtraction_CISS(EPS eps,EPSCISSExtraction extraction)
1146: {
1147:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1150:   if (ctx->extraction != extraction) {
1151:     ctx->extraction = extraction;
1152:     eps->state      = EPS_STATE_INITIAL;
1153:   }
1154:   return(0);
1155: }

1157: /*@
1158:    EPSCISSSetExtraction - Sets the extraction technique used in the CISS solver.

1160:    Logically Collective on eps

1162:    Input Parameters:
1163: +  eps        - the eigenproblem solver context
1164: -  extraction - the extraction technique

1166:    Options Database Key:
1167: .  -eps_ciss_extraction - Sets the extraction technique (either 'ritz' or
1168:                            'hankel')

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

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

1176:    Level: advanced

1178: .seealso: EPSCISSGetExtraction(), EPSCISSExtraction
1179: @*/
1180: PetscErrorCode EPSCISSSetExtraction(EPS eps,EPSCISSExtraction extraction)
1181: {

1187:   PetscTryMethod(eps,"EPSCISSSetExtraction_C",(EPS,EPSCISSExtraction),(eps,extraction));
1188:   return(0);
1189: }

1191: static PetscErrorCode EPSCISSGetExtraction_CISS(EPS eps,EPSCISSExtraction *extraction)
1192: {
1193:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1196:   *extraction = ctx->extraction;
1197:   return(0);
1198: }

1200: /*@
1201:    EPSCISSGetExtraction - Gets the extraction technique used in the CISS solver.

1203:    Not Collective

1205:    Input Parameter:
1206: .  eps - the eigenproblem solver context

1208:    Output Parameters:
1209: .  extraction - extraction technique

1211:    Level: advanced

1213: .seealso: EPSCISSSetExtraction() EPSCISSExtraction
1214: @*/
1215: PetscErrorCode EPSCISSGetExtraction(EPS eps,EPSCISSExtraction *extraction)
1216: {

1222:   PetscUseMethod(eps,"EPSCISSGetExtraction_C",(EPS,EPSCISSExtraction*),(eps,extraction));
1223:   return(0);
1224: }

1226: static PetscErrorCode EPSCISSGetKSPs_CISS(EPS eps,PetscInt *nsolve,KSP **ksp)
1227: {
1228:   PetscErrorCode   ierr;
1229:   EPS_CISS         *ctx = (EPS_CISS*)eps->data;
1230:   SlepcContourData contour;
1231:   PetscInt         i;
1232:   PC               pc;

1235:   if (!ctx->contour) {  /* initialize contour data structure first */
1236:     RGCanUseConjugates(eps->rg,ctx->isreal,&ctx->useconj);
1237:     SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)eps,&ctx->contour);
1238:   }
1239:   contour = ctx->contour;
1240:   if (!contour->ksp) {
1241:     PetscMalloc1(contour->npoints,&contour->ksp);
1242:     for (i=0;i<contour->npoints;i++) {
1243:       KSPCreate(PetscSubcommChild(contour->subcomm),&contour->ksp[i]);
1244:       PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)eps,1);
1245:       KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)eps)->prefix);
1246:       KSPAppendOptionsPrefix(contour->ksp[i],"eps_ciss_");
1247:       PetscLogObjectParent((PetscObject)eps,(PetscObject)contour->ksp[i]);
1248:       PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)eps)->options);
1249:       KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE);
1250:       KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(eps->tol),PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
1251:       KSPGetPC(contour->ksp[i],&pc);
1252:       KSPSetType(contour->ksp[i],KSPPREONLY);
1253:       PCSetType(pc,PCLU);
1254:     }
1255:   }
1256:   if (nsolve) *nsolve = contour->npoints;
1257:   if (ksp)    *ksp    = contour->ksp;
1258:   return(0);
1259: }

1261: /*@C
1262:    EPSCISSGetKSPs - Retrieve the array of linear solver objects associated with
1263:    the CISS solver.

1265:    Not Collective

1267:    Input Parameter:
1268: .  eps - the eigenproblem solver solver

1270:    Output Parameters:
1271: +  nsolve - number of solver objects
1272: -  ksp - array of linear solver object

1274:    Notes:
1275:    The number of KSP solvers is equal to the number of integration points divided by
1276:    the number of partitions. This value is halved in the case of real matrices with
1277:    a region centered at the real axis.

1279:    Level: advanced

1281: .seealso: EPSCISSSetSizes()
1282: @*/
1283: PetscErrorCode EPSCISSGetKSPs(EPS eps,PetscInt *nsolve,KSP **ksp)
1284: {

1289:   PetscUseMethod(eps,"EPSCISSGetKSPs_C",(EPS,PetscInt*,KSP**),(eps,nsolve,ksp));
1290:   return(0);
1291: }

1293: PetscErrorCode EPSReset_CISS(EPS eps)
1294: {
1296:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

1299:   BVDestroy(&ctx->S);
1300:   BVDestroy(&ctx->V);
1301:   BVDestroy(&ctx->Y);
1302:   if (!ctx->usest) {
1303:     SlepcContourDataReset(ctx->contour);
1304:   }
1305:   BVDestroy(&ctx->pV);
1306:   return(0);
1307: }

1309: PetscErrorCode EPSSetFromOptions_CISS(PetscOptionItems *PetscOptionsObject,EPS eps)
1310: {
1311:   PetscErrorCode    ierr;
1312:   PetscReal         r3,r4;
1313:   PetscInt          i,i1,i2,i3,i4,i5,i6,i7;
1314:   PetscBool         b1,b2,flg,flg2,flg3,flg4,flg5,flg6;
1315:   EPS_CISS          *ctx = (EPS_CISS*)eps->data;
1316:   EPSCISSQuadRule   quad;
1317:   EPSCISSExtraction extraction;

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

1322:     EPSCISSGetSizes(eps,&i1,&i2,&i3,&i4,&i5,&b1);
1323:     PetscOptionsInt("-eps_ciss_integration_points","Number of integration points","EPSCISSSetSizes",i1,&i1,&flg);
1324:     PetscOptionsInt("-eps_ciss_blocksize","Block size","EPSCISSSetSizes",i2,&i2,&flg2);
1325:     PetscOptionsInt("-eps_ciss_moments","Moment size","EPSCISSSetSizes",i3,&i3,&flg3);
1326:     PetscOptionsInt("-eps_ciss_partitions","Number of partitions","EPSCISSSetSizes",i4,&i4,&flg4);
1327:     PetscOptionsInt("-eps_ciss_maxblocksize","Maximum block size","EPSCISSSetSizes",i5,&i5,&flg5);
1328:     PetscOptionsBool("-eps_ciss_realmats","True if A and B are real","EPSCISSSetSizes",b1,&b1,&flg6);
1329:     if (flg || flg2 || flg3 || flg4 || flg5 || flg6) { EPSCISSSetSizes(eps,i1,i2,i3,i4,i5,b1); }

1331:     EPSCISSGetThreshold(eps,&r3,&r4);
1332:     PetscOptionsReal("-eps_ciss_delta","Threshold for numerical rank","EPSCISSSetThreshold",r3,&r3,&flg);
1333:     PetscOptionsReal("-eps_ciss_spurious_threshold","Threshold for the spurious eigenpairs","EPSCISSSetThreshold",r4,&r4,&flg2);
1334:     if (flg || flg2) { EPSCISSSetThreshold(eps,r3,r4); }

1336:     EPSCISSGetRefinement(eps,&i6,&i7);
1337:     PetscOptionsInt("-eps_ciss_refine_inner","Number of inner iterative refinement iterations","EPSCISSSetRefinement",i6,&i6,&flg);
1338:     PetscOptionsInt("-eps_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","EPSCISSSetRefinement",i7,&i7,&flg2);
1339:     if (flg || flg2) { EPSCISSSetRefinement(eps,i6,i7); }

1341:     EPSCISSGetUseST(eps,&b2);
1342:     PetscOptionsBool("-eps_ciss_usest","Use ST for linear solves","EPSCISSSetUseST",b2,&b2,&flg);
1343:     if (flg) { EPSCISSSetUseST(eps,b2); }

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

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

1351:   PetscOptionsTail();

1353:   if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
1354:   RGSetFromOptions(eps->rg); /* this is necessary here to set useconj */
1355:   if (!ctx->contour || !ctx->contour->ksp) { EPSCISSGetKSPs(eps,NULL,NULL); }
1356:   for (i=0;i<ctx->contour->npoints;i++) {
1357:     KSPSetFromOptions(ctx->contour->ksp[i]);
1358:   }
1359:   PetscSubcommSetFromOptions(ctx->contour->subcomm);
1360:   return(0);
1361: }

1363: PetscErrorCode EPSDestroy_CISS(EPS eps)
1364: {
1366:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

1369:   SlepcContourDataDestroy(&ctx->contour);
1370:   PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma);
1371:   PetscFree(eps->data);
1372:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",NULL);
1373:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",NULL);
1374:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",NULL);
1375:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",NULL);
1376:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",NULL);
1377:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",NULL);
1378:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetUseST_C",NULL);
1379:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetUseST_C",NULL);
1380:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetQuadRule_C",NULL);
1381:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetQuadRule_C",NULL);
1382:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetExtraction_C",NULL);
1383:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetExtraction_C",NULL);
1384:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetKSPs_C",NULL);
1385:   return(0);
1386: }

1388: PetscErrorCode EPSView_CISS(EPS eps,PetscViewer viewer)
1389: {
1391:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1392:   PetscBool      isascii;
1393:   PetscViewer    sviewer;

1396:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1397:   if (isascii) {
1398:     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);
1399:     if (ctx->isreal) {
1400:       PetscViewerASCIIPrintf(viewer,"  exploiting symmetry of integration points\n");
1401:     }
1402:     PetscViewerASCIIPrintf(viewer,"  threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold);
1403:     PetscViewerASCIIPrintf(viewer,"  iterative refinement { inner: %D, blocksize: %D }\n",ctx->refine_inner, ctx->refine_blocksize);
1404:     PetscViewerASCIIPrintf(viewer,"  extraction: %s\n",EPSCISSExtractions[ctx->extraction]);
1405:     PetscViewerASCIIPrintf(viewer,"  quadrature rule: %s\n",EPSCISSQuadRules[ctx->quad]);
1406:     if (ctx->usest) {
1407:       PetscViewerASCIIPrintf(viewer,"  using ST for linear solves\n");
1408:     } else {
1409:       if (!ctx->contour || !ctx->contour->ksp) { EPSCISSGetKSPs(eps,NULL,NULL); }
1410:       PetscViewerASCIIPushTab(viewer);
1411:       if (ctx->npart>1 && ctx->contour->subcomm) {
1412:         PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1413:         if (!ctx->contour->subcomm->color) {
1414:           KSPView(ctx->contour->ksp[0],sviewer);
1415:         }
1416:         PetscViewerFlush(sviewer);
1417:         PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer);
1418:         PetscViewerFlush(viewer);
1419:         /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
1420:         PetscViewerASCIIPopSynchronized(viewer);
1421:       } else {
1422:         KSPView(ctx->contour->ksp[0],viewer);
1423:       }
1424:       PetscViewerASCIIPopTab(viewer);
1425:     }
1426:   }
1427:   return(0);
1428: }

1430: PetscErrorCode EPSSetDefaultST_CISS(EPS eps)
1431: {
1433:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1434:   PetscBool      usest = ctx->usest;

1437:   if (!((PetscObject)eps->st)->type_name) {
1438:     if (!ctx->usest_set) usest = (ctx->npart>1)? PETSC_FALSE: PETSC_TRUE;
1439:     if (usest) {
1440:       STSetType(eps->st,STSINVERT);
1441:     } else {
1442:       /* we are not going to use ST, so avoid factorizing the matrix */
1443:       STSetType(eps->st,STSHIFT);
1444:     }
1445:   }
1446:   return(0);
1447: }

1449: SLEPC_EXTERN PetscErrorCode EPSCreate_CISS(EPS eps)
1450: {
1452:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

1455:   PetscNewLog(eps,&ctx);
1456:   eps->data = ctx;

1458:   eps->useds = PETSC_TRUE;
1459:   eps->categ = EPS_CATEGORY_CONTOUR;

1461:   eps->ops->solve          = EPSSolve_CISS;
1462:   eps->ops->setup          = EPSSetUp_CISS;
1463:   eps->ops->setupsort      = EPSSetUpSort_CISS;
1464:   eps->ops->setfromoptions = EPSSetFromOptions_CISS;
1465:   eps->ops->destroy        = EPSDestroy_CISS;
1466:   eps->ops->reset          = EPSReset_CISS;
1467:   eps->ops->view           = EPSView_CISS;
1468:   eps->ops->computevectors = EPSComputeVectors_CISS;
1469:   eps->ops->setdefaultst   = EPSSetDefaultST_CISS;

1471:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",EPSCISSSetSizes_CISS);
1472:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",EPSCISSGetSizes_CISS);
1473:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",EPSCISSSetThreshold_CISS);
1474:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",EPSCISSGetThreshold_CISS);
1475:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",EPSCISSSetRefinement_CISS);
1476:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",EPSCISSGetRefinement_CISS);
1477:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetUseST_C",EPSCISSSetUseST_CISS);
1478:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetUseST_C",EPSCISSGetUseST_CISS);
1479:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetQuadRule_C",EPSCISSSetQuadRule_CISS);
1480:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetQuadRule_C",EPSCISSGetQuadRule_CISS);
1481:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetExtraction_C",EPSCISSSetExtraction_CISS);
1482:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetExtraction_C",EPSCISSGetExtraction_CISS);
1483:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetKSPs_C",EPSCISSGetKSPs_CISS);

1485:   /* set default values of parameters */
1486:   ctx->N                  = 32;
1487:   ctx->L                  = 16;
1488:   ctx->M                  = ctx->N/4;
1489:   ctx->delta              = SLEPC_DEFAULT_TOL*1e-4;
1490:   ctx->L_max              = 64;
1491:   ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
1492:   ctx->usest              = PETSC_TRUE;
1493:   ctx->usest_set          = PETSC_FALSE;
1494:   ctx->isreal             = PETSC_FALSE;
1495:   ctx->refine_inner       = 0;
1496:   ctx->refine_blocksize   = 0;
1497:   ctx->npart              = 1;
1498:   ctx->quad               = (EPSCISSQuadRule)0;
1499:   ctx->extraction         = EPS_CISS_EXTRACTION_RITZ;
1500:   return(0);
1501: }