Actual source code: cheby.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2:  #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h>

  4: static PetscErrorCode KSPReset_Chebyshev(KSP ksp)
  5: {
  6:   KSP_Chebyshev  *cheb = (KSP_Chebyshev*)ksp->data;

 10:   KSPReset(cheb->kspest);
 11:   return(0);
 12: }

 14: static PetscErrorCode KSPSetUp_Chebyshev(KSP ksp)
 15: {
 16:   KSP_Chebyshev  *cheb = (KSP_Chebyshev*)ksp->data;

 20:   KSPSetWorkVecs(ksp,3);
 21:   if ((cheb->emin == 0. || cheb->emax == 0.) && !cheb->kspest) { /* We need to estimate eigenvalues */
 22:     KSPChebyshevEstEigSet(ksp,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);
 23:   }
 24:   return(0);
 25: }

 27: static PetscErrorCode KSPChebyshevSetEigenvalues_Chebyshev(KSP ksp,PetscReal emax,PetscReal emin)
 28: {
 29:   KSP_Chebyshev  *chebyshevP = (KSP_Chebyshev*)ksp->data;

 33:   if (emax <= emin) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin);
 34:   if (emax*emin <= 0.0) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin);
 35:   chebyshevP->emax = emax;
 36:   chebyshevP->emin = emin;

 38:   KSPChebyshevEstEigSet(ksp,0.,0.,0.,0.); /* Destroy any estimation setup */
 39:   return(0);
 40: }

 42: static PetscErrorCode KSPChebyshevEstEigSet_Chebyshev(KSP ksp,PetscReal a,PetscReal b,PetscReal c,PetscReal d)
 43: {
 44:   KSP_Chebyshev  *cheb = (KSP_Chebyshev*)ksp->data;

 48:   if (a != 0.0 || b != 0.0 || c != 0.0 || d != 0.0) {
 49:     if (!cheb->kspest) { /* should this block of code be moved to KSPSetUp_Chebyshev()? */
 50:       KSPCreate(PetscObjectComm((PetscObject)ksp),&cheb->kspest);
 51:       PetscObjectIncrementTabLevel((PetscObject)cheb->kspest,(PetscObject)ksp,1);
 52:       KSPSetPC(cheb->kspest,ksp->pc);
 53:       /* use PetscObjectSet/AppendOptionsPrefix() instead of KSPSet/AppendOptionsPrefix() so that the PC prefix is not changed */
 54:       PetscObjectSetOptionsPrefix((PetscObject)cheb->kspest,((PetscObject)ksp)->prefix);
 55:       PetscObjectAppendOptionsPrefix((PetscObject)cheb->kspest,"esteig_");
 56:       KSPSetSkipPCSetFromOptions(cheb->kspest,PETSC_TRUE);

 58:       KSPSetComputeEigenvalues(cheb->kspest,PETSC_TRUE);

 60:       /* We cannot turn off convergence testing because GMRES will break down if you attempt to keep iterating after a zero norm is obtained */
 61:       KSPSetTolerances(cheb->kspest,1.e-12,PETSC_DEFAULT,PETSC_DEFAULT,cheb->eststeps);
 62:     }
 63:     if (a >= 0) cheb->tform[0] = a;
 64:     if (b >= 0) cheb->tform[1] = b;
 65:     if (c >= 0) cheb->tform[2] = c;
 66:     if (d >= 0) cheb->tform[3] = d;
 67:     cheb->amatid    = 0;
 68:     cheb->pmatid    = 0;
 69:     cheb->amatstate = -1;
 70:     cheb->pmatstate = -1;
 71:   } else {
 72:     KSPDestroy(&cheb->kspest);
 73:   }
 74:   return(0);
 75: }

 77: static PetscErrorCode KSPChebyshevEstEigSetUseNoisy_Chebyshev(KSP ksp,PetscBool use)
 78: {
 79:   KSP_Chebyshev  *cheb = (KSP_Chebyshev*)ksp->data;

 82:   cheb->usenoisy = use;
 83:   return(0);
 84: }

 86: /*@
 87:    KSPChebyshevSetEigenvalues - Sets estimates for the extreme eigenvalues
 88:    of the preconditioned problem.

 90:    Logically Collective on KSP

 92:    Input Parameters:
 93: +  ksp - the Krylov space context
 94: -  emax, emin - the eigenvalue estimates

 96:   Options Database:
 97: .  -ksp_chebyshev_eigenvalues emin,emax

 99:    Note: Call KSPChebyshevEstEigSet() or use the option -ksp_chebyshev_esteig a,b,c,d to have the KSP
100:          estimate the eigenvalues and use these estimated values automatically

102:    Level: intermediate

104: .keywords: KSP, Chebyshev, set, eigenvalues
105: @*/
106: PetscErrorCode  KSPChebyshevSetEigenvalues(KSP ksp,PetscReal emax,PetscReal emin)
107: {

114:   PetscTryMethod(ksp,"KSPChebyshevSetEigenvalues_C",(KSP,PetscReal,PetscReal),(ksp,emax,emin));
115:   return(0);
116: }

