Actual source code: cheby.c

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

  4: static const char *const KSPChebyshevKinds[] = {"FIRST", "FOURTH", "OPT_FOURTH", "KSPChebyshevKinds", "KSP_CHEBYSHEV_", NULL};

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

 10:   PetscFunctionBegin;
 11:   if (cheb->kspest) PetscCall(KSPReset(cheb->kspest));
 12:   PetscFunctionReturn(PETSC_SUCCESS);
 13: }

 15: /*
 16:  * Must be passed a KSP solver that has "converged", with KSPSetComputeEigenvalues() called before the solve
 17:  */
 18: static PetscErrorCode KSPChebyshevComputeExtremeEigenvalues_Private(KSP kspest, PetscReal *emin, PetscReal *emax)
 19: {
 20:   PetscInt   n, neig;
 21:   PetscReal *re, *im, min, max;

 23:   PetscFunctionBegin;
 24:   PetscCall(KSPGetIterationNumber(kspest, &n));
 25:   PetscCall(PetscMalloc2(n, &re, n, &im));
 26:   PetscCall(KSPComputeEigenvalues(kspest, n, re, im, &neig));
 27:   min = PETSC_MAX_REAL;
 28:   max = PETSC_MIN_REAL;
 29:   for (n = 0; n < neig; n++) {
 30:     min = PetscMin(min, re[n]);
 31:     max = PetscMax(max, re[n]);
 32:   }
 33:   PetscCall(PetscFree2(re, im));
 34:   *emax = max;
 35:   *emin = min;
 36:   PetscFunctionReturn(PETSC_SUCCESS);
 37: }

 39: static PetscErrorCode KSPChebyshevGetEigenvalues_Chebyshev(KSP ksp, PetscReal *emax, PetscReal *emin)
 40: {
 41:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

 43:   PetscFunctionBegin;
 44:   *emax = 0;
 45:   *emin = 0;
 46:   if (cheb->emax != 0.) {
 47:     *emax = cheb->emax;
 48:   } else if (cheb->emax_computed != 0.) {
 49:     *emax = cheb->tform[2] * cheb->emin_computed + cheb->tform[3] * cheb->emax_computed;
 50:   } else if (cheb->emax_provided != 0.) {
 51:     *emax = cheb->tform[2] * cheb->emin_provided + cheb->tform[3] * cheb->emax_provided;
 52:   }
 53:   if (cheb->emin != 0.) {
 54:     *emin = cheb->emin;
 55:   } else if (cheb->emin_computed != 0.) {
 56:     *emin = cheb->tform[0] * cheb->emin_computed + cheb->tform[1] * cheb->emax_computed;
 57:   } else if (cheb->emin_provided != 0.) {
 58:     *emin = cheb->tform[0] * cheb->emin_provided + cheb->tform[1] * cheb->emax_provided;
 59:   }
 60:   PetscFunctionReturn(PETSC_SUCCESS);
 61: }

 63: static PetscErrorCode KSPChebyshevSetEigenvalues_Chebyshev(KSP ksp, PetscReal emax, PetscReal emin)
 64: {
 65:   KSP_Chebyshev *chebyshevP = (KSP_Chebyshev *)ksp->data;

 67:   PetscFunctionBegin;
 68:   PetscCheck(emax > emin, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "Maximum eigenvalue must be larger than minimum: max %g min %g", (double)emax, (double)emin);
 69:   PetscCheck(emax * emin > 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "Both eigenvalues must be of the same sign: max %g min %g", (double)emax, (double)emin);
 70:   chebyshevP->emax = emax;
 71:   chebyshevP->emin = emin;

 73:   PetscCall(KSPChebyshevEstEigSet(ksp, 0., 0., 0., 0.)); /* Destroy any estimation setup */
 74:   PetscFunctionReturn(PETSC_SUCCESS);
 75: }

 77: static PetscErrorCode KSPChebyshevEstEigSet_Chebyshev(KSP ksp, PetscReal a, PetscReal b, PetscReal c, PetscReal d)
 78: {
 79:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

 81:   PetscFunctionBegin;
 82:   if (a != 0.0 || b != 0.0 || c != 0.0 || d != 0.0) {
 83:     if ((cheb->emin_provided == 0. || cheb->emax_provided == 0.) && !cheb->kspest) { /* should this block of code be moved to KSPSetUp_Chebyshev()? */
 84:       PetscCall(KSPCreate(PetscObjectComm((PetscObject)ksp), &cheb->kspest));
 85:       PetscCall(KSPSetErrorIfNotConverged(cheb->kspest, ksp->errorifnotconverged));
 86:       PetscCall(PetscObjectIncrementTabLevel((PetscObject)cheb->kspest, (PetscObject)ksp, 1));
 87:       /* use PetscObjectSet/AppendOptionsPrefix() instead of KSPSet/AppendOptionsPrefix() so that the PC prefix is not changed */
 88:       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)cheb->kspest, ((PetscObject)ksp)->prefix));
 89:       PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)cheb->kspest, "esteig_"));
 90:       PetscCall(KSPSetSkipPCSetFromOptions(cheb->kspest, PETSC_TRUE));

 92:       PetscCall(KSPSetComputeEigenvalues(cheb->kspest, PETSC_TRUE));

 94:       /* We cannot turn off convergence testing because GMRES will break down if you attempt to keep iterating after a zero norm is obtained */
 95:       PetscCall(KSPSetTolerances(cheb->kspest, 1.e-12, PETSC_DEFAULT, PETSC_DEFAULT, cheb->eststeps));
 96:     }
 97:     if (a >= 0) cheb->tform[0] = a;
 98:     if (b >= 0) cheb->tform[1] = b;
 99:     if (c >= 0) cheb->tform[2] = c;
