Actual source code: ciss.c

  1: /*

  3:    SLEPc eigensolver: "ciss"

  5:    Method: Contour Integral Spectral Slicing

  7:    Algorithm:

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

 13:    Based on code contributed by Tetsuya Sakurai.

 15:    References:

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

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

 24:    Last update: Jun 2013

 26:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 27:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 28:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

 30:    This file is part of SLEPc.

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

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

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

 46: #include <slepc-private/epsimpl.h>                /*I "slepceps.h" I*/
 47: #include <slepcblaslapack.h>

 49: PetscErrorCode EPSSolve_CISS(EPS);

 51: typedef struct {
 52:   /* parameters */
 53:   PetscScalar center;     /* center of the region where to find eigenpairs (default: 0.0) */
 54:   PetscReal   radius;     /* radius of the region (1.0) */
 55:   PetscReal   vscale;     /* vertical scale of the region (1.0; 0.1 if spectrum real) */
 56:   PetscInt    N;          /* number of integration points (32) */
 57:   PetscInt    L;          /* block size (16) */
 58:   PetscInt    M;          /* moment degree (N/4 = 4) */
 59:   PetscReal   delta;      /* threshold of singular value (1e-12) */
 60:   PetscInt    npart;      /* number of partitions of the matrix (1) */
 61:   PetscReal   *sigma;     /* threshold for numerical rank */
 62:   PetscInt    L_max;      /* maximum number of columns of the source matrix V */
 63:   PetscReal   spurious_threshold; /* discard spurious eigenpairs */
 64:   PetscBool   isreal;     /* A and B are real */
 65:   PetscInt    refine_inner;
 66:   PetscInt    refine_outer;
 67:   PetscInt    refine_blocksize;
 68:   /* private data */
 69:   PetscInt    solver_comm_id;
 70:   PetscInt    num_solve_point;
 71:   PetscScalar *weight;
 72:   PetscScalar *omega;
 73:   PetscScalar *pp;
 74:   Vec         *V;
 75:   Vec         *Y;
 76:   Vec         *S;
 77:   KSP         *ksp;
 78:   PetscBool   useconj;
 79:   PetscReal   est_eig;
 80: } EPS_CISS;

 84: static PetscErrorCode SetSolverComm(EPS eps)
 85: {
 86:   EPS_CISS *ctx = (EPS_CISS*)eps->data;
 87:   PetscInt N = ctx->N;

 90:   if (ctx->useconj) N = N/2;
 91:   ctx->solver_comm_id = 0;
 92:   ctx->num_solve_point = N;
 93:   return(0);
 94: }

 98: static PetscErrorCode SetPathParameter(EPS eps)
 99: {
100:   EPS_CISS  *ctx = (EPS_CISS*)eps->data;
101:   PetscInt  i;
102:   PetscReal theta;

105:   for (i=0;i<ctx->N;i++){
106:     theta = ((2*PETSC_PI)/ctx->N)*(i+0.5);
107:     ctx->pp[i] = cos(theta) + PETSC_i*ctx->vscale*sin(theta);
108:     ctx->omega[i] = ctx->center + ctx->radius*ctx->pp[i];
109:     ctx->weight[i] = ctx->vscale*cos(theta) + PETSC_i*sin(theta);
110:   }
111:   return(0);
112: }

116: static PetscErrorCode CISSVecSetRandom(Vec x,PetscRandom rctx)
117: {
119:   PetscInt       j,nlocal;
120:   PetscScalar    *vdata;

123:   SlepcVecSetRandom(x,rctx);
124:   VecGetLocalSize(x,&nlocal);
125:   VecGetArray(x,&vdata);
126:   for (j=0;j<nlocal;j++) {
127:     vdata[j] = PetscRealPart(vdata[j]);
128:     if (PetscRealPart(vdata[j]) < 0.5) vdata[j] = -1.0;
129:     else vdata[j] = 1.0;
130:   }
131:   VecRestoreArray(x,&vdata);
132:   return(0);
133: }

137: static PetscErrorCode SolveLinearSystem(EPS eps)
138: {
140:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
141:   PetscInt       i,j,nmat,p_id;
142:   Mat            A,B,Fz;
143:   PC             pc;
144:   Vec            BV;

147:   STGetNumMatrices(eps->st,&nmat);
148:   STGetOperators(eps->st,0,&A);
149:   if (nmat>1) { STGetOperators(eps->st,1,&B); }
150:   else B = NULL;
151:   MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Fz);
152:   VecDuplicate(ctx->V[0],&BV);

154:   for (i=0;i<ctx->num_solve_point;i++) {
155:     p_id = ctx->solver_comm_id * ctx->num_solve_point + i;
156:     MatCopy(A,Fz,DIFFERENT_NONZERO_PATTERN);
157:     if (nmat>1) {
158:       MatAXPY(Fz,-ctx->omega[p_id],B,DIFFERENT_NONZERO_PATTERN);
159:     } else {
160:       MatShift(Fz,-ctx->omega[p_id]);
161:     }
162:     KSPSetOperators(ctx->ksp[i],Fz,Fz,SAME_NONZERO_PATTERN);
163:     KSPSetType(ctx->ksp[i],KSPPREONLY);
164:     KSPGetPC(ctx->ksp[i],&pc);
165:     PCSetType(pc,PCREDUNDANT);
166:     KSPSetFromOptions(ctx->ksp[i]);
167:     for (j=0;j<ctx->L;j++) {
168:       VecDuplicate(ctx->V[0],&ctx->Y[i*ctx->L_max+j]);
169:       PetscLogObjectParent(eps,ctx->Y[i*ctx->L_max+j]);
170:       if (nmat==2) {
171:         MatMult(B,ctx->V[j],BV);
172:         KSPSolve(ctx->ksp[i],BV,ctx->Y[i*ctx->L_max+j]);
173:       } else {
174:         KSPSolve(ctx->ksp[i],ctx->V[j],ctx->Y[i*ctx->L_max+j]);
175:       }
176:     }
177:   }
178:   MatDestroy(&Fz);
179:   VecDestroy(&BV);
180:   return(0);
181: }