118: /*@
119:    KSPChebyshevEstEigSet - Automatically estimate the eigenvalues to use for Chebyshev

121:    Logically Collective on KSP

123:    Input Parameters:
124: +  ksp - the Krylov space context
125: .  a - multiple of min eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE)
126: .  b - multiple of max eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE)
127: .  c - multiple of min eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE)
128: -  d - multiple of max eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE)

130:   Options Database:
131: .  -ksp_chebyshev_esteig a,b,c,d

133:    Notes:
134:    The Chebyshev bounds are set using
135: .vb
136:    minbound = a*minest + b*maxest
137:    maxbound = c*minest + d*maxest
138: .ve
139:    The default configuration targets the upper part of the spectrum for use as a multigrid smoother, so only the maximum eigenvalue estimate is used.
140:    The minimum eigenvalue estimate obtained by Krylov iteration is typically not accurate until the method has converged.

142:    If 0.0 is passed for all transform arguments (a,b,c,d), eigenvalue estimation is disabled.

144:    The default transform is (0,0.1; 0,1.1) which targets the "upper" part of the spectrum, as desirable for use with multigrid.

146:    The eigenvalues are estimated using the Lanczo (KSPCG) or Arnoldi (KSPGMRES) process using a noisy right hand side vector.

148:    Level: intermediate

150: .keywords: KSP, Chebyshev, set, eigenvalues, PCMG
151: @*/
152: PetscErrorCode KSPChebyshevEstEigSet(KSP ksp,PetscReal a,PetscReal b,PetscReal c,PetscReal d)
153: {

162:   PetscTryMethod(ksp,"KSPChebyshevEstEigSet_C",(KSP,PetscReal,PetscReal,PetscReal,PetscReal),(ksp,a,b,c,d));
163:   return(0);
164: }

166: /*@
167:    KSPChebyshevEstEigSetUseNoisy - use a noisy right hand side in order to do the estimate instead of the given right hand side

169:    Logically Collective

171:    Input Arguments:
172: +  ksp - linear solver context
173: -  use - PETSC_TRUE to use noisy

175:    Options Database:
176: +  -ksp_chebyshev_esteig_noisy <true,false>

178:   Notes:
179:     This alledgely works better for multigrid smoothers

181:   Level: intermediate

183: .seealso: KSPChebyshevEstEigSet()
184: @*/
185: PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP ksp,PetscBool use)
186: {

191:   PetscTryMethod(ksp,"KSPChebyshevEstEigSetUseNoisy_C",(KSP,PetscBool),(ksp,use));
192:   return(0);
193: }

195: /*@
196:   KSPChebyshevEstEigGetKSP - Get the Krylov method context used to estimate eigenvalues for the Chebyshev method.  If
197:   a Krylov method is not being used for this purpose, NULL is returned.  The reference count of the returned KSP is
198:   not incremented: it should not be destroyed by the user.

200:   Input Parameters:
201: . ksp - the Krylov space context

203:   Output Parameters:
204: . kspest the eigenvalue estimation Krylov space context

206:   Level: intermediate

208: .seealso: KSPChebyshevEstEigSet()
209: @*/
210: PetscErrorCode KSPChebyshevEstEigGetKSP(KSP ksp, KSP *kspest)
211: {

216:   *kspest = NULL;
217:   PetscTryMethod(ksp,"KSPChebyshevEstEigGetKSP_C",(KSP,KSP*),(ksp,kspest));
218:   return(0);
219: }

221: static PetscErrorCode KSPChebyshevEstEigGetKSP_Chebyshev(KSP ksp, KSP *kspest)
222: {
223:   KSP_Chebyshev *cheb = (KSP_Chebyshev*)ksp->data;

226:   *kspest = cheb->kspest;
227:   return(0);
228: }