100:     if (d >= 0) cheb->tform[3] = d;
101:     cheb->amatid    = 0;
102:     cheb->pmatid    = 0;
103:     cheb->amatstate = -1;
104:     cheb->pmatstate = -1;
105:   } else {
106:     PetscCall(KSPDestroy(&cheb->kspest));
107:   }
108:   PetscFunctionReturn(PETSC_SUCCESS);
109: }

111: static PetscErrorCode KSPChebyshevEstEigSetUseNoisy_Chebyshev(KSP ksp, PetscBool use)
112: {
113:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

115:   PetscFunctionBegin;
116:   cheb->usenoisy = use;
117:   PetscFunctionReturn(PETSC_SUCCESS);
118: }

120: static PetscErrorCode KSPChebyshevSetKind_Chebyshev(KSP ksp, KSPChebyshevKind kind)
121: {
122:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

124:   PetscFunctionBegin;
125:   cheb->chebykind = kind;
126:   PetscFunctionReturn(PETSC_SUCCESS);
127: }

129: /*@
130:    KSPChebyshevSetEigenvalues - Sets estimates for the extreme eigenvalues
131:    of the preconditioned problem.

133:    Logically Collective

135:    Input Parameters:
136: +  ksp - the Krylov space context
137:    emax - the eigenvalue maximum estimate
138: -  emin - the eigenvalue minimum estimate

140:   Options Database Key:
141: .  -ksp_chebyshev_eigenvalues emin,emax - extreme eigenvalues

143:    Notes:
144:    Call `KSPChebyshevEstEigSet()` or use the option -ksp_chebyshev_esteig a,b,c,d to have the KSP
145:    estimate the eigenvalues and use these estimated values automatically.

147:    When `KSPCHEBYSHEV` is used as a smoother, one often wants to target a portion of the spectrum rather than the entire
148:    spectrum. This function takes the range of target eigenvalues for Chebyshev, which will often slightly over-estimate
149:    the largest eigenvalue of the actual operator (for safety) and greatly overestimate the smallest eigenvalue to
150:    improve the smoothing properties of Chebyshev iteration on the higher frequencies in the spectrum.

152:    Level: intermediate

154: .seealso: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevEstEigSet()`,
155: @*/
156: PetscErrorCode KSPChebyshevSetEigenvalues(KSP ksp, PetscReal emax, PetscReal emin)
157: {
158:   PetscFunctionBegin;
162:   PetscTryMethod(ksp, "KSPChebyshevSetEigenvalues_C", (KSP, PetscReal, PetscReal), (ksp, emax, emin));
163:   PetscFunctionReturn(PETSC_SUCCESS);
164: }

166: /*@
167:    KSPChebyshevEstEigSet - Automatically estimate the eigenvalues to use for Chebyshev

169:    Logically Collective

171:    Input Parameters:
172: +  ksp - the Krylov space context
173: .  a - multiple of min eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE)
174: .  b - multiple of max eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE)
175: .  c - multiple of min eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE)
176: -  d - multiple of max eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE)

178:   Options Database Key:
179: .  -ksp_chebyshev_esteig a,b,c,d - estimate eigenvalues using a Krylov method, then use this transform for Chebyshev eigenvalue bounds

181:    Notes:
182:    The Chebyshev bounds are set using
183: .vb
184:    minbound = a*minest + b*maxest
185:    maxbound = c*minest + d*maxest
186: .ve
187:    The default configuration targets the upper part of the spectrum for use as a multigrid smoother, so only the maximum eigenvalue estimate is used.
188:    The minimum eigenvalue estimate obtained by Krylov iteration is typically not accurate until the method has converged.

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

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

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

196:    Level: intermediate

198: .seealso: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevEstEigSet()`, `KSPChebyshevEstEigSetUseNoisy()`, `KSPChebyshevEstEigGetKSP()`
199: @*/
200: PetscErrorCode KSPChebyshevEstEigSet(KSP ksp, PetscReal a, PetscReal b, PetscReal c, PetscReal d)
201: {
202:   PetscFunctionBegin;
208:   PetscTryMethod(ksp, "KSPChebyshevEstEigSet_C", (KSP, PetscReal, PetscReal, PetscReal, PetscReal), (ksp, a, b, c, d));
209:   PetscFunctionReturn(PETSC_SUCCESS);
210: }

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

215:    Logically Collective

217:    Input Parameters:
218: +  ksp - linear solver context
219: -  use - `PETSC_TRUE` to use noisy

221:    Options Database Key:
222: .  -ksp_chebyshev_esteig_noisy <true,false> - Use noisy right hand side for estimate

224:    Note:
225:     This allegedly works better for multigrid smoothers

227:   Level: intermediate

229: .seealso: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevEstEigSet()`, `KSPChebyshevEstEigGetKSP()`
230: @*/
231: PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP ksp, PetscBool use)
232: {
233:   PetscFunctionBegin;
236:   PetscTryMethod(ksp, "KSPChebyshevEstEigSetUseNoisy_C", (KSP, PetscBool), (ksp, use));
237:   PetscFunctionReturn(PETSC_SUCCESS);
238: }

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

245:   Input Parameter:
246: . ksp - the Krylov space context

248:   Output Parameter:
249: . kspest - the eigenvalue estimation Krylov space context

251:   Level: advanced

253: .seealso: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevEstEigSet()`
254: @*/
255: PetscErrorCode KSPChebyshevEstEigGetKSP(KSP ksp, KSP *kspest)
256: {
257:   PetscFunctionBegin;
260:   *kspest = NULL;
261:   PetscTryMethod(ksp, "KSPChebyshevEstEigGetKSP_C", (KSP, KSP *), (ksp, kspest));
262:   PetscFunctionReturn(PETSC_SUCCESS);
263: }