185: static PetscErrorCode ConstructS(EPS eps,PetscInt M,Vec **S)
186: {
188:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
189:   PetscInt       i,j,k;
190:   Vec            v;
191:   PetscScalar    *ppk;

194:   VecDuplicateVecs(ctx->Y[0],M*ctx->L,S);
195:   PetscMalloc(ctx->num_solve_point*sizeof(PetscScalar),&ppk);
196:   for (i=0;i<ctx->num_solve_point;i++) ppk[i] = 1;
197:   VecDuplicate(ctx->Y[0],&v);
198:   for (k=0;k<M;k++) {
199:     for (j=0;j<ctx->L;j++) {
200:       VecSet(v,0);
201:       for (i=0;i<ctx->num_solve_point; i++) {
202:         VecAXPY(v,ppk[i]*ctx->weight[ctx->solver_comm_id*ctx->num_solve_point+i]/(PetscReal)ctx->N,ctx->Y[i*ctx->L_max+j]);
203:       }
204:       VecCopy(v,(*S)[k*ctx->L+j]);
205:     }
206:     for (i=0;i<ctx->num_solve_point;i++) {
207:       ppk[i] *= ctx->pp[ctx->solver_comm_id*ctx->num_solve_point+i];
208:     }
209:   }
210:   PetscFree(ppk);
211:   VecDestroy(&v);
212:   return(0);
213: }

217: static PetscErrorCode EstimateNumberEigs(EPS eps,Vec *S1,PetscInt *L_add)
218: {
220:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
221:   PetscInt       i,j,istart,p_start,p_end;
222:   PetscScalar    *data,*p_data,tmp,sum = 0.0;
223:   Vec            V_p;
224:   PetscReal      eta;

227:   VecGetOwnershipRange(ctx->V[0],&istart,NULL);
228:   VecGetOwnershipRange(S1[0],&p_start,&p_end);

230:   VecDuplicate(S1[0],&V_p);
231:   for (i=0;i<ctx->L;i++) {
232:     VecGetArray(ctx->V[i],&data);
233:     VecGetArray(V_p,&p_data);
234:     for (j=p_start;j<p_end;j++) p_data[j-p_start] = data[j-istart];
235:     VecRestoreArray(ctx->V[i],&data);
236:     VecRestoreArray(V_p,&p_data);
237:     VecDot(V_p,S1[i],&tmp);
238:     sum += tmp;
239:   }
240:   VecDestroy(&V_p);
241:   ctx->est_eig = PetscAbsScalar(ctx->radius*sum/(PetscReal)ctx->L);
242:   eta = PetscPowReal(10,-log10(eps->tol)/ctx->N);
243:   PetscInfo1(eps,"Estimation_#Eig %F\n",ctx->est_eig);
244:   *L_add = (PetscInt)ceil((ctx->est_eig*eta)/ctx->M) - ctx->L;
245:   if (*L_add < 0) *L_add = 0;
246:   if (*L_add>ctx->L_max-ctx->L) {
247:     PetscInfo(eps,"Number of eigenvalues around the contour path may be too large\n");
248:     *L_add = ctx->L_max-ctx->L;
249:   }
250:   return(0);
251: }

255: static PetscErrorCode SetAddVector(EPS eps,PetscInt Ladd_end)
256: {
258:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
259:   PetscInt       i,j,nlocal,Ladd_start=ctx->L;
260:   Vec            *newV;
261:   PetscScalar    *vdata;

264:   PetscMalloc(Ladd_end*sizeof(Vec*),&newV);
265:   for (i=0;i<ctx->L;i++) { newV[i] = ctx->V[i]; }
266:   PetscFree(ctx->V);
267:   ctx->V = newV;
268:   VecGetLocalSize(ctx->V[0],&nlocal);
269:   for (i=Ladd_start;i<Ladd_end;i++) {
270:     VecDuplicate(ctx->V[0],&ctx->V[i]);
271:     PetscLogObjectParent(eps,ctx->V[i]);
272:     CISSVecSetRandom(ctx->V[i],eps->rand);
273:     VecGetArray(ctx->V[i],&vdata);
274:     for (j=0;j<nlocal;j++) {
275:       vdata[j] = PetscRealPart(vdata[j]);
276:       if (PetscRealPart(vdata[j]) < 0.5) vdata[j] = -1.0;
277:       else vdata[j] = 1.0;
278:     }
279:     VecRestoreArray(ctx->V[i],&vdata);
280:   }
281:   return(0);
282: }

286: static PetscErrorCode SolveAddLinearSystem(EPS eps,PetscInt Ladd_start,PetscInt Ladd_end)
287: {
289:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
290:   PetscInt       i,j;

293:   for (i=0;i<ctx->num_solve_point;i++) {
294:     for (j=Ladd_start;j<Ladd_end;j++) {
295:       VecDestroy(&ctx->Y[i*ctx->L_max+j]);
296:       VecDuplicate(ctx->V[0],&ctx->Y[i*ctx->L_max+j]);
297:       PetscLogObjectParent(eps,ctx->Y[i*ctx->L_max+j]);
298:       KSPSolve(ctx->ksp[i],ctx->V[j],ctx->Y[i*ctx->L_max+j]);
299:     }
300:   }
301:   return(0);
302: }

306: static PetscErrorCode CalcMu(EPS eps,PetscScalar *Mu)
307: {
309:   PetscInt       i,j,k,s;
310:   PetscInt       rank_region,icolor,ikey;
311:   PetscScalar    *temp,*temp2,*ppk,alp;
312:   MPI_Comm       Row_Comm;
313:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

316:   MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank_region);
317:   icolor = rank_region % ctx->npart;
318:   ikey = rank_region / ctx->npart;
319:   MPI_Comm_split(PetscObjectComm((PetscObject)eps),icolor,ikey,&Row_Comm);

321:   PetscMalloc(ctx->num_solve_point*ctx->L*ctx->L*sizeof(PetscScalar),&temp);
322:   PetscMalloc(2*ctx->M*ctx->L*ctx->L*sizeof(PetscScalar),&temp2);
323:   PetscMalloc(ctx->num_solve_point*sizeof(PetscScalar),&ppk);
324:   for (i=0;i<2*ctx->M*ctx->L*ctx->L;i++) temp2[i] = 0;
325:   for (i=0; i<ctx->num_solve_point;i++) {
326:     for (j=0;j<ctx->L;j++) {
327:       VecMDot(ctx->Y[i*ctx->L_max+j],ctx->L,ctx->V,&temp[(j+i*ctx->L)*ctx->L]);
328:     }
329:   }

331:   for (i=0;i<ctx->num_solve_point;i++) ppk[i] = 1;
332:   for (k=0;k<2*ctx->M;k++) {
333:     for (j=0;j<ctx->L;j++) {
334:       for (i=0;i<ctx->num_solve_point;i++) {
335:           alp = ppk[i]*ctx->weight[ctx->solver_comm_id*ctx->num_solve_point+i]/(PetscReal)ctx->N;
336:         for (s=0;s<ctx->L;s++) {
337:           if (ctx->useconj) temp2[s+(j+k*ctx->L)*ctx->L] += PetscRealPart(alp*temp[s+(j+i*ctx->L)*ctx->L])*2;
338:           else temp2[s+(j+k*ctx->L)*ctx->L] += alp*temp[s+(j+i*ctx->L)*ctx->L];
339:         }
340:       }
341:     }
342:     for (i=0;i<ctx->num_solve_point;i++)
343:       ppk[i] *= ctx->pp[ctx->solver_comm_id*ctx->num_solve_point+i];
344:   }
345:   MPI_Allreduce(temp2,Mu,2*ctx->M*ctx->L*ctx->L,MPIU_SCALAR,MPIU_SUM,Row_Comm);