230: static PetscErrorCode KSPSetFromOptions_Chebyshev(PetscOptionItems *PetscOptionsObject,KSP ksp)
231: {
232:   KSP_Chebyshev  *cheb = (KSP_Chebyshev*)ksp->data;
234:   PetscInt       neigarg = 2, nestarg = 4;
235:   PetscReal      eminmax[2] = {0., 0.};
236:   PetscReal      tform[4] = {PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE};
237:   PetscBool      flgeig, flgest;

240:   PetscOptionsHead(PetscOptionsObject,"KSP Chebyshev Options");
241:   PetscOptionsInt("-ksp_chebyshev_esteig_steps","Number of est steps in Chebyshev","",cheb->eststeps,&cheb->eststeps,NULL);
242:   PetscOptionsRealArray("-ksp_chebyshev_eigenvalues","extreme eigenvalues","KSPChebyshevSetEigenvalues",eminmax,&neigarg,&flgeig);
243:   if (flgeig) {
244:     if (neigarg != 2) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"-ksp_chebyshev_eigenvalues: must specify 2 parameters, min and max eigenvalues");
245:     KSPChebyshevSetEigenvalues(ksp, eminmax[1], eminmax[0]);
246:   }
247:   PetscOptionsRealArray("-ksp_chebyshev_esteig","estimate eigenvalues using a Krylov method, then use this transform for Chebyshev eigenvalue bounds","KSPChebyshevEstEigSet",tform,&nestarg,&flgest);
248:   if (flgest) {
249:     switch (nestarg) {
250:     case 0:
251:       KSPChebyshevEstEigSet(ksp,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);
252:       break;
253:     case 2:                     /* Base everything on the max eigenvalues */
254:       KSPChebyshevEstEigSet(ksp,PETSC_DECIDE,tform[0],PETSC_DECIDE,tform[1]);
255:       break;
256:     case 4:                     /* Use the full 2x2 linear transformation */
257:       KSPChebyshevEstEigSet(ksp,tform[0],tform[1],tform[2],tform[3]);
258:       break;
259:     default: SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"Must specify either 0, 2, or 4 parameters for eigenvalue estimation");
260:     }
261:   }

263:   /* We need to estimate eigenvalues; need to set this here so that KSPSetFromOptions() is called on the estimator */
264:   if ((cheb->emin == 0. || cheb->emax == 0.) && !cheb->kspest) {
265:    KSPChebyshevEstEigSet(ksp,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);
266:   }

268:   if (cheb->kspest) {
269:     PetscOptionsBool("-ksp_chebyshev_esteig_noisy","Use noisy right hand side for estimate","KSPChebyshevEstEigSetUseNoisy",cheb->usenoisy,&cheb->usenoisy,NULL);
270:     KSPSetFromOptions(cheb->kspest);
271:   }
272:   PetscOptionsTail();
273:   return(0);
274: }

276: /*
277:  * Must be passed a KSP solver that has "converged", with KSPSetComputeEigenvalues() called before the solve
278:  */
279: static PetscErrorCode KSPChebyshevComputeExtremeEigenvalues_Private(KSP kspest,PetscReal *emin,PetscReal *emax)
280: {
282:   PetscInt       n,neig;
283:   PetscReal      *re,*im,min,max;

286:   KSPGetIterationNumber(kspest,&n);
287:   PetscMalloc2(n,&re,n,&im);
288:   KSPComputeEigenvalues(kspest,n,re,im,&neig);
289:   min  = PETSC_MAX_REAL;
290:   max  = PETSC_MIN_REAL;
291:   for (n=0; n<neig; n++) {
292:     min = PetscMin(min,re[n]);
293:     max = PetscMax(max,re[n]);
294:   }
295:   PetscFree2(re,im);
296:   *emax = max;
297:   *emin = min;
298:   return(0);
299: }

301: PETSC_STATIC_INLINE PetscScalar chebyhash(PetscInt xx)
302: {
303:   unsigned int x = xx;
304:   x = ((x >> 16) ^ x) * 0x45d9f3b;
305:   x = ((x >> 16) ^ x) * 0x45d9f3b;
306:   x = ((x >> 16) ^ x);
307:   return (PetscScalar)((PetscInt64)x-2147483648)*5.e-10; /* center around zero, scaled about -1. to 1.*/
308: }