265: /*@
266:   KSPChebyshevSetKind - set the kind of Chebyshev polynomial to use

268:   Logically Collective

270:   Input Parameters:
271: + ksp  - Linear solver context
272: - kind - The kind of Chebyshev polynomial to use

274:   Options Database Key:
275: . -ksp_chebyshev_kind <kind> - which kind of Chebyshev polynomial to use

277:   Note:
278:   When using multigrid methods for problems with a poor quality coarse space (e.g., due to anisotropy or aggressive
279:   coarsening), it is necessary for the smoother to handle smaller eigenvalues. With first-kind Chebyshev smoothing, this
280:   requires using higher degree Chebyhev polynomials and reducing the lower end of the target spectrum, at which point
281:   the whole target spectrum experiences about the same damping. Fourth kind Chebyshev polynomials (and the "optimized"
282:   fourth kind) avoid the ad-hoc choice of lower bound and extend smoothing to smaller eigenvalues while preferentially
283:   smoothing higher modes faster as needed to minimize the energy norm of the error.

285:   References:
286: +  * - Malachi Phillips and Paul Fischer, Optimal Chebyshev Smoothers and One-sided V-cycles, https://arxiv.org/abs/2210.03179.
287: -  * - James Lottes, Optimal Polynomial Smoothers for Multigrid V-cycles, https://arxiv.org/abs/2202.08830.

289:   Level: intermediate

291: .seealso: [](chapter_ksp), `KSPCHEBYSHEV` `KSPChebyshevKind`
292: @*/
293: PetscErrorCode KSPChebyshevSetKind(KSP ksp, KSPChebyshevKind kind)
294: {
295:   PetscFunctionBegin;
298:   PetscUseMethod(ksp, "KSPChebyshevSetKind_C", (KSP, KSPChebyshevKind), (ksp, kind));
299:   PetscFunctionReturn(PETSC_SUCCESS);
300: }

302: static PetscErrorCode KSPChebyshevEstEigGetKSP_Chebyshev(KSP ksp, KSP *kspest)
303: {
304:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

306:   PetscFunctionBegin;
307:   *kspest = cheb->kspest;
308:   PetscFunctionReturn(PETSC_SUCCESS);
309: }