347:   PetscFree(ppk);
348:   PetscFree(temp);
349:   PetscFree(temp2);
350:   return(0);
351: }

355: static PetscErrorCode BlockHankel(EPS eps,PetscScalar *Mu,PetscInt s,Vec *H)
356: {
357:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
358:   PetscInt       i,j,k,L=ctx->L,M=ctx->M;
359:   PetscScalar    *H_data;

363:   for (k=0;k<L*M;k++) {
364:     VecGetArray(H[k],&H_data);
365:     for (j=0;j<M;j++)
366:       for (i=0;i<L;i++)
367:         H_data[j*L+i] = Mu[i+k*L+(j+s)*L*L];
368:     VecRestoreArray(H[k],&H_data);
369:   }
370:   return(0);
371: }

375: static PetscErrorCode SVD(EPS eps,Vec *Q,PetscInt *K,PetscBool isqr)
376: {
378:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
379:   PetscInt       i,j,ld,ml=ctx->L*ctx->M,n=eps->n;
380:   PetscScalar    *R,*w,*s;
381:   DS             ds;

384:   if (isqr) {
385:     PetscMalloc(ml*ml*sizeof(PetscScalar),&s);
386:     IPQRDecomposition(eps->ip,Q,0,ml,s,ml);
387:   }

389:   DSCreate(PETSC_COMM_WORLD,&ds);
390:   DSSetType(ds,DSSVD);
391:   DSSetFromOptions(ds);
392:   ld = ml;
393:   DSAllocate(ds,ld);
394:   DSSetDimensions(ds,PetscMin(n,ml),ml,0,0);
395:   DSGetArray(ds,DS_MAT_A,&R);

397:   if (isqr) {
398:     for (i=0;i<ml;i++)
399:       for (j=0;j<PetscMin(n,ml);j++)
400:         R[i*ld+j] = s[i*ml+j];
401:   } else {
402:     for (i=0;i<ml;i++) {
403:       VecGetArray(Q[i],&s);
404:       for (j=0;j<PetscMin(n,ml);j++) {
405:         R[i*ld+j] = s[j];
406:       }
407:       VecRestoreArray(Q[i],&s);
408:     }
409:   }
410:   DSRestoreArray(ds,DS_MAT_A,&R);
411:   if (isqr) { PetscFree(s); }
412:   DSSetState(ds,DS_STATE_RAW);
413:   PetscMalloc(ml*sizeof(PetscScalar),&w);
414:   DSSetEigenvalueComparison(ds,SlepcCompareLargestReal,NULL);
415:   DSSolve(ds,w,NULL);
416:   DSSort(ds,w,NULL,NULL,NULL,NULL);
417:   (*K) = 0;
418:   for (i=0;i<ml;i++) {
419:     ctx->sigma[i] = PetscRealPart(w[i]);
420:     if (ctx->sigma[i]/PetscMax(ctx->sigma[0],1)>ctx->delta) (*K)++;
421:   }
422:   PetscFree(w);
423:   DSDestroy(&ds);
424:   return(0);
425: }

429: static PetscErrorCode ProjectMatrix(Mat A,PetscInt nv,PetscInt ld,Vec *Q,PetscScalar *H,Vec w,PetscBool isherm)
430: {
432:   PetscInt       i,j;

435:   if (isherm) {
436:     for (j=0;j<nv;j++) {
437:       MatMult(A,Q[j],w);
438:       VecMDot(w,j+1,Q,H+j*ld);
439:       for (i=0;i<j;i++)
440:         H[j+i*ld] = PetscConj(H[i+j*ld]);
441:     }
442:   } else {
443:     for (j=0;j<nv;j++) {
444:       MatMult(A,Q[j],w);
445:       VecMDot(w,nv,Q,H+j*ld);
446:     }
447:   }
448:   return(0);
449: }

453: static PetscErrorCode isInsideGamma(EPS eps,PetscInt nv,PetscBool *fl)
454: {
455:   EPS_CISS    *ctx = (EPS_CISS*)eps->data;
456:   PetscInt    i;
457:   PetscScalar d;
458:   PetscReal   dx,dy;
459:   for (i=0;i<nv;i++) {
460:     d = (eps->eigr[i]-ctx->center)/ctx->radius;
461:     dx = PetscRealPart(d);
462:     dy = PetscImaginaryPart(d);
463:     if ((dx*dx+(dy*dy)/(ctx->vscale*ctx->vscale))<=1) fl[i] = PETSC_TRUE;
464:     else fl[i] = PETSC_FALSE;
465:   }
466:   return(0);
467: }