310: static PetscErrorCode KSPSolve_Chebyshev(KSP ksp)
311: {
312:   KSP_Chebyshev  *cheb = (KSP_Chebyshev*)ksp->data;
314:   PetscInt       k,kp1,km1,maxit,ktmp,i;
315:   PetscScalar    alpha,omegaprod,mu,omega,Gamma,c[3],scale;
316:   PetscReal      rnorm = 0.0;
317:   Vec            sol_orig,b,p[3],r;
318:   Mat            Amat,Pmat;
319:   PetscBool      diagonalscale;

322:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
323:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

325:   PCGetOperators(ksp->pc,&Amat,&Pmat);
326:   if (cheb->kspest) {
327:     PetscObjectId    amatid,    pmatid;
328:     PetscObjectState amatstate, pmatstate;

330:     PetscObjectGetId((PetscObject)Amat,&amatid);
331:     PetscObjectGetId((PetscObject)Pmat,&pmatid);
332:     PetscObjectStateGet((PetscObject)Amat,&amatstate);
333:     PetscObjectStateGet((PetscObject)Pmat,&pmatstate);
334:     if (amatid != cheb->amatid || pmatid != cheb->pmatid || amatstate != cheb->amatstate || pmatstate != cheb->pmatstate) {
335:       PetscReal          max=0.0,min=0.0;
336:       Vec                B;
337:       KSPConvergedReason reason;

339:       if (cheb->usenoisy) {
340:         B  = ksp->work[1];
341:         {
343:           PetscInt       n,i,istart;
344:           PetscScalar    *xx;
345:           VecGetOwnershipRange(B,&istart,NULL);
346:           VecGetLocalSize(B,&n);
347:           VecGetArray(B,&xx);
348:           for (i=0; i<n; i++) {
349:             PetscScalar v = chebyhash(i+istart);
350:             xx[i] = v;
351:           }
352:           VecRestoreArray(B,&xx);
353:         }
354:       } else {
355:         PC        pc;
356:         PetscBool change;

358:         KSPGetPC(cheb->kspest,&pc);
359:         PCPreSolveChangeRHS(pc,&change);
360:         if (change) {
361:           B = ksp->work[1];
362:           VecCopy(ksp->vec_rhs,B);
363:         } else {
364:           B = ksp->vec_rhs;
365:         }
366:       }
367:       KSPSolve(cheb->kspest,B,ksp->work[0]);

369:       KSPGetConvergedReason(cheb->kspest,&reason);
370:       if (reason < 0) {
371:         if (reason == KSP_DIVERGED_ITS) {
372:           PetscInfo(ksp,"Eigen estimator ran for prescribed number of iterations\n");
373:         } else {
374:           PetscInt its;
375:           KSPGetIterationNumber(cheb->kspest,&its);
376:           if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
377:             PetscInfo(ksp,"Eigen estimator KSP_DIVERGED_PCSETUP_FAILED\n");
378:             ksp->reason = KSP_DIVERGED_PCSETUP_FAILED;
379:             VecSetInf(ksp->vec_sol);
380:           } else {
381:             SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Eigen estimator failed: %s at iteration %D",KSPConvergedReasons[reason],its);
382:           }
383:         }
384:       } else if (reason==KSP_CONVERGED_RTOL ||reason==KSP_CONVERGED_ATOL) {
385:         PetscInfo(ksp,"Eigen estimator converged prematurely. Should not happen except for small or low rank problem\n");
386:       } else {
387:         PetscInfo1(ksp,"Eigen estimator did not converge by iteration: %s\n",KSPConvergedReasons[reason]);
388:       }

390:       KSPChebyshevComputeExtremeEigenvalues_Private(cheb->kspest,&min,&max);

392:       cheb->emin_computed = min;
393:       cheb->emax_computed = max;
394:       cheb->emin = cheb->tform[0]*min + cheb->tform[1]*max;
395:       cheb->emax = cheb->tform[2]*min + cheb->tform[3]*max;

397:       cheb->amatid    = amatid;
398:       cheb->pmatid    = pmatid;
399:       cheb->amatstate = amatstate;
400:       cheb->pmatstate = pmatstate;
401:     }
402:   }

404:   ksp->its = 0;
405:   maxit    = ksp->max_it;

407:   /* These three point to the three active solutions, we
408:      rotate these three at each solution update */
409:   km1      = 0; k = 1; kp1 = 2;
410:   sol_orig = ksp->vec_sol; /* ksp->vec_sol will be asigned to rotating vector p[k], thus save its address */
411:   b        = ksp->vec_rhs;
412:   p[km1]   = sol_orig;
413:   p[k]     = ksp->work[0];
414:   p[kp1]   = ksp->work[1];
415:   r        = ksp->work[2];

417:   /* use scale*B as our preconditioner */
418:   scale = 2.0/(cheb->emax + cheb->emin);