311: static PetscErrorCode KSPSetFromOptions_Chebyshev(KSP ksp, PetscOptionItems *PetscOptionsObject)
312: {
313:   KSP_Chebyshev *cheb    = (KSP_Chebyshev *)ksp->data;
314:   PetscInt       neigarg = 2, nestarg = 4;
315:   PetscReal      eminmax[2] = {0., 0.};
316:   PetscReal      tform[4]   = {PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE};
317:   PetscBool      flgeig, flgest;

319:   PetscFunctionBegin;
320:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP Chebyshev Options");
321:   PetscCall(PetscOptionsInt("-ksp_chebyshev_esteig_steps", "Number of est steps in Chebyshev", "", cheb->eststeps, &cheb->eststeps, NULL));
322:   PetscCall(PetscOptionsRealArray("-ksp_chebyshev_eigenvalues", "extreme eigenvalues", "KSPChebyshevSetEigenvalues", eminmax, &neigarg, &flgeig));
323:   if (flgeig) {
324:     PetscCheck(neigarg == 2, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "-ksp_chebyshev_eigenvalues: must specify 2 parameters, min and max eigenvalues");
325:     PetscCall(KSPChebyshevSetEigenvalues(ksp, eminmax[1], eminmax[0]));
326:   }
327:   PetscCall(PetscOptionsRealArray("-ksp_chebyshev_esteig", "estimate eigenvalues using a Krylov method, then use this transform for Chebyshev eigenvalue bounds", "KSPChebyshevEstEigSet", tform, &nestarg, &flgest));
328:   if (flgest) {
329:     switch (nestarg) {
330:     case 0:
331:       PetscCall(KSPChebyshevEstEigSet(ksp, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
332:       break;
333:     case 2: /* Base everything on the max eigenvalues */
334:       PetscCall(KSPChebyshevEstEigSet(ksp, PETSC_DECIDE, tform[0], PETSC_DECIDE, tform[1]));
335:       break;
336:     case 4: /* Use the full 2x2 linear transformation */
337:       PetscCall(KSPChebyshevEstEigSet(ksp, tform[0], tform[1], tform[2], tform[3]));
338:       break;
339:     default:
340:       SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "Must specify either 0, 2, or 4 parameters for eigenvalue estimation");
341:     }
342:   }

344:   cheb->chebykind = KSP_CHEBYSHEV_FIRST; /* Default to 1st-kind Chebyshev polynomial */
345:   PetscCall(PetscOptionsEnum("-ksp_chebyshev_kind", "Type of Chebyshev polynomial", "KSPChebyshevKind", KSPChebyshevKinds, (PetscEnum)cheb->chebykind, (PetscEnum *)&cheb->chebykind, NULL));

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

350:   if (cheb->kspest) {
351:     PetscCall(PetscOptionsBool("-ksp_chebyshev_esteig_noisy", "Use noisy right hand side for estimate", "KSPChebyshevEstEigSetUseNoisy", cheb->usenoisy, &cheb->usenoisy, NULL));
352:     PetscCall(KSPSetFromOptions(cheb->kspest));
353:   }
354:   PetscOptionsHeadEnd();
355:   PetscFunctionReturn(PETSC_SUCCESS);
356: }

358: static PetscErrorCode KSPSolve_Chebyshev_FirstKind(KSP ksp)
359: {
360:   PetscInt    k, kp1, km1, ktmp, i;
361:   PetscScalar alpha, omegaprod, mu, omega, Gamma, c[3], scale;
362:   PetscReal   rnorm = 0.0, emax, emin;
363:   Vec         sol_orig, b, p[3], r;
364:   Mat         Amat, Pmat;
365:   PetscBool   diagonalscale;

367:   PetscFunctionBegin;
368:   PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
369:   PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);

371:   PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
372:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
373:   ksp->its = 0;
374:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
375:   /* These three point to the three active solutions, we
376:      rotate these three at each solution update */
377:   km1      = 0;
378:   k        = 1;
379:   kp1      = 2;
380:   sol_orig = ksp->vec_sol; /* ksp->vec_sol will be assigned to rotating vector p[k], thus save its address */
381:   b        = ksp->vec_rhs;
382:   p[km1]   = sol_orig;
383:   p[k]     = ksp->work[0];
384:   p[kp1]   = ksp->work[1];
385:   r        = ksp->work[2];

387:   PetscCall(KSPChebyshevGetEigenvalues_Chebyshev(ksp, &emax, &emin));
388:   /* use scale*B as our preconditioner */
389:   scale = 2.0 / (emax + emin);

391:   /*   -alpha <=  scale*lambda(B^{-1}A) <= alpha   */
392:   alpha     = 1.0 - scale * emin;
393:   Gamma     = 1.0;
394:   mu        = 1.0 / alpha;
395:   omegaprod = 2.0 / alpha;

397:   c[km1] = 1.0;
398:   c[k]   = mu;

400:   if (!ksp->guess_zero) {
401:     PetscCall(KSP_MatMult(ksp, Amat, sol_orig, r)); /*  r = b - A*p[km1] */
402:     PetscCall(VecAYPX(r, -1.0, b));
403:   } else {
404:     PetscCall(VecCopy(b, r));
405:   }

407:   /* calculate residual norm if requested, we have done one iteration */
408:   if (ksp->normtype) {
409:     switch (ksp->normtype) {
410:     case KSP_NORM_PRECONDITIONED:
411:       PetscCall(KSP_PCApply(ksp, r, p[k])); /* p[k] = B^{-1}r */
412:       PetscCall(VecNorm(p[k], NORM_2, &rnorm));
413:       break;
414:     case KSP_NORM_UNPRECONDITIONED:
415:     case KSP_NORM_NATURAL:
416:       PetscCall(VecNorm(r, NORM_2, &rnorm));
417:       break;
418:     default:
419:       SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
420:     }
421:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
422:     ksp->rnorm = rnorm;
423:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
424:     PetscCall(KSPLogResidualHistory(ksp, rnorm));
425:     PetscCall(KSPLogErrorHistory(ksp));
426:     PetscCall(KSPMonitor(ksp, 0, rnorm));
427:     PetscCall((*ksp->converged)(ksp, 0, rnorm, &ksp->reason, ksp->cnvP));
428:   } else ksp->reason = KSP_CONVERGED_ITERATING;
429:   if (ksp->reason || ksp->max_it == 0) {
430:     if (ksp->max_it == 0) ksp->reason = KSP_DIVERGED_ITS; /* This for a V(0,x) cycle */
431:     PetscFunctionReturn(PETSC_SUCCESS);
432:   }
433:   if (ksp->normtype != KSP_NORM_PRECONDITIONED) { PetscCall(KSP_PCApply(ksp, r, p[k])); /* p[k] = B^{-1}r */ }
434:   PetscCall(VecAYPX(p[k], scale, p[km1])); /* p[k] = scale B^{-1}r + p[km1] */
435:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
436:   ksp->its = 1;
437:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));

439:   for (i = 1; i < ksp->max_it; i++) {
440:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
441:     ksp->its++;
442:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));

444:     PetscCall(KSP_MatMult(ksp, Amat, p[k], r)); /*  r = b - Ap[k]    */
445:     PetscCall(VecAYPX(r, -1.0, b));
446:     /* calculate residual norm if requested */
447:     if (ksp->normtype) {
448:       switch (ksp->normtype) {
449:       case KSP_NORM_PRECONDITIONED:
450:         PetscCall(KSP_PCApply(ksp, r, p[kp1])); /*  p[kp1] = B^{-1}r  */
451:         PetscCall(VecNorm(p[kp1], NORM_2, &rnorm));
452:         break;
453:       case KSP_NORM_UNPRECONDITIONED:
454:       case KSP_NORM_NATURAL:
455:         PetscCall(VecNorm(r, NORM_2, &rnorm));
456:         break;
457:       default:
458:         rnorm = 0.0;
459:         break;
460:       }
461:       KSPCheckNorm(ksp, rnorm);
462:       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
463:       ksp->rnorm = rnorm;
464:       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
465:       PetscCall(KSPLogResidualHistory(ksp, rnorm));
466:       PetscCall(KSPMonitor(ksp, i, rnorm));
467:       PetscCall((*ksp->converged)(ksp, i, rnorm, &ksp->reason, ksp->cnvP));
468:       if (ksp->reason) break;
469:       if (ksp->normtype != KSP_NORM_PRECONDITIONED) { PetscCall(KSP_PCApply(ksp, r, p[kp1])); /*  p[kp1] = B^{-1}r  */ }
470:     } else {
471:       PetscCall(KSP_PCApply(ksp, r, p[kp1])); /*  p[kp1] = B^{-1}r  */
472:     }
473:     ksp->vec_sol = p[k];
474:     PetscCall(KSPLogErrorHistory(ksp));