471: PetscErrorCode EPSSetUp_CISS(EPS eps)
472: {
474:   PetscInt       i;
475:   Vec            stemp;
476:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
477:   const char     *prefix;

480: #if !defined(PETSC_USE_COMPLEX)
481:   SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"CISS only works for complex scalars");
482: #endif
483:   eps->ncv = PetscMin(eps->n,ctx->L*ctx->M);
484:   if (!eps->mpd) eps->mpd = eps->ncv;
485:   if (!eps->which) eps->which = EPS_ALL;
486:   if (!eps->extraction) {
487:     EPSSetExtraction(eps,EPS_RITZ);
488:   } else if (eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
489:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");

491:   if (ctx->isreal && PetscImaginaryPart(ctx->center) == 0.0) ctx->useconj = PETSC_TRUE;
492:   else ctx->useconj = PETSC_FALSE;

494:   if (!ctx->delta) ctx->delta = PetscMin((eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL*1e-1:eps->tol*1e-1),1e-12);

496:   if (!ctx->vscale) {
497:     if (eps->ishermitian && (eps->ispositive || !eps->isgeneralized) && PetscImaginaryPart(ctx->center) == 0.0) ctx->vscale = 0.1;
498:     else ctx->vscale = 1.0;
499:   }

501:   /* create split comm */
502:   SetSolverComm(eps);

504:   EPSAllocateSolution(eps);
505:   PetscMalloc(ctx->N*sizeof(PetscScalar),&ctx->weight);
506:   PetscMalloc(ctx->N*sizeof(PetscScalar),&ctx->omega);
507:   PetscMalloc(ctx->N*sizeof(PetscScalar),&ctx->pp);
508:   PetscLogObjectMemory(eps,3*ctx->N*sizeof(PetscScalar));
509:   PetscMalloc(ctx->L*ctx->M*sizeof(PetscReal),&ctx->sigma);

511:   /* create a template vector for Vecs on solver communicator */
512:   VecCreateMPI(PetscObjectComm((PetscObject)eps),PETSC_DECIDE,eps->n,&stemp);
513:   VecDuplicateVecs(stemp,ctx->L,&ctx->V);
514:   PetscLogObjectParents(eps,ctx->L,ctx->V);
515:   VecDestroy(&stemp);

517:   PetscMalloc(ctx->num_solve_point*sizeof(KSP),&ctx->ksp);
518:   PetscLogObjectMemory(eps,ctx->num_solve_point*sizeof(KSP));
519:   for (i=0;i<ctx->num_solve_point;i++) {
520:     KSPCreate(PetscObjectComm((PetscObject)eps),&ctx->ksp[i]);
521:     PetscObjectIncrementTabLevel((PetscObject)ctx->ksp[i],(PetscObject)eps,1);
522:     PetscLogObjectParent(eps,ctx->ksp[i]);
523:     KSPAppendOptionsPrefix(ctx->ksp[i],"eps_ciss_");
524:     EPSGetOptionsPrefix(eps,&prefix);
525:     KSPAppendOptionsPrefix(ctx->ksp[i],prefix);
526:   }
527:   PetscMalloc(ctx->num_solve_point*ctx->L_max*sizeof(Vec),&ctx->Y);
528:   PetscMemzero(ctx->Y,ctx->num_solve_point*ctx->L_max*sizeof(Vec));
529:   PetscLogObjectMemory(eps,ctx->num_solve_point*ctx->L_max*sizeof(Vec));

531:   if (eps->isgeneralized) {
532:     if (eps->ishermitian && eps->ispositive) {
533:       DSSetType(eps->ds,DSGHEP);
534:     } else {
535:       DSSetType(eps->ds,DSGNHEP);
536:     }
537:   } else {
538:     if (eps->ishermitian) {
539:       DSSetType(eps->ds,DSHEP);
540:     } else {
541:       DSSetType(eps->ds,DSNHEP);
542:     }
543:   }
544:   DSAllocate(eps->ds,eps->ncv);
545:   EPSSetWorkVecs(eps,2);

547:   /* dispatch solve method */
548:   if (eps->leftvecs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Left vectors not supported in this solver");
549:   eps->ops->solve = EPSSolve_CISS;
550:   return(0);
551: }

555: PetscErrorCode EPSSolve_CISS(EPS eps)
556: {
558:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
559:   PetscInt       i,j,k,ld,nv,nmat,nvecs,L_add=0,inner,outer,L_base=ctx->L;
560:   PetscScalar    *H,*rr,*pX,*tdata,*vdata,*Mu;
561:   PetscReal      *tau,s1,s2,tau_max=0.0,*temp,error,max_error;
562:   PetscBool      *fl;
563:   Mat            A,B;
564:   Vec            w=eps->work[0],tempv=eps->work[1],*H0;

567:   DSGetLeadingDimension(eps->ds,&ld);
568:   STGetNumMatrices(eps->st,&nmat);
569:   STGetOperators(eps->st,0,&A);
570:   if (nmat>1) { STGetOperators(eps->st,1,&B); }
571:   else B = NULL;

573:   SetPathParameter(eps);
574:   for (i=0;i<ctx->L;i++) {
575:     CISSVecSetRandom(ctx->V[i],eps->rand);
576:   }
577:   SolveLinearSystem(eps);
578:   ConstructS(eps,1,&ctx->S);
579:   nvecs = ctx->L;
580:   EstimateNumberEigs(eps,ctx->S,&L_add);

582:   if (L_add>0) {
583:     PetscInfo2(eps,"Changing L %d -> %d by Estimate #Eig\n",ctx->L,ctx->L+L_add);
584:     SetAddVector(eps,ctx->L+L_add);
585:     SolveAddLinearSystem(eps,ctx->L,ctx->L+L_add);
586:     ctx->L += L_add;
587:     PetscFree(ctx->sigma);
588:     PetscMalloc(ctx->L*ctx->M*sizeof(PetscReal),&ctx->sigma);
589:   }

591:   for (i=0;i<ctx->refine_blocksize;i++) {
592:     PetscMalloc(ctx->L*ctx->L*ctx->M*2*sizeof(PetscScalar),&Mu);
593:     CalcMu(eps,Mu);
594:     VecDuplicateVecs(eps->t,ctx->L*ctx->M,&H0);
595:     BlockHankel(eps,Mu,0,H0);
596:     SVD(eps,H0,&nv,PETSC_FALSE);
597:     PetscFree(Mu);
598:     VecDestroyVecs(ctx->L*ctx->M,&H0);
599:     if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M) break;
600:     L_add = L_base;
601:     PetscInfo2(eps,"Changing L %d -> %d by SVD(H0)\n",ctx->L,ctx->L+L_add);
602:     SetAddVector(eps,ctx->L+L_add);
603:     SolveAddLinearSystem(eps,ctx->L,ctx->L+L_add);
604:     ctx->L += L_add;
605:     PetscFree(ctx->sigma);
606:     PetscMalloc(ctx->L*ctx->M*sizeof(PetscReal),&ctx->sigma);
607:   }

609:   if (ctx->L != L_base) {
610:     eps->ncv = PetscMin(eps->n,ctx->L*ctx->M);
611:     eps->mpd = eps->ncv;
612:     EPSAllocateSolution(eps);
613:     DSReset(eps->ds);
614:     DSSetEigenvalueComparison(eps->ds,eps->comparison,eps->comparisonctx);
615:     if (eps->isgeneralized) {
616:       if (eps->ishermitian && eps->ispositive) {
617:         DSSetType(eps->ds,DSGHEP);
618:       } else {
619:         DSSetType(eps->ds,DSGNHEP);
620:       }
621:     } else {
622:       if (eps->ishermitian) {
623:         DSSetType(eps->ds,DSHEP);
624:       } else {
625:         DSSetType(eps->ds,DSNHEP);
626:       }
627:     }
628:     DSAllocate(eps->ds,eps->ncv);
629:     DSGetLeadingDimension(eps->ds,&ld);
630:   }

632:   for (outer=0;outer<=ctx->refine_outer;outer++) {
633:     for (inner=0;inner<=ctx->refine_inner;inner++) {
634:       VecDestroyVecs(nvecs,&ctx->S);
635:       ConstructS(eps,ctx->M,&ctx->S);
636:       nvecs = ctx->M*ctx->L;
637:       SVD(eps,ctx->S,&nv,PETSC_TRUE);
638:       if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
639:         for (i=0;i<ctx->L;i++) {
640:           VecCopy(ctx->S[i],ctx->V[i]);
641:         }
642:         SolveAddLinearSystem(eps,0,ctx->L);
643:       } else break;
644:     }
645:     eps->nconv = 0;
646:     if (nv == 0) break;
647:     DSSetDimensions(eps->ds,nv,0,0,0);
648:     DSSetState(eps->ds,DS_STATE_RAW);

650:     DSGetArray(eps->ds,DS_MAT_A,&H);
651:     ProjectMatrix(A,nv,ld,ctx->S,H,w,eps->ishermitian);
652:     DSRestoreArray(eps->ds,DS_MAT_A,&H);
653: 
654:     if (nmat>1) {
655:       DSGetArray(eps->ds,DS_MAT_B,&H);
656:       ProjectMatrix(B,nv,ld,ctx->S,H,w,eps->ishermitian);
657:       DSRestoreArray(eps->ds,DS_MAT_B,&H);
658:     }
659: 
660:     DSSolve(eps->ds,eps->eigr,NULL);

662:     PetscMalloc(nv*sizeof(PetscReal),&tau);
663:     DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
664:     DSGetArray(eps->ds,DS_MAT_X,&pX);
665:     for (i=0;i<nv;i++) {
666:       s1 = 0;
667:       s2 = 0;
668:       for (j=0;j<nv;j++) {
669:         s1 += PetscAbsScalar(PetscPowScalar(pX[i*ld+j],2));
670:         s2 += PetscPowScalar(PetscAbsScalar(pX[i*ld+j]),2)/ctx->sigma[j];
671:       }
672:       tau[i] = s1/s2;
673:       tau_max = PetscMax(tau_max,tau[i]);
674:     }
675:     tau_max /= ctx->sigma[0];
676:     DSRestoreArray(eps->ds,DS_MAT_X,&pX);
677:     for (i=0;i<nv;i++) tau[i] /= tau_max;
678:     PetscMalloc(nv*sizeof(PetscBool),&fl);
679:     isInsideGamma(eps,nv,fl);
680:     PetscMalloc(nv*sizeof(PetscScalar),&rr);
681:     for (i=0;i<nv;i++) {
682:       if (fl[i] && tau[i]>=ctx->spurious_threshold*tau_max) {
683:         rr[i] = 1.0;
684:         eps->nconv++;
685:       } else rr[i] = 0.0;
686:     }

688:     PetscFree(tau);
689:     PetscFree(fl);
690:     DSSetEigenvalueComparison(eps->ds,SlepcCompareLargestMagnitude,NULL);
691:     DSSort(eps->ds,eps->eigr,NULL,rr,NULL,&eps->nconv);
692:     DSSetEigenvalueComparison(eps->ds,eps->comparison,eps->comparisonctx);
693:     PetscFree(rr);
694:     for (i=0;i<nv;i++) {
695:       VecCopy(ctx->S[i],eps->V[i]);
696:     }
697: 
698:     DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
699:     DSGetArray(eps->ds,DS_MAT_X,&pX);
700:     SlepcUpdateVectors(nv,ctx->S,0,eps->nconv,pX,ld,PETSC_FALSE);
701:     if (eps->ishermitian) {  /* compute eigenvectors */
702:       SlepcUpdateVectors(nv,eps->V,0,eps->nconv,pX,ld,PETSC_FALSE);
703:     }
704:     DSRestoreArray(eps->ds,DS_MAT_X,&pX);

706:     max_error = 0.0;
707:     for (i=0;i<eps->nconv;i++) {
708:       VecNormalize(eps->V[i],NULL);
709:       VecNormalize(ctx->S[i],NULL);
710:       EPSComputeRelativeError_Private(eps,eps->eigr[i],0,ctx->S[i],NULL,&error);
711:       max_error = PetscMax(max_error,error);
712:     }
713:     if (max_error <= eps->tol || outer == ctx->refine_outer) break;
714:     PetscMalloc(ctx->L*eps->nconv*sizeof(PetscReal),&temp);
715:     for (i=0;i<ctx->L*eps->nconv;i++) {
716:       PetscRandomGetValueReal(eps->rand,&temp[i]);
717:       temp[i] = 2*temp[i]-1;
718:     }
719: 
720:     for (k=0;k<ctx->L;k++) {
721:       VecGetArray(tempv,&tdata);
722:       for (j=0;j<eps->nconv;j++) {
723:         VecGetArray(ctx->S[j],&vdata);
724:         for (i=0;i<eps->n;i++) {
725:           if (j==0) tdata[i] = vdata[i]*temp[j+eps->nconv*k];
726:           else tdata[i] = tdata[i]+vdata[i]*temp[j+eps->nconv*k];
727:         }
728:         VecRestoreArray(ctx->S[j],&vdata);
729:       }
730:       VecRestoreArray(tempv,&tdata);
731:       VecCopy(tempv,ctx->V[k]);
732:     }
733: 
734:     PetscFree(temp);
735:     SolveAddLinearSystem(eps,0,ctx->L);
736:   }
737:   eps->reason = EPS_CONVERGED_TOL;
738:   return(0);
739: }