420:   /*   -alpha <=  scale*lambda(B^{-1}A) <= alpha   */
421:   alpha     = 1.0 - scale*(cheb->emin);
422:   Gamma     = 1.0;
423:   mu        = 1.0/alpha;
424:   omegaprod = 2.0/alpha;

426:   c[km1] = 1.0;
427:   c[k]   = mu;

429:   if (!ksp->guess_zero) {
430:     KSP_MatMult(ksp,Amat,p[km1],r);     /*  r = b - A*p[km1] */
431:     VecAYPX(r,-1.0,b);
432:   } else {
433:     VecCopy(b,r);
434:   }

436:   KSP_PCApply(ksp,r,p[k]);  /* p[k] = scale B^{-1}r + p[km1] */
437:   VecAYPX(p[k],scale,p[km1]);

439:   for (i=0; i<maxit; i++) {
440:     PetscObjectSAWsTakeAccess((PetscObject)ksp);

442:     ksp->its++;
443:     PetscObjectSAWsGrantAccess((PetscObject)ksp);
444:     c[kp1] = 2.0*mu*c[k] - c[km1];
445:     omega  = omegaprod*c[k]/c[kp1];

447:     KSP_MatMult(ksp,Amat,p[k],r);          /*  r = b - Ap[k]    */
448:     VecAYPX(r,-1.0,b);
449:     KSP_PCApply(ksp,r,p[kp1]);             /*  p[kp1] = B^{-1}r  */
450:     ksp->vec_sol = p[k];

452:     /* calculate residual norm if requested */
453:     if (ksp->normtype != KSP_NORM_NONE || ksp->numbermonitors) {
454:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
455:         VecNorm(r,NORM_2,&rnorm);
456:       } else {
457:         VecNorm(p[kp1],NORM_2,&rnorm);
458:       }
459:       PetscObjectSAWsTakeAccess((PetscObject)ksp);
460:       ksp->rnorm   = rnorm;
461:       PetscObjectSAWsGrantAccess((PetscObject)ksp);
462:       KSPLogResidualHistory(ksp,rnorm);
463:       KSPMonitor(ksp,i,rnorm);
464:       (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
465:       if (ksp->reason) break;
466:     }

468:     /* y^{k+1} = omega(y^{k} - y^{k-1} + Gamma*r^{k}) + y^{k-1} */
469:     VecAXPBYPCZ(p[kp1],1.0-omega,omega,omega*Gamma*scale,p[km1],p[k]);

471:     ktmp = km1;
472:     km1  = k;
473:     k    = kp1;
474:     kp1  = ktmp;
475:   }
476:   if (!ksp->reason) {
477:     if (ksp->normtype != KSP_NORM_NONE) {
478:       KSP_MatMult(ksp,Amat,p[k],r);       /*  r = b - Ap[k]    */
479:       VecAYPX(r,-1.0,b);
480:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
481:         VecNorm(r,NORM_2,&rnorm);
482:       } else {
483:         KSP_PCApply(ksp,r,p[kp1]); /* p[kp1] = B^{-1}r */
484:         VecNorm(p[kp1],NORM_2,&rnorm);
485:       }
486:       PetscObjectSAWsTakeAccess((PetscObject)ksp);
487:       ksp->rnorm   = rnorm;
488:       PetscObjectSAWsGrantAccess((PetscObject)ksp);
489:       ksp->vec_sol = p[k];
490:       KSPLogResidualHistory(ksp,rnorm);
491:       KSPMonitor(ksp,i,rnorm);
492:     }
493:     if (ksp->its >= ksp->max_it) {
494:       if (ksp->normtype != KSP_NORM_NONE) {
495:         (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
496:         if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
497:       } else ksp->reason = KSP_CONVERGED_ITS;
498:     }
499:   }

501:   /* make sure solution is in vector x */
502:   ksp->vec_sol = sol_orig;
503:   if (k) {
504:     VecCopy(p[k],sol_orig);
505:   }
506:   return(0);
507: }