476:     c[kp1] = 2.0 * mu * c[k] - c[km1];
477:     omega  = omegaprod * c[k] / c[kp1];

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

482:     ktmp = km1;
483:     km1  = k;
484:     k    = kp1;
485:     kp1  = ktmp;
486:   }
487:   if (!ksp->reason) {
488:     if (ksp->normtype) {
489:       PetscCall(KSP_MatMult(ksp, Amat, p[k], r)); /*  r = b - Ap[k]    */
490:       PetscCall(VecAYPX(r, -1.0, b));
491:       switch (ksp->normtype) {
492:       case KSP_NORM_PRECONDITIONED:
493:         PetscCall(KSP_PCApply(ksp, r, p[kp1])); /* p[kp1] = B^{-1}r */
494:         PetscCall(VecNorm(p[kp1], NORM_2, &rnorm));
495:         break;
496:       case KSP_NORM_UNPRECONDITIONED:
497:       case KSP_NORM_NATURAL:
498:         PetscCall(VecNorm(r, NORM_2, &rnorm));
499:         break;
500:       default:
501:         rnorm = 0.0;
502:         break;
503:       }
504:       KSPCheckNorm(ksp, rnorm);
505:       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
506:       ksp->rnorm = rnorm;
507:       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
508:       PetscCall(KSPLogResidualHistory(ksp, rnorm));
509:       PetscCall(KSPMonitor(ksp, i, rnorm));
510:     }
511:     if (ksp->its >= ksp->max_it) {
512:       if (ksp->normtype != KSP_NORM_NONE) {
513:         PetscCall((*ksp->converged)(ksp, i, rnorm, &ksp->reason, ksp->cnvP));
514:         if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
515:       } else ksp->reason = KSP_CONVERGED_ITS;
516:     }
517:   }

519:   /* make sure solution is in vector x */
520:   ksp->vec_sol = sol_orig;
521:   if (k) PetscCall(VecCopy(p[k], sol_orig));
522:   if (ksp->reason == KSP_CONVERGED_ITS) PetscCall(KSPLogErrorHistory(ksp));
523:   PetscFunctionReturn(PETSC_SUCCESS);
524: }

526: static PetscErrorCode KSPSolve_Chebyshev_FourthKind(KSP ksp)
527: {
528:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
529:   PetscInt       i;
530:   PetscScalar    scale, rScale, dScale;
531:   PetscReal      rnorm = 0.0, emax, emin;
532:   Vec            x, b, d, r, Br;
533:   Mat            Amat, Pmat;
534:   PetscBool      diagonalscale;
535:   PetscReal     *betas = cheb->betas;

537:   PetscFunctionBegin;
538:   PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
539:   PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);

541:   PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
542:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
543:   ksp->its = 0;
544:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));

546:   x  = ksp->vec_sol;
547:   b  = ksp->vec_rhs;
548:   r  = ksp->work[0];
549:   d  = ksp->work[1];
550:   Br = ksp->work[2];

552:   PetscCall(KSPChebyshevGetEigenvalues_Chebyshev(ksp, &emax, &emin));
553:   /* use scale*B as our preconditioner */
554:   scale = 1.0 / emax;

556:   if (!ksp->guess_zero) {
557:     PetscCall(KSP_MatMult(ksp, Amat, x, r)); /*  r = b - A*x */
558:     PetscCall(VecAYPX(r, -1.0, b));
559:   } else {
560:     PetscCall(VecCopy(b, r));
561:   }

563:   /* calculate residual norm if requested, we have done one iteration */
564:   if (ksp->normtype) {
565:     switch (ksp->normtype) {
566:     case KSP_NORM_PRECONDITIONED:
567:       PetscCall(KSP_PCApply(ksp, r, Br)); /* Br = B^{-1}r */
568:       PetscCall(VecNorm(Br, NORM_2, &rnorm));
569:       break;
570:     case KSP_NORM_UNPRECONDITIONED:
571:     case KSP_NORM_NATURAL:
572:       PetscCall(VecNorm(r, NORM_2, &rnorm));
573:       break;
574:     default:
575:       SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
576:     }
577:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
578:     ksp->rnorm = rnorm;
579:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
580:     PetscCall(KSPLogResidualHistory(ksp, rnorm));
581:     PetscCall(KSPLogErrorHistory(ksp));
582:     PetscCall(KSPMonitor(ksp, 0, rnorm));
583:     PetscCall((*ksp->converged)(ksp, 0, rnorm, &ksp->reason, ksp->cnvP));
584:   } else ksp->reason = KSP_CONVERGED_ITERATING;
585:   if (ksp->reason || ksp->max_it == 0) {
586:     if (ksp->max_it == 0) ksp->reason = KSP_DIVERGED_ITS; /* This for a V(0,x) cycle */
587:     PetscFunctionReturn(PETSC_SUCCESS);
588:   }
589:   if (ksp->normtype != KSP_NORM_PRECONDITIONED) { PetscCall(KSP_PCApply(ksp, r, Br)); /* Br = B^{-1}r */ }
590:   PetscCall(VecAXPBY(d, 4.0 / 3.0 * scale, 0.0, Br)); /* d = 4/3 * scale B^{-1}r */
591:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
592:   ksp->its = 1;
593:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));