743: static PetscErrorCode EPSCISSSetRegion_CISS(EPS eps,PetscScalar center,PetscReal radius,PetscReal vscale)
744: {
745:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

748:   ctx->center = center;
749:   if (radius) {
750:     if (radius == PETSC_DEFAULT) {
751:       ctx->radius = 1.0;
752:     } else {
753:       if (radius<0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The radius argument must be > 0.0");
754:       ctx->radius = radius;
755:     }
756:   }
757:   if (vscale) {
758:     if (vscale<0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The vscale argument must be > 0.0");
759:     ctx->vscale = vscale;
760:   }
761:   return(0);
762: }

766: /*@
767:    EPSCISSSetRegion - Sets the parameters defining the region where eigenvalues
768:    must be computed.

770:    Logically Collective on EPS

772:    Input Parameters:
773: +  eps - the eigenproblem solver context
774: .  center - center of the region
775: .  radius - radius of the region
776: -  vscale - vertical scale of the region

778:    Options Database Keys:
779: +  -eps_ciss_center - Sets the center
780: .  -eps_ciss_radius - Sets the radius
781: -  -eps_ciss_vscale - Sets the vertical scale

783:    Level: advanced

785: .seealso: EPSCISSGetRegion()
786: @*/
787: PetscErrorCode EPSCISSSetRegion(EPS eps,PetscScalar center,PetscReal radius,PetscReal vscale)
788: {

796:   PetscTryMethod(eps,"EPSCISSSetRegion_C",(EPS,PetscScalar,PetscReal,PetscReal),(eps,center,radius,vscale));
797:   return(0);
798: }

802: static PetscErrorCode EPSCISSGetRegion_CISS(EPS eps,PetscScalar *center,PetscReal *radius,PetscReal *vscale)
803: {
804:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

807:   if (center) *center = ctx->center;
808:   if (radius) *radius = ctx->radius;
809:   if (vscale) *vscale = ctx->vscale;
810:   return(0);
811: }

815: /*@
816:    EPSCISSGetRegion - Gets the parameters that define the region where eigenvalues
817:    must be computed.

819:    Not Collective

821:    Input Parameter:
822: .  eps - the eigenproblem solver context

824:    Output Parameters:
825: +  center - center of the region
826: .  radius - radius of the region
827: -  vscale - vertical scale of the region

829:    Level: advanced

831: .seealso: EPSCISSSetRegion()
832: @*/
833: PetscErrorCode EPSCISSGetRegion(EPS eps,PetscScalar *center,PetscReal *radius,PetscReal *vscale)
834: {

839:   PetscTryMethod(eps,"EPSCISSGetRegion_C",(EPS,PetscScalar*,PetscReal*,PetscReal*),(eps,center,radius,vscale));
840:   return(0);
841: }

845: static PetscErrorCode EPSCISSSetSizes_CISS(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool isreal)
846: {
848:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

851:   if (ip) {
852:     if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) {
853:       if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
854:     } else {
855:       if (ip<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
856:       if (ip%2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
857:       if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
858:     }
859:   }
860:   if (bs) {
861:     if (bs == PETSC_DECIDE || bs == PETSC_DEFAULT) {
862:       ctx->L = 16;
863:     } else {
864:       if (bs<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
865:       if (bs>ctx->L_max) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be less than or equal to the maximum number of block size");
866:       ctx->L = bs;
867:     }
868:   }
869:   if (ms) {
870:     if (ms == PETSC_DECIDE || ms == PETSC_DEFAULT) {
871:       ctx->M = ctx->N/4;
872:     } else {
873:       if (ms<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
874:       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");
875:       ctx->M = ms;
876:     }
877:   }
878:   if (npart) {
879:     if (npart == PETSC_DECIDE || npart == PETSC_DEFAULT) {
880:       ctx->npart = 1;
881:     } else {
882:       if (npart<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The npart argument must be > 0");
883:       ctx->npart = npart;
884:     }
885:   }
886:   if (bsmax) {
887:     if (bsmax == PETSC_DECIDE || bsmax == PETSC_DEFAULT) {
888:       ctx->L = 256;
889:     } else {
890:       if (bsmax<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
891:       if (bsmax<ctx->L) ctx->L_max = ctx->L;
892:       else ctx->L_max = bsmax;
893:     }
894:   }
895:   ctx->isreal = isreal;
896:   EPSReset(eps);   /* clean allocated arrays and force new setup */
897:   return(0);
898: }

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

905:    Logically Collective on EPS

907:    Input Parameters:
908: +  eps   - the eigenproblem solver context
909: .  ip    - number of integration points
910: .  bs    - block size
911: .  ms    - moment size
912: .  npart - number of partitions when splitting the communicator
913: .  bsmax - max block size
914: -  isreal - A and B are real

916:    Options Database Keys:
917: +  -eps_ciss_integration_points - Sets the number of integration points
918: .  -eps_ciss_blocksize - Sets the block size
919: .  -eps_ciss_moments - Sets the moment size
920: .  -eps_ciss_partitions - Sets the number of partitions
921: .  -eps_ciss_maxblocksize - Sets the maximum block size
922: -  -eps_ciss_realmats - A and B are real

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

929:    Level: advanced

931: .seealso: EPSCISSGetSizes()
932: @*/
933: PetscErrorCode EPSCISSSetSizes(EPS eps,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool isreal)
934: {

945:   PetscTryMethod(eps,"EPSCISSSetSizes_C",(EPS,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(eps,ip,bs,ms,npart,bsmax,isreal));
946:   return(0);
947: }

951: static PetscErrorCode EPSCISSGetSizes_CISS(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *isreal)
952: {
953:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

956:   if (ip) *ip = ctx->N;
957:   if (bs) *bs = ctx->L;
958:   if (ms) *ms = ctx->M;
959:   if (npart) *npart = ctx->npart;
960:   if (bsmax) *bsmax = ctx->L_max;
961:   if (isreal) *isreal = ctx->isreal;
962:   return(0);
963: }

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

970:    Not Collective

972:    Input Parameter:
973: .  eps - the eigenproblem solver context

975:    Output Parameters:
976: +  ip    - number of integration points
977: .  bs    - block size
978: .  ms    - moment size
979: .  npart - number of partitions when splitting the communicator
980: .  bsmax - max block size
981: -  isreal - A and B are real

983:    Level: advanced

985: .seealso: EPSCISSSetSizes()
986: @*/
987: PetscErrorCode EPSCISSGetSizes(EPS eps,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *isreal)
988: {

993:   PetscTryMethod(eps,"EPSCISSGetSizes_C",(EPS,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(eps,ip,bs,ms,npart,bsmax,isreal));
994:   return(0);
995: }

999: static PetscErrorCode EPSCISSSetThreshold_CISS(EPS eps,PetscReal delta,PetscReal spur)
1000: {
1001:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1004:   if (delta) {
1005:     if (delta == PETSC_DEFAULT) {
1006:       ctx->delta = 1e-12;
1007:     } else {
1008:       if (delta<=0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
1009:       ctx->delta = delta;
1010:     }
1011:   }
1012:   if (spur) {
1013:     if (spur == PETSC_DEFAULT) {
1014:       ctx->spurious_threshold = 1e-4;
1015:     } else {
1016:       if (spur<=0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
1017:       ctx->spurious_threshold = spur;
1018:     }
1019:   }
1020:   return(0);
1021: }

1025: /*@
1026:    EPSCISSSetThreshold - Sets the values of various threshold parameters in
1027:    the CISS solver.

1029:    Logically Collective on EPS

1031:    Input Parameters:
1032: +  eps   - the eigenproblem solver context
1033: .  delta - threshold for numerical rank
1034: -  spur  - spurious threshold (to discard spurious eigenpairs)

1036:    Options Database Keys:
1037: +  -eps_ciss_delta - Sets the delta
1038: -  -eps_ciss_spurious_threshold - Sets the spurious threshold

1040:    Level: advanced

1042: .seealso: EPSCISSGetThreshold()
1043: @*/
1044: PetscErrorCode EPSCISSSetThreshold(EPS eps,PetscReal delta,PetscReal spur)
1045: {

1052:   PetscTryMethod(eps,"EPSCISSSetThreshold_C",(EPS,PetscReal,PetscReal),(eps,delta,spur));
1053:   return(0);
1054: }

1058: static PetscErrorCode EPSCISSGetThreshold_CISS(EPS eps,PetscReal *delta,PetscReal *spur)
1059: {
1060:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1063:   if (delta) *delta = ctx->delta;
1064:   if (spur)  *spur = ctx->spurious_threshold;
1065:   return(0);
1066: }

1070: /*@
1071:    EPSCISSGetThreshold - Gets the values of various threshold parameters
1072:    in the CISS solver.

1074:    Not Collective

1076:    Input Parameter:
1077: .  eps - the eigenproblem solver context

1079:    Output Parameters:
1080: +  delta - threshold for numerical rank
1081: -  spur  - spurious threshold (to discard spurious eigenpairs)

1083:    Level: advanced

1085: .seealso: EPSCISSSetThreshold()
1086: @*/
1087: PetscErrorCode EPSCISSGetThreshold(EPS eps,PetscReal *delta,PetscReal *spur)
1088: {

1093:   PetscTryMethod(eps,"EPSCISSGetThreshold_C",(EPS,PetscReal*,PetscReal*),(eps,delta,spur));
1094:   return(0);
1095: }

1099: static PetscErrorCode EPSCISSSetRefinement_CISS(EPS eps,PetscInt inner,PetscInt outer,PetscInt blsize)
1100: {
1101:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1104:   if (inner == PETSC_DEFAULT) {
1105:     ctx->refine_inner = 0;
1106:   } else {
1107:     if (inner<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
1108:     ctx->refine_inner = inner;
1109:   }
1110:   if (outer == PETSC_DEFAULT) {
1111:     ctx->refine_outer = 0;
1112:   } else {
1113:     if (outer<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine outer argument must be >= 0");
1114:     ctx->refine_outer = outer;
1115:   }
1116:   if (blsize == PETSC_DEFAULT) {
1117:     ctx->refine_blocksize = 0;
1118:   } else {
1119:     if (blsize<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
1120:     ctx->refine_blocksize = blsize;
1121:   }
1122:   return(0);
1123: }

1127: /*@
1128:    EPSCISSSetRefinement - Sets the values of various refinement parameters
1129:    in the CISS solver.

1131:    Logically Collective on EPS

1133:    Input Parameters:
1134: +  eps    - the eigenproblem solver context
1135: .  inner  - number of iterative refinement iterations (inner loop)
1136: .  outer  - number of iterative refinement iterations (outer loop)
1137: -  blsize - number of iterative refinement iterations (blocksize loop)

1139:    Options Database Keys:
1140: +  -eps_ciss_refine_inner - Sets number of inner iterations
1141: .  -eps_ciss_refine_outer - Sets number of outer iterations
1142: -  -eps_ciss_refine_blocksize - Sets number of blocksize iterations

1144:    Level: advanced

1146: .seealso: EPSCISSGetRefinement()
1147: @*/
1148: PetscErrorCode EPSCISSSetRefinement(EPS eps,PetscInt inner,PetscInt outer,PetscInt blsize)
1149: {

1157:   PetscTryMethod(eps,"EPSCISSSetRefinement_C",(EPS,PetscInt,PetscInt,PetscInt),(eps,inner,outer,blsize));
1158:   return(0);
1159: }

1163: static PetscErrorCode EPSCISSGetRefinement_CISS(EPS eps,PetscInt *inner,PetscInt *outer,PetscInt *blsize)
1164: {
1165:   EPS_CISS *ctx = (EPS_CISS*)eps->data;

1168:   if (inner)  *inner = ctx->refine_inner;
1169:   if (outer)  *outer = ctx->refine_outer;
1170:   if (blsize) *blsize = ctx->refine_blocksize;
1171:   return(0);
1172: }

1176: /*@
1177:    EPSCISSGetRefinement - Gets the values of various refinement parameters
1178:    in the CISS solver.

1180:    Not Collective

1182:    Input Parameter:
1183: .  eps - the eigenproblem solver context

1185:    Output Parameters:
1186: +  inner  - number of iterative refinement iterations (inner loop)
1187: .  outer  - number of iterative refinement iterations (outer loop)
1188: -  blsize - number of iterative refinement iterations (blocksize loop)

1190:    Level: advanced

1192: .seealso: EPSCISSSetRefinement()
1193: @*/
1194: PetscErrorCode EPSCISSGetRefinement(EPS eps, PetscInt *inner, PetscInt *outer,PetscInt *blsize)
1195: {

1200:   PetscTryMethod(eps,"EPSCISSGetRefinement_C",(EPS,PetscInt*,PetscInt*,PetscInt*),(eps,inner,outer,blsize));
1201:   return(0);
1202: }

1206: PetscErrorCode EPSReset_CISS(EPS eps)
1207: {
1209:   PetscInt       i;
1210:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

1213:   PetscFree(ctx->weight);
1214:   PetscFree(ctx->omega);
1215:   PetscFree(ctx->pp);
1216:   VecDestroyVecs(ctx->L,&ctx->V);
1217:   for (i=0;i<ctx->num_solve_point;i++) {
1218:     KSPDestroy(&ctx->ksp[i]);
1219:   }
1220:   PetscFree(ctx->ksp);
1221:   PetscFree(ctx->sigma);
1222:   for (i=0;i<ctx->num_solve_point*ctx->L_max;i++) {
1223:     VecDestroy(&ctx->Y[i]);
1224:   }
1225:   PetscFree(ctx->Y);
1226:   VecDestroyVecs(ctx->M*ctx->L,&ctx->S);
1227:   EPSReset_Default(eps);
1228:   return(0);
1229: }

1233: PetscErrorCode EPSSetFromOptions_CISS(EPS eps)
1234: {
1236:   PetscScalar    s;
1237:   PetscReal      r1,r2,r3,r4;
1238:   PetscInt       i1=0,i2=0,i3=0,i4=0,i5=0,i6=0,i7=0,i8=0;
1239:   PetscBool      b1=PETSC_FALSE;

1242:   PetscOptionsHead("EPS CISS Options");
1243:   EPSCISSGetRegion(eps,&s,&r1,&r2);
1244:   PetscOptionsReal("-eps_ciss_radius","CISS radius of region","EPSCISSSetRegion",r1,&r1,NULL);
1245:   PetscOptionsScalar("-eps_ciss_center","CISS center of region","EPSCISSSetRegion",s,&s,NULL);
1246:   PetscOptionsReal("-eps_ciss_vscale","CISS vertical scale of region","EPSCISSSetRegion",r2,&r2,NULL);
1247:   EPSCISSSetRegion(eps,s,r1,r2);

1249:   PetscOptionsInt("-eps_ciss_integration_points","CISS number of integration points","EPSCISSSetSizes",i1,&i1,NULL);
1250:   PetscOptionsInt("-eps_ciss_blocksize","CISS block size","EPSCISSSetSizes",i2,&i2,NULL);
1251:   PetscOptionsInt("-eps_ciss_moments","CISS moment size","EPSCISSSetSizes",i3,&i3,NULL);
1252:   PetscOptionsInt("-eps_ciss_partitions","CISS number of partitions","EPSCISSSetSizes",i4,&i4,NULL);
1253:   PetscOptionsInt("-eps_ciss_maxblocksize","CISS maximum block size","EPSCISSSetSizes",i5,&i5,NULL);
1254:   PetscOptionsBool("-eps_ciss_realmats","CISS A and B are real","EPSCISSSetSizes",b1,&b1,NULL);
1255:   EPSCISSSetSizes(eps,i1,i2,i3,i4,i5,b1);

1257:   EPSCISSGetThreshold(eps,&r3,&r4);
1258:   PetscOptionsReal("-eps_ciss_delta","CISS threshold for numerical rank","EPSCISSSetThreshold",r3,&r3,NULL);
1259:   PetscOptionsReal("-eps_ciss_spurious_threshold","CISS threshold for the spurious eigenpairs","EPSCISSSetThreshold",r4,&r4,NULL);
1260:   EPSCISSSetThreshold(eps,r3,r4);

1262:   EPSCISSGetRefinement(eps,&i6,&i7,&i8);
1263:   PetscOptionsInt("-eps_ciss_refine_inner","CISS number of inner iterative refinement iterations","EPSCISSSetRefinement",i6,&i6,NULL);
1264:   PetscOptionsInt("-eps_ciss_refine_outer","CISS number of outer iterative refinement iterations","EPSCISSSetRefinement",i7,&i7,NULL);
1265:   PetscOptionsInt("-eps_ciss_refine_blocksize","CISS number of blocksize iterative refinement iterations","EPSCISSSetRefinement",i8,&i8,NULL);
1266:   EPSCISSSetRefinement(eps,i6,i7,i8);

1268:   PetscOptionsTail();
1269:   return(0);
1270: }

1274: PetscErrorCode EPSDestroy_CISS(EPS eps)
1275: {

1279:   PetscFree(eps->data);
1280:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRegion_C",NULL);
1281:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRegion_C",NULL);
1282:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",NULL);
1283:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",NULL);
1284:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",NULL);
1285:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",NULL);
1286:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",NULL);
1287:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",NULL);
1288:   return(0);
1289: }

1293: PetscErrorCode EPSView_CISS(EPS eps,PetscViewer viewer)
1294: {
1296:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;
1297:   PetscBool      isascii;
1298:   char           str[50];

1301:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1302:   if (isascii) {
1303:     SlepcSNPrintfScalar(str,50,ctx->center,PETSC_FALSE);
1304:     PetscViewerASCIIPrintf(viewer,"  CISS: region { center: %s, radius: %G, vscale: %G }\n",str,ctx->radius,ctx->vscale);
1305:     PetscViewerASCIIPrintf(viewer,"  CISS: sizes { integration points: %D, block size: %D, moment size: %D, partitions: %D, maximum block size: %D }\n",ctx->N,ctx->L,ctx->M,ctx->npart,ctx->L_max);
1306:     if (ctx->isreal) {
1307:       PetscViewerASCIIPrintf(viewer,"  CISS: exploiting symmetry of integration points\n");
1308:     }
1309:     PetscViewerASCIIPrintf(viewer,"  CISS: threshold { delta: %G, spurious threshold: %G }\n",ctx->delta,ctx->spurious_threshold);
1310:     PetscViewerASCIIPrintf(viewer,"  CISS: iterative refinement  { inner: %D, outer: %D, blocksize: %D }\n",ctx->refine_inner,ctx->refine_outer, ctx->refine_blocksize);
1311:     PetscViewerASCIIPushTab(viewer);
1312:     KSPView(ctx->ksp[0],viewer);
1313:     PetscViewerASCIIPopTab(viewer);

1315:   }
1316:   return(0);
1317: }

1321: PETSC_EXTERN PetscErrorCode EPSCreate_CISS(EPS eps)
1322: {
1324:   EPS_CISS       *ctx = (EPS_CISS*)eps->data;

1327:   PetscNewLog(eps,EPS_CISS,&ctx);
1328:   eps->data = ctx;
1329:   eps->ops->setup          = EPSSetUp_CISS;
1330:   eps->ops->setfromoptions = EPSSetFromOptions_CISS;
1331:   eps->ops->destroy        = EPSDestroy_CISS;
1332:   eps->ops->reset          = EPSReset_CISS;
1333:   eps->ops->view           = EPSView_CISS;
1334:   eps->ops->backtransform  = PETSC_NULL;
1335:   eps->ops->computevectors = EPSComputeVectors_Schur;
1336:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRegion_C",EPSCISSSetRegion_CISS);
1337:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRegion_C",EPSCISSGetRegion_CISS);
1338:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetSizes_C",EPSCISSSetSizes_CISS);
1339:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetSizes_C",EPSCISSGetSizes_CISS);
1340:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetThreshold_C",EPSCISSSetThreshold_CISS);
1341:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetThreshold_C",EPSCISSGetThreshold_CISS);
1342:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSSetRefinement_C",EPSCISSSetRefinement_CISS);
1343:   PetscObjectComposeFunction((PetscObject)eps,"EPSCISSGetRefinement_C",EPSCISSGetRefinement_CISS);
1344:   /* set default values of parameters */
1345:   ctx->center  = 0.0;
1346:   ctx->radius  = 1.0;
1347:   ctx->vscale  = 0.0;
1348:   ctx->N       = 32;
1349:   ctx->L       = 16;
1350:   ctx->M       = ctx->N/4;
1351:   ctx->delta   = 0;
1352:   ctx->npart   = 1;
1353:   ctx->L_max   = 128;
1354:   ctx->spurious_threshold = 1e-4;
1355:   ctx->isreal  = PETSC_FALSE;
1356:   ctx->refine_outer = 1;
1357:   ctx->refine_inner = 1;
1358:   ctx->refine_blocksize = 1;
1359:   return(0);
1360: }