509: static  PetscErrorCode KSPView_Chebyshev(KSP ksp,PetscViewer viewer)
510: {
511:   KSP_Chebyshev  *cheb = (KSP_Chebyshev*)ksp->data;
513:   PetscBool      iascii;

516:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
517:   if (iascii) {
518:     PetscViewerASCIIPrintf(viewer,"  eigenvalue estimates used:  min = %g, max = %g\n",(double)cheb->emin,(double)cheb->emax);
519:     if (cheb->kspest) {
520:       PetscViewerASCIIPrintf(viewer,"  eigenvalues estimate via %s min %g, max %g\n",((PetscObject)(cheb->kspest))->type_name,(double)cheb->emin_computed,(double)cheb->emax_computed);
521:       PetscViewerASCIIPrintf(viewer,"  eigenvalues estimated using %s with translations  [%g %g; %g %g]\n",((PetscObject) cheb->kspest)->type_name,(double)cheb->tform[0],(double)cheb->tform[1],(double)cheb->tform[2],(double)cheb->tform[3]);
522:       PetscViewerASCIIPushTab(viewer);
523:       KSPView(cheb->kspest,viewer);
524:       PetscViewerASCIIPopTab(viewer);
525:       if (cheb->usenoisy) {
526:         PetscViewerASCIIPrintf(viewer,"  estimating eigenvalues using noisy right hand side\n");
527:       }
528:     }
529:   }
530:   return(0);
531: }

533: static PetscErrorCode KSPDestroy_Chebyshev(KSP ksp)
534: {
535:   KSP_Chebyshev  *cheb = (KSP_Chebyshev*)ksp->data;

539:   KSPDestroy(&cheb->kspest);
540:   PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevSetEigenvalues_C",NULL);
541:   PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigSet_C",NULL);
542:   PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigSetUseNoisy_C",NULL);
543:   PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigGetKSP_C",NULL);
544:   KSPDestroyDefault(ksp);
545:   return(0);
546: }

548: /*MC
549:      KSPCHEBYSHEV - The preconditioned Chebyshev iterative method

551:    Options Database Keys:
552: +   -ksp_chebyshev_eigenvalues <emin,emax> - set approximations to the smallest and largest eigenvalues
553:                   of the preconditioned operator. If these are accurate you will get much faster convergence.
554: .   -ksp_chebyshev_esteig <a,b,c,d> - estimate eigenvalues using a Krylov method, then use this
555:                          transform for Chebyshev eigenvalue bounds (KSPChebyshevEstEigSet())
556: .   -ksp_chebyshev_esteig_steps - number of estimation steps
557: -   -ksp_chebyshev_esteig_noisy - use noisy number generator to create right hand side for eigenvalue estimator

559:    Level: beginner

561:    Notes:
562:     The Chebyshev method requires both the matrix and preconditioner to
563:           be symmetric positive (semi) definite.
564:           Only support for left preconditioning.

566:           Chebyshev is configured as a smoother by default, targetting the "upper" part of the spectrum.
567:           The user should call KSPChebyshevSetEigenvalues() if they have eigenvalue estimates.

569: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
570:            KSPChebyshevSetEigenvalues(), KSPChebyshevEstEigSet(), KSPChebyshevEstEigSetUseNoisy()
571:            KSPRICHARDSON, KSPCG, PCMG

573: M*/

575: PETSC_EXTERN PetscErrorCode KSPCreate_Chebyshev(KSP ksp)
576: {
578:   KSP_Chebyshev  *chebyshevP;

581:   PetscNewLog(ksp,&chebyshevP);

583:   ksp->data = (void*)chebyshevP;
584:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
585:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
586:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
587:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);

589:   chebyshevP->emin = 0.;
590:   chebyshevP->emax = 0.;

592:   chebyshevP->tform[0] = 0.0;
593:   chebyshevP->tform[1] = 0.1;
594:   chebyshevP->tform[2] = 0;
595:   chebyshevP->tform[3] = 1.1;
596:   chebyshevP->eststeps = 10;
597:   chebyshevP->usenoisy = PETSC_TRUE;

599:   ksp->ops->setup          = KSPSetUp_Chebyshev;
600:   ksp->ops->solve          = KSPSolve_Chebyshev;
601:   ksp->ops->destroy        = KSPDestroy_Chebyshev;
602:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
603:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
604:   ksp->ops->setfromoptions = KSPSetFromOptions_Chebyshev;
605:   ksp->ops->view           = KSPView_Chebyshev;
606:   ksp->ops->reset          = KSPReset_Chebyshev;

608:   PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevSetEigenvalues_C",KSPChebyshevSetEigenvalues_Chebyshev);
609:   PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigSet_C",KSPChebyshevEstEigSet_Chebyshev);
610:   PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigSetUseNoisy_C",KSPChebyshevEstEigSetUseNoisy_Chebyshev);
611:   PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigGetKSP_C",KSPChebyshevEstEigGetKSP_Chebyshev);
612:   return(0);
613: }