595:   for (i = 1; i < ksp->max_it; i++) {
596:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
597:     ksp->its++;
598:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));

600:     PetscCall(VecAXPBY(x, betas[i - 1], 1.0, d)); /* x = x + \beta_k d */

602:     PetscCall(KSP_MatMult(ksp, Amat, d, Br)); /*  r = r - Ad */
603:     PetscCall(VecAXPBY(r, -1.0, 1.0, Br));

605:     /* calculate residual norm if requested */
606:     if (ksp->normtype) {
607:       switch (ksp->normtype) {
608:       case KSP_NORM_PRECONDITIONED:
609:         PetscCall(KSP_PCApply(ksp, r, Br)); /*  Br = B^{-1}r  */
610:         PetscCall(VecNorm(Br, NORM_2, &rnorm));
611:         break;
612:       case KSP_NORM_UNPRECONDITIONED:
613:       case KSP_NORM_NATURAL:
614:         PetscCall(VecNorm(r, NORM_2, &rnorm));
615:         break;
616:       default:
617:         rnorm = 0.0;
618:         break;
619:       }
620:       KSPCheckNorm(ksp, rnorm);
621:       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
622:       ksp->rnorm = rnorm;
623:       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
624:       PetscCall(KSPLogResidualHistory(ksp, rnorm));
625:       PetscCall(KSPMonitor(ksp, i, rnorm));
626:       PetscCall((*ksp->converged)(ksp, i, rnorm, &ksp->reason, ksp->cnvP));
627:       if (ksp->reason) break;
628:       if (ksp->normtype != KSP_NORM_PRECONDITIONED) { PetscCall(KSP_PCApply(ksp, r, Br)); /*  Br = B^{-1}r  */ }
629:     } else {
630:       PetscCall(KSP_PCApply(ksp, r, Br)); /*  Br = B^{-1}r  */
631:     }
632:     PetscCall(KSPLogErrorHistory(ksp));

634:     rScale = scale * (8.0 * i + 4.0) / (2.0 * i + 3.0);
635:     dScale = (2.0 * i - 1.0) / (2.0 * i + 3.0);

637:     /* d_k+1 = \dfrac{2k-1}{2k+3} d_k + \dfrac{8k+4}{2k+3} \dfrac{1}{\rho(SA)} Br */
638:     PetscCall(VecAXPBY(d, rScale, dScale, Br));
639:   }

641:   /* on last pass, update solution vector */
642:   PetscCall(VecAXPBY(x, betas[ksp->max_it - 1], 1.0, d)); /* x = x + d */

644:   if (!ksp->reason) {
645:     if (ksp->normtype) {
646:       PetscCall(KSP_MatMult(ksp, Amat, x, r)); /*  r = b - Ax    */
647:       PetscCall(VecAYPX(r, -1.0, b));
648:       switch (ksp->normtype) {
649:       case KSP_NORM_PRECONDITIONED:
650:         PetscCall(KSP_PCApply(ksp, r, Br)); /* Br= B^{-1}r */
651:         PetscCall(VecNorm(Br, NORM_2, &rnorm));
652:         break;
653:       case KSP_NORM_UNPRECONDITIONED:
654:       case KSP_NORM_NATURAL:
655:         PetscCall(VecNorm(r, NORM_2, &rnorm));
656:         break;
657:       default:
658:         rnorm = 0.0;
659:         break;
660:       }
661:       KSPCheckNorm(ksp, rnorm);
662:       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
663:       ksp->rnorm = rnorm;
664:       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
665:       PetscCall(KSPLogResidualHistory(ksp, rnorm));
666:       PetscCall(KSPMonitor(ksp, i, rnorm));
667:     }
668:     if (ksp->its >= ksp->max_it) {
669:       if (ksp->normtype != KSP_NORM_NONE) {
670:         PetscCall((*ksp->converged)(ksp, i, rnorm, &ksp->reason, ksp->cnvP));
671:         if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
672:       } else ksp->reason = KSP_CONVERGED_ITS;
673:     }
674:   }

676:   if (ksp->reason == KSP_CONVERGED_ITS) PetscCall(KSPLogErrorHistory(ksp));
677:   PetscFunctionReturn(PETSC_SUCCESS);
678: }

680: static PetscErrorCode KSPView_Chebyshev(KSP ksp, PetscViewer viewer)
681: {
682:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
683:   PetscBool      iascii;

685:   PetscFunctionBegin;
686:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
687:   if (iascii) {
688:     switch (cheb->chebykind) {
689:     case KSP_CHEBYSHEV_FIRST:
690:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Chebyshev polynomial of first kind\n"));
691:       break;
692:     case KSP_CHEBYSHEV_FOURTH:
693:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Chebyshev polynomial of fourth kind\n"));
694:       break;
695:     case KSP_CHEBYSHEV_OPT_FOURTH:
696:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Chebyshev polynomial of opt. fourth kind\n"));
697:       break;
698:     }
699:     PetscReal emax, emin;
700:     PetscCall(KSPChebyshevGetEigenvalues_Chebyshev(ksp, &emax, &emin));
701:     PetscCall(PetscViewerASCIIPrintf(viewer, "  eigenvalue targets used: min %g, max %g\n", (double)emin, (double)emax));
702:     if (cheb->kspest) {
703:       PetscCall(PetscViewerASCIIPrintf(viewer, "  eigenvalues estimated via %s: min %g, max %g\n", ((PetscObject)(cheb->kspest))->type_name, (double)cheb->emin_computed, (double)cheb->emax_computed));
704:       PetscCall(PetscViewerASCIIPrintf(viewer, "  eigenvalues estimated using %s with transform: [%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]));
705:       PetscCall(PetscViewerASCIIPushTab(viewer));
706:       PetscCall(KSPView(cheb->kspest, viewer));
707:       PetscCall(PetscViewerASCIIPopTab(viewer));
708:       if (cheb->usenoisy) PetscCall(PetscViewerASCIIPrintf(viewer, "  estimating eigenvalues using noisy right hand side\n"));
709:     } else if (cheb->emax_provided != 0.) {
710:       PetscCall(PetscViewerASCIIPrintf(viewer, "  eigenvalues provided (min %g, max %g) with transform: [%g %g; %g %g]\n", (double)cheb->emin_provided, (double)cheb->emax_provided, (double)cheb->tform[0], (double)cheb->tform[1], (double)cheb->tform[2],
711:                                        (double)cheb->tform[3]));
712:     }
713:   }
714:   PetscFunctionReturn(PETSC_SUCCESS);
715: }

717: static PetscErrorCode KSPSetUp_Chebyshev(KSP ksp)
718: {
719:   KSP_Chebyshev   *cheb = (KSP_Chebyshev *)ksp->data;
720:   PetscBool        isset, flg;
721:   Mat              Pmat, Amat;
722:   PetscObjectId    amatid, pmatid;
723:   PetscObjectState amatstate, pmatstate;

725:   PetscFunctionBegin;
726:   switch (cheb->chebykind) {
727:   case KSP_CHEBYSHEV_FIRST:
728:     ksp->ops->solve = KSPSolve_Chebyshev_FirstKind;
729:     break;
730:   case KSP_CHEBYSHEV_FOURTH:
731:   case KSP_CHEBYSHEV_OPT_FOURTH:
732:     ksp->ops->solve = KSPSolve_Chebyshev_FourthKind;
733:     break;
734:   }

736:   if (ksp->max_it > cheb->num_betas_alloc) {
737:     PetscCall(PetscFree(cheb->betas));
738:     PetscCall(PetscMalloc1(ksp->max_it, &cheb->betas));
739:     cheb->num_betas_alloc = ksp->max_it;
740:   }

742:   // coefficients for 4th-kind Chebyshev
743:   for (PetscInt i = 0; i < ksp->max_it; i++) cheb->betas[i] = 1.0;

745:   // coefficients for optimized 4th-kind Chebyshev
746:   if (cheb->chebykind == KSP_CHEBYSHEV_OPT_FOURTH) PetscCall(KSPChebyshevGetBetas_Private(ksp));

748:   PetscCall(KSPSetWorkVecs(ksp, 3));
749:   if (cheb->emin == 0. || cheb->emax == 0.) { // User did not specify eigenvalues
750:     PC pc;

752:     PetscCall(KSPGetPC(ksp, &pc));
753:     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &flg));
754:     if (!flg) { // Provided estimates are only relevant for Jacobi
755:       cheb->emax_provided = 0;
756:       cheb->emin_provided = 0;
757:     }
758:     /* We need to estimate eigenvalues */
759:     if (!cheb->kspest) PetscCall(KSPChebyshevEstEigSet(ksp, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
760:   }
761:   if (cheb->kspest) {
762:     PetscCall(KSPGetOperators(ksp, &Amat, &Pmat));
763:     PetscCall(MatIsSPDKnown(Pmat, &isset, &flg));
764:     if (isset && flg) {
765:       const char *prefix;

767:       PetscCall(KSPGetOptionsPrefix(cheb->kspest, &prefix));
768:       PetscCall(PetscOptionsHasName(NULL, prefix, "-ksp_type", &flg));
769:       if (!flg) PetscCall(KSPSetType(cheb->kspest, KSPCG));
770:     }
771:     PetscCall(PetscObjectGetId((PetscObject)Amat, &amatid));
772:     PetscCall(PetscObjectGetId((PetscObject)Pmat, &pmatid));
773:     PetscCall(PetscObjectStateGet((PetscObject)Amat, &amatstate));
774:     PetscCall(PetscObjectStateGet((PetscObject)Pmat, &pmatstate));
775:     if (amatid != cheb->amatid || pmatid != cheb->pmatid || amatstate != cheb->amatstate || pmatstate != cheb->pmatstate) {
776:       PetscReal          max = 0.0, min = 0.0;
777:       Vec                B;
778:       KSPConvergedReason reason;

780:       PetscCall(KSPSetPC(cheb->kspest, ksp->pc));
781:       if (cheb->usenoisy) {
782:         B = ksp->work[1];
783:         PetscCall(KSPSetNoisy_Private(B));
784:       } else {
785:         PetscBool change;

787:         PetscCheck(ksp->vec_rhs, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Chebyshev must use a noisy right hand side to estimate the eigenvalues when no right hand side is available");
788:         PetscCall(PCPreSolveChangeRHS(ksp->pc, &change));
789:         if (change) {
790:           B = ksp->work[1];
791:           PetscCall(VecCopy(ksp->vec_rhs, B));
792:         } else B = ksp->vec_rhs;
793:       }
794:       PetscCall(KSPSolve(cheb->kspest, B, ksp->work[0]));
795:       PetscCall(KSPGetConvergedReason(cheb->kspest, &reason));
796:       if (reason == KSP_DIVERGED_ITS) {
797:         PetscCall(PetscInfo(ksp, "Eigen estimator ran for prescribed number of iterations\n"));
798:       } else if (reason == KSP_DIVERGED_PC_FAILED) {
799:         PetscInt       its;
800:         PCFailedReason pcreason;

802:         PetscCall(KSPGetIterationNumber(cheb->kspest, &its));
803:         if (ksp->normtype == KSP_NORM_NONE) {
804:           PetscInt sendbuf, recvbuf;

806:           PetscCall(PCGetFailedReasonRank(ksp->pc, &pcreason));
807:           sendbuf = (PetscInt)pcreason;
808:           PetscCallMPI(MPI_Allreduce(&sendbuf, &recvbuf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ksp)));
809:           PetscCall(PCSetFailedReason(ksp->pc, (PCFailedReason)recvbuf));
810:         }
811:         PetscCall(PCGetFailedReason(ksp->pc, &pcreason));
812:         ksp->reason = KSP_DIVERGED_PC_FAILED;
813:         PetscCall(PetscInfo(ksp, "Eigen estimator failed: %s %s at iteration %" PetscInt_FMT, KSPConvergedReasons[reason], PCFailedReasons[pcreason], its));
814:         PetscFunctionReturn(PETSC_SUCCESS);
815:       } else if (reason == KSP_CONVERGED_RTOL || reason == KSP_CONVERGED_ATOL) {
816:         PetscCall(PetscInfo(ksp, "Eigen estimator converged prematurely. Should not happen except for small or low rank problem\n"));
817:       } else if (reason < 0) {
818:         PetscCall(PetscInfo(ksp, "Eigen estimator failed %s, using estimates anyway\n", KSPConvergedReasons[reason]));
819:       }

821:       PetscCall(KSPChebyshevComputeExtremeEigenvalues_Private(cheb->kspest, &min, &max));
822:       PetscCall(KSPSetPC(cheb->kspest, NULL));

824:       cheb->emin_computed = min;
825:       cheb->emax_computed = max;

827:       cheb->amatid    = amatid;
828:       cheb->pmatid    = pmatid;
829:       cheb->amatstate = amatstate;
830:       cheb->pmatstate = pmatstate;
831:     }
832:   }
833:   PetscFunctionReturn(PETSC_SUCCESS);
834: }

836: static PetscErrorCode KSPDestroy_Chebyshev(KSP ksp)
837: {
838:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

840:   PetscFunctionBegin;
841:   PetscCall(PetscFree(cheb->betas));
842:   PetscCall(KSPDestroy(&cheb->kspest));
843:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevSetEigenvalues_C", NULL));
844:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSet_C", NULL));
845:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSetUseNoisy_C", NULL));
846:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevSetKind_C", NULL));
847:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigGetKSP_C", NULL));
848:   PetscCall(KSPDestroyDefault(ksp));
849:   PetscFunctionReturn(PETSC_SUCCESS);
850: }

852: /*MC
853:      KSPCHEBYSHEV - The preconditioned Chebyshev iterative method

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

863:    Level: beginner

865:    Notes:
866:    The Chebyshev method requires both the matrix and preconditioner to be symmetric positive (semi) definite, but it can work as a smoother in other situations

868:    Only support for left preconditioning.

870:    Chebyshev is configured as a smoother by default, targeting the "upper" part of the spectrum.

872:    The user should call `KSPChebyshevSetEigenvalues()` to get eigenvalue estimates.

874: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`,
875:           `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`, `KSPChebyshevEstEigSetUseNoisy()`
876:           `KSPRICHARDSON`, `KSPCG`, `PCMG`
877: M*/

879: PETSC_EXTERN PetscErrorCode KSPCreate_Chebyshev(KSP ksp)
880: {
881:   KSP_Chebyshev *chebyshevP;

883:   PetscFunctionBegin;
884:   PetscCall(PetscNew(&chebyshevP));

886:   ksp->data = (void *)chebyshevP;
887:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3));
888:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 2));
889:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1));
890:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1));

892:   chebyshevP->emin = 0.;
893:   chebyshevP->emax = 0.;

895:   chebyshevP->tform[0] = 0.0;
896:   chebyshevP->tform[1] = 0.1;
897:   chebyshevP->tform[2] = 0;
898:   chebyshevP->tform[3] = 1.1;
899:   chebyshevP->eststeps = 10;
900:   chebyshevP->usenoisy = PETSC_TRUE;
901:   ksp->setupnewmatrix  = PETSC_TRUE;

903:   ksp->ops->setup          = KSPSetUp_Chebyshev;
904:   ksp->ops->destroy        = KSPDestroy_Chebyshev;
905:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
906:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
907:   ksp->ops->setfromoptions = KSPSetFromOptions_Chebyshev;
908:   ksp->ops->view           = KSPView_Chebyshev;
909:   ksp->ops->reset          = KSPReset_Chebyshev;

911:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevSetEigenvalues_C", KSPChebyshevSetEigenvalues_Chebyshev));
912:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSet_C", KSPChebyshevEstEigSet_Chebyshev));
913:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSetUseNoisy_C", KSPChebyshevEstEigSetUseNoisy_Chebyshev));
914:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevSetKind_C", KSPChebyshevSetKind_Chebyshev));
915:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigGetKSP_C", KSPChebyshevEstEigGetKSP_Chebyshev));
916:   PetscFunctionReturn(PETSC_SUCCESS);
917: }