Actual source code: cheby.c


  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;

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

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

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

 37: static PetscErrorCode KSPSetUp_Chebyshev(KSP ksp)
 38: {
 39:   KSP_Chebyshev   *cheb = (KSP_Chebyshev *)ksp->data;
 40:   PetscBool        isset, flg;
 41:   Mat              Pmat, Amat;
 42:   PetscObjectId    amatid, pmatid;
 43:   PetscObjectState amatstate, pmatstate;

 45:   PetscFunctionBegin;
 46:   PetscCall(KSPSetWorkVecs(ksp, 3));
 47:   if (cheb->emin == 0. || cheb->emax == 0.) { // User did not specify eigenvalues
 48:     PC pc;
 49:     PetscCall(KSPGetPC(ksp, &pc));
 50:     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &flg));
 51:     if (!flg) { // Provided estimates are only relevant for Jacobi
 52:       cheb->emax_provided = 0;
 53:       cheb->emin_provided = 0;
 54:     }
 55:     if (!cheb->kspest) { /* We need to estimate eigenvalues */
 56:       PetscCall(KSPChebyshevEstEigSet(ksp, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
 57:     }
 58:   }
 59:   if (cheb->kspest) {
 60:     PetscCall(KSPGetOperators(ksp, &Amat, &Pmat));
 61:     PetscCall(MatIsSPDKnown(Pmat, &isset, &flg));
 62:     if (isset && flg) {
 63:       const char *prefix;
 64:       PetscCall(KSPGetOptionsPrefix(cheb->kspest, &prefix));
 65:       PetscCall(PetscOptionsHasName(NULL, prefix, "-ksp_type", &flg));
 66:       if (!flg) PetscCall(KSPSetType(cheb->kspest, KSPCG));
 67:     }
 68:     PetscCall(PetscObjectGetId((PetscObject)Amat, &amatid));
 69:     PetscCall(PetscObjectGetId((PetscObject)Pmat, &pmatid));
 70:     PetscCall(PetscObjectStateGet((PetscObject)Amat, &amatstate));
 71:     PetscCall(PetscObjectStateGet((PetscObject)Pmat, &pmatstate));
 72:     if (amatid != cheb->amatid || pmatid != cheb->pmatid || amatstate != cheb->amatstate || pmatstate != cheb->pmatstate) {
 73:       PetscReal          max = 0.0, min = 0.0;
 74:       Vec                B;
 75:       KSPConvergedReason reason;
 76:       PetscCall(KSPSetPC(cheb->kspest, ksp->pc));
 77:       if (cheb->usenoisy) {
 78:         B = ksp->work[1];
 79:         PetscCall(KSPSetNoisy_Private(B));
 80:       } else {
 81:         PetscBool change;

 83:         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");
 84:         PetscCall(PCPreSolveChangeRHS(ksp->pc, &change));
 85:         if (change) {
 86:           B = ksp->work[1];
 87:           PetscCall(VecCopy(ksp->vec_rhs, B));
 88:         } else B = ksp->vec_rhs;
 89:       }
 90:       PetscCall(KSPSolve(cheb->kspest, B, ksp->work[0]));
 91:       PetscCall(KSPGetConvergedReason(cheb->kspest, &reason));
 92:       if (reason == KSP_DIVERGED_ITS) {
 93:         PetscCall(PetscInfo(ksp, "Eigen estimator ran for prescribed number of iterations\n"));
 94:       } else if (reason == KSP_DIVERGED_PC_FAILED) {
 95:         PetscInt       its;
 96:         PCFailedReason pcreason;

 98:         PetscCall(KSPGetIterationNumber(cheb->kspest, &its));
 99:         if (ksp->normtype == KSP_NORM_NONE) {
100:           PetscInt sendbuf, recvbuf;
101:           PetscCall(PCGetFailedReasonRank(ksp->pc, &pcreason));
102:           sendbuf = (PetscInt)pcreason;
103:           PetscCallMPI(MPI_Allreduce(&sendbuf, &recvbuf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ksp)));
104:           PetscCall(PCSetFailedReason(ksp->pc, (PCFailedReason)recvbuf));
105:         }
106:         PetscCall(PCGetFailedReason(ksp->pc, &pcreason));
107:         ksp->reason = KSP_DIVERGED_PC_FAILED;
108:         PetscCall(PetscInfo(ksp, "Eigen estimator failed: %s %s at iteration %" PetscInt_FMT, KSPConvergedReasons[reason], PCFailedReasons[pcreason], its));
109:         PetscFunctionReturn(PETSC_SUCCESS);
110:       } else if (reason == KSP_CONVERGED_RTOL || reason == KSP_CONVERGED_ATOL) {
111:         PetscCall(PetscInfo(ksp, "Eigen estimator converged prematurely. Should not happen except for small or low rank problem\n"));
112:       } else if (reason < 0) {
113:         PetscCall(PetscInfo(ksp, "Eigen estimator failed %s, using estimates anyway\n", KSPConvergedReasons[reason]));
114:       }

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

119:       cheb->emin_computed = min;
120:       cheb->emax_computed = max;

122:       cheb->amatid    = amatid;
123:       cheb->pmatid    = pmatid;
124:       cheb->amatstate = amatstate;
125:       cheb->pmatstate = pmatstate;
126:     }
127:   }
128:   PetscFunctionReturn(PETSC_SUCCESS);
129: }

131: static PetscErrorCode KSPChebyshevGetEigenvalues_Chebyshev(KSP ksp, PetscReal *emax, PetscReal *emin)
132: {
133:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

135:   PetscFunctionBegin;
136:   *emax = 0;
137:   *emin = 0;
138:   if (cheb->emax != 0.) {
139:     *emax = cheb->emax;
140:   } else if (cheb->emax_computed != 0.) {
141:     *emax = cheb->tform[2] * cheb->emin_computed + cheb->tform[3] * cheb->emax_computed;
142:   } else if (cheb->emax_provided != 0.) {
143:     *emax = cheb->tform[2] * cheb->emin_provided + cheb->tform[3] * cheb->emax_provided;
144:   }
145:   if (cheb->emin != 0.) {
146:     *emin = cheb->emin;
147:   } else if (cheb->emin_computed != 0.) {
148:     *emin = cheb->tform[0] * cheb->emin_computed + cheb->tform[1] * cheb->emax_computed;
149:   } else if (cheb->emin_provided != 0.) {
150:     *emin = cheb->tform[0] * cheb->emin_provided + cheb->tform[1] * cheb->emax_provided;
151:   }
152:   PetscFunctionReturn(PETSC_SUCCESS);
153: }

155: static PetscErrorCode KSPChebyshevSetEigenvalues_Chebyshev(KSP ksp, PetscReal emax, PetscReal emin)
156: {
157:   KSP_Chebyshev *chebyshevP = (KSP_Chebyshev *)ksp->data;

159:   PetscFunctionBegin;
160:   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);
161:   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);
162:   chebyshevP->emax = emax;
163:   chebyshevP->emin = emin;

165:   PetscCall(KSPChebyshevEstEigSet(ksp, 0., 0., 0., 0.)); /* Destroy any estimation setup */
166:   PetscFunctionReturn(PETSC_SUCCESS);
167: }

169: static PetscErrorCode KSPChebyshevEstEigSet_Chebyshev(KSP ksp, PetscReal a, PetscReal b, PetscReal c, PetscReal d)
170: {
171:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

173:   PetscFunctionBegin;
174:   if (a != 0.0 || b != 0.0 || c != 0.0 || d != 0.0) {
175:     if ((cheb->emin_provided == 0. || cheb->emax_provided == 0.) && !cheb->kspest) { /* should this block of code be moved to KSPSetUp_Chebyshev()? */
176:       PetscCall(KSPCreate(PetscObjectComm((PetscObject)ksp), &cheb->kspest));
177:       PetscCall(KSPSetErrorIfNotConverged(cheb->kspest, ksp->errorifnotconverged));
178:       PetscCall(PetscObjectIncrementTabLevel((PetscObject)cheb->kspest, (PetscObject)ksp, 1));
179:       /* use PetscObjectSet/AppendOptionsPrefix() instead of KSPSet/AppendOptionsPrefix() so that the PC prefix is not changed */
180:       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)cheb->kspest, ((PetscObject)ksp)->prefix));
181:       PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)cheb->kspest, "esteig_"));
182:       PetscCall(KSPSetSkipPCSetFromOptions(cheb->kspest, PETSC_TRUE));

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

186:       /* We cannot turn off convergence testing because GMRES will break down if you attempt to keep iterating after a zero norm is obtained */
187:       PetscCall(KSPSetTolerances(cheb->kspest, 1.e-12, PETSC_DEFAULT, PETSC_DEFAULT, cheb->eststeps));
188:     }
189:     if (a >= 0) cheb->tform[0] = a;
190:     if (b >= 0) cheb->tform[1] = b;
191:     if (c >= 0) cheb->tform[2] = c;
192:     if (d >= 0) cheb->tform[3] = d;
193:     cheb->amatid    = 0;
194:     cheb->pmatid    = 0;
195:     cheb->amatstate = -1;
196:     cheb->pmatstate = -1;
197:   } else {
198:     PetscCall(KSPDestroy(&cheb->kspest));
199:   }
200:   PetscFunctionReturn(PETSC_SUCCESS);
201: }

203: static PetscErrorCode KSPChebyshevEstEigSetUseNoisy_Chebyshev(KSP ksp, PetscBool use)
204: {
205:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

207:   PetscFunctionBegin;
208:   cheb->usenoisy = use;
209:   PetscFunctionReturn(PETSC_SUCCESS);
210: }

212: /*@
213:    KSPChebyshevSetEigenvalues - Sets estimates for the extreme eigenvalues
214:    of the preconditioned problem.

216:    Logically Collective

218:    Input Parameters:
219: +  ksp - the Krylov space context
220: -  emax, emin - the eigenvalue estimates

222:   Options Database Key:
223: .  -ksp_chebyshev_eigenvalues emin,emax - extreme eigenvalues

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

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

234:    Level: intermediate

236: .seealso: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevEstEigSet()`,
237: @*/
238: PetscErrorCode KSPChebyshevSetEigenvalues(KSP ksp, PetscReal emax, PetscReal emin)
239: {
240:   PetscFunctionBegin;
244:   PetscTryMethod(ksp, "KSPChebyshevSetEigenvalues_C", (KSP, PetscReal, PetscReal), (ksp, emax, emin));
245:   PetscFunctionReturn(PETSC_SUCCESS);
246: }

248: /*@
249:    KSPChebyshevEstEigSet - Automatically estimate the eigenvalues to use for Chebyshev

251:    Logically Collective

253:    Input Parameters:
254: +  ksp - the Krylov space context
255: .  a - multiple of min eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE)
256: .  b - multiple of max eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE)
257: .  c - multiple of min eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE)
258: -  d - multiple of max eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE)

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

263:    Notes:
264:    The Chebyshev bounds are set using
265: .vb
266:    minbound = a*minest + b*maxest
267:    maxbound = c*minest + d*maxest
268: .ve
269:    The default configuration targets the upper part of the spectrum for use as a multigrid smoother, so only the maximum eigenvalue estimate is used.
270:    The minimum eigenvalue estimate obtained by Krylov iteration is typically not accurate until the method has converged.

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

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

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

278:    Level: intermediate

280: .seealso: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevEstEigSet()`, `KSPChebyshevEstEigSetUseNoisy()`, `KSPChebyshevEstEigGetKSP()`
281: @*/
282: PetscErrorCode KSPChebyshevEstEigSet(KSP ksp, PetscReal a, PetscReal b, PetscReal c, PetscReal d)
283: {
284:   PetscFunctionBegin;
290:   PetscTryMethod(ksp, "KSPChebyshevEstEigSet_C", (KSP, PetscReal, PetscReal, PetscReal, PetscReal), (ksp, a, b, c, d));
291:   PetscFunctionReturn(PETSC_SUCCESS);
292: }

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

297:    Logically Collective

299:    Input Parameters:
300: +  ksp - linear solver context
301: -  use - `PETSC_TRUE` to use noisy

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

306:    Note:
307:     This allegedly works better for multigrid smoothers

309:   Level: intermediate

311: .seealso: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevEstEigSet()`, `KSPChebyshevEstEigGetKSP()`
312: @*/
313: PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP ksp, PetscBool use)
314: {
315:   PetscFunctionBegin;
317:   PetscTryMethod(ksp, "KSPChebyshevEstEigSetUseNoisy_C", (KSP, PetscBool), (ksp, use));
318:   PetscFunctionReturn(PETSC_SUCCESS);
319: }

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

326:   Input Parameters:
327: . ksp - the Krylov space context

329:   Output Parameters:
330: . kspest - the eigenvalue estimation Krylov space context

332:   Level: advanced

334: .seealso: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevEstEigSet()`
335: @*/
336: PetscErrorCode KSPChebyshevEstEigGetKSP(KSP ksp, KSP *kspest)
337: {
338:   PetscFunctionBegin;
341:   *kspest = NULL;
342:   PetscTryMethod(ksp, "KSPChebyshevEstEigGetKSP_C", (KSP, KSP *), (ksp, kspest));
343:   PetscFunctionReturn(PETSC_SUCCESS);
344: }

346: static PetscErrorCode KSPChebyshevEstEigGetKSP_Chebyshev(KSP ksp, KSP *kspest)
347: {
348:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

350:   PetscFunctionBegin;
351:   *kspest = cheb->kspest;
352:   PetscFunctionReturn(PETSC_SUCCESS);
353: }

355: static PetscErrorCode KSPSetFromOptions_Chebyshev(KSP ksp, PetscOptionItems *PetscOptionsObject)
356: {
357:   KSP_Chebyshev *cheb    = (KSP_Chebyshev *)ksp->data;
358:   PetscInt       neigarg = 2, nestarg = 4;
359:   PetscReal      eminmax[2] = {0., 0.};
360:   PetscReal      tform[4]   = {PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE};
361:   PetscBool      flgeig, flgest;

363:   PetscFunctionBegin;
364:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP Chebyshev Options");
365:   PetscCall(PetscOptionsInt("-ksp_chebyshev_esteig_steps", "Number of est steps in Chebyshev", "", cheb->eststeps, &cheb->eststeps, NULL));
366:   PetscCall(PetscOptionsRealArray("-ksp_chebyshev_eigenvalues", "extreme eigenvalues", "KSPChebyshevSetEigenvalues", eminmax, &neigarg, &flgeig));
367:   if (flgeig) {
368:     PetscCheck(neigarg == 2, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "-ksp_chebyshev_eigenvalues: must specify 2 parameters, min and max eigenvalues");
369:     PetscCall(KSPChebyshevSetEigenvalues(ksp, eminmax[1], eminmax[0]));
370:   }
371:   PetscCall(PetscOptionsRealArray("-ksp_chebyshev_esteig", "estimate eigenvalues using a Krylov method, then use this transform for Chebyshev eigenvalue bounds", "KSPChebyshevEstEigSet", tform, &nestarg, &flgest));
372:   if (flgest) {
373:     switch (nestarg) {
374:     case 0:
375:       PetscCall(KSPChebyshevEstEigSet(ksp, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
376:       break;
377:     case 2: /* Base everything on the max eigenvalues */
378:       PetscCall(KSPChebyshevEstEigSet(ksp, PETSC_DECIDE, tform[0], PETSC_DECIDE, tform[1]));
379:       break;
380:     case 4: /* Use the full 2x2 linear transformation */
381:       PetscCall(KSPChebyshevEstEigSet(ksp, tform[0], tform[1], tform[2], tform[3]));
382:       break;
383:     default:
384:       SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "Must specify either 0, 2, or 4 parameters for eigenvalue estimation");
385:     }
386:   }

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

391:   if (cheb->kspest) {
392:     PetscCall(PetscOptionsBool("-ksp_chebyshev_esteig_noisy", "Use noisy right hand side for estimate", "KSPChebyshevEstEigSetUseNoisy", cheb->usenoisy, &cheb->usenoisy, NULL));
393:     PetscCall(KSPSetFromOptions(cheb->kspest));
394:   }
395:   PetscOptionsHeadEnd();
396:   PetscFunctionReturn(PETSC_SUCCESS);
397: }

399: static PetscErrorCode KSPSolve_Chebyshev(KSP ksp)
400: {
401:   PetscInt    k, kp1, km1, ktmp, i;
402:   PetscScalar alpha, omegaprod, mu, omega, Gamma, c[3], scale;
403:   PetscReal   rnorm = 0.0, emax, emin;
404:   Vec         sol_orig, b, p[3], r;
405:   Mat         Amat, Pmat;
406:   PetscBool   diagonalscale;

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

412:   PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
413:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
414:   ksp->its = 0;
415:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
416:   /* These three point to the three active solutions, we
417:      rotate these three at each solution update */
418:   km1      = 0;
419:   k        = 1;
420:   kp1      = 2;
421:   sol_orig = ksp->vec_sol; /* ksp->vec_sol will be assigned to rotating vector p[k], thus save its address */
422:   b        = ksp->vec_rhs;
423:   p[km1]   = sol_orig;
424:   p[k]     = ksp->work[0];
425:   p[kp1]   = ksp->work[1];
426:   r        = ksp->work[2];

428:   PetscCall(KSPChebyshevGetEigenvalues_Chebyshev(ksp, &emax, &emin));
429:   /* use scale*B as our preconditioner */
430:   scale = 2.0 / (emax + emin);

432:   /*   -alpha <=  scale*lambda(B^{-1}A) <= alpha   */
433:   alpha     = 1.0 - scale * emin;
434:   Gamma     = 1.0;
435:   mu        = 1.0 / alpha;
436:   omegaprod = 2.0 / alpha;

438:   c[km1] = 1.0;
439:   c[k]   = mu;

441:   if (!ksp->guess_zero) {
442:     PetscCall(KSP_MatMult(ksp, Amat, sol_orig, r)); /*  r = b - A*p[km1] */
443:     PetscCall(VecAYPX(r, -1.0, b));
444:   } else {
445:     PetscCall(VecCopy(b, r));
446:   }

448:   /* calculate residual norm if requested, we have done one iteration */
449:   if (ksp->normtype) {
450:     switch (ksp->normtype) {
451:     case KSP_NORM_PRECONDITIONED:
452:       PetscCall(KSP_PCApply(ksp, r, p[k])); /* p[k] = B^{-1}r */
453:       PetscCall(VecNorm(p[k], NORM_2, &rnorm));
454:       break;
455:     case KSP_NORM_UNPRECONDITIONED:
456:     case KSP_NORM_NATURAL:
457:       PetscCall(VecNorm(r, NORM_2, &rnorm));
458:       break;
459:     default:
460:       SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
461:     }
462:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
463:     ksp->rnorm = rnorm;
464:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
465:     PetscCall(KSPLogResidualHistory(ksp, rnorm));
466:     PetscCall(KSPLogErrorHistory(ksp));
467:     PetscCall(KSPMonitor(ksp, 0, rnorm));
468:     PetscCall((*ksp->converged)(ksp, 0, rnorm, &ksp->reason, ksp->cnvP));
469:   } else ksp->reason = KSP_CONVERGED_ITERATING;
470:   if (ksp->reason || ksp->max_it == 0) {
471:     if (ksp->max_it == 0) ksp->reason = KSP_DIVERGED_ITS; /* This for a V(0,x) cycle */
472:     PetscFunctionReturn(PETSC_SUCCESS);
473:   }
474:   if (ksp->normtype != KSP_NORM_PRECONDITIONED) { PetscCall(KSP_PCApply(ksp, r, p[k])); /* p[k] = B^{-1}r */ }
475:   PetscCall(VecAYPX(p[k], scale, p[km1])); /* p[k] = scale B^{-1}r + p[km1] */
476:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
477:   ksp->its = 1;
478:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));

480:   for (i = 1; i < ksp->max_it; i++) {
481:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
482:     ksp->its++;
483:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));

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

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

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

523:     ktmp = km1;
524:     km1  = k;
525:     k    = kp1;
526:     kp1  = ktmp;
527:   }
528:   if (!ksp->reason) {
529:     if (ksp->normtype) {
530:       PetscCall(KSP_MatMult(ksp, Amat, p[k], r)); /*  r = b - Ap[k]    */
531:       PetscCall(VecAYPX(r, -1.0, b));
532:       switch (ksp->normtype) {
533:       case KSP_NORM_PRECONDITIONED:
534:         PetscCall(KSP_PCApply(ksp, r, p[kp1])); /* p[kp1] = B^{-1}r */
535:         PetscCall(VecNorm(p[kp1], NORM_2, &rnorm));
536:         break;
537:       case KSP_NORM_UNPRECONDITIONED:
538:       case KSP_NORM_NATURAL:
539:         PetscCall(VecNorm(r, NORM_2, &rnorm));
540:         break;
541:       default:
542:         rnorm = 0.0;
543:         break;
544:       }
545:       KSPCheckNorm(ksp, rnorm);
546:       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
547:       ksp->rnorm = rnorm;
548:       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
549:       PetscCall(KSPLogResidualHistory(ksp, rnorm));
550:       PetscCall(KSPMonitor(ksp, i, rnorm));
551:     }
552:     if (ksp->its >= ksp->max_it) {
553:       if (ksp->normtype != KSP_NORM_NONE) {
554:         PetscCall((*ksp->converged)(ksp, i, rnorm, &ksp->reason, ksp->cnvP));
555:         if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
556:       } else ksp->reason = KSP_CONVERGED_ITS;
557:     }
558:   }

560:   /* make sure solution is in vector x */
561:   ksp->vec_sol = sol_orig;
562:   if (k) PetscCall(VecCopy(p[k], sol_orig));
563:   if (ksp->reason == KSP_CONVERGED_ITS) PetscCall(KSPLogErrorHistory(ksp));
564:   PetscFunctionReturn(PETSC_SUCCESS);
565: }

567: static PetscErrorCode KSPView_Chebyshev(KSP ksp, PetscViewer viewer)
568: {
569:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
570:   PetscBool      iascii;

572:   PetscFunctionBegin;
573:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
574:   if (iascii) {
575:     PetscReal emax, emin;
576:     PetscCall(KSPChebyshevGetEigenvalues_Chebyshev(ksp, &emax, &emin));
577:     PetscCall(PetscViewerASCIIPrintf(viewer, "  eigenvalue targets used: min %g, max %g\n", (double)emin, (double)emax));
578:     if (cheb->kspest) {
579:       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));
580:       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]));
581:       PetscCall(PetscViewerASCIIPushTab(viewer));
582:       PetscCall(KSPView(cheb->kspest, viewer));
583:       PetscCall(PetscViewerASCIIPopTab(viewer));
584:       if (cheb->usenoisy) PetscCall(PetscViewerASCIIPrintf(viewer, "  estimating eigenvalues using noisy right hand side\n"));
585:     } else if (cheb->emax_provided != 0.) {
586:       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],
587:                                        (double)cheb->tform[3]));
588:     }
589:   }
590:   PetscFunctionReturn(PETSC_SUCCESS);
591: }

593: static PetscErrorCode KSPDestroy_Chebyshev(KSP ksp)
594: {
595:   KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;

597:   PetscFunctionBegin;
598:   PetscCall(KSPDestroy(&cheb->kspest));
599:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevSetEigenvalues_C", NULL));
600:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSet_C", NULL));
601:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSetUseNoisy_C", NULL));
602:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigGetKSP_C", NULL));
603:   PetscCall(KSPDestroyDefault(ksp));
604:   PetscFunctionReturn(PETSC_SUCCESS);
605: }

607: /*MC
608:      KSPCHEBYSHEV - The preconditioned Chebyshev iterative method

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

618:    Level: beginner

620:    Notes:
621:    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

623:    Only support for left preconditioning.

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

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

629: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`,
630:           `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`, `KSPChebyshevEstEigSetUseNoisy()`
631:           `KSPRICHARDSON`, `KSPCG`, `PCMG`
632: M*/

634: PETSC_EXTERN PetscErrorCode KSPCreate_Chebyshev(KSP ksp)
635: {
636:   KSP_Chebyshev *chebyshevP;

638:   PetscFunctionBegin;
639:   PetscCall(PetscNew(&chebyshevP));

641:   ksp->data = (void *)chebyshevP;
642:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3));
643:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 2));
644:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1));
645:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1));

647:   chebyshevP->emin = 0.;
648:   chebyshevP->emax = 0.;

650:   chebyshevP->tform[0] = 0.0;
651:   chebyshevP->tform[1] = 0.1;
652:   chebyshevP->tform[2] = 0;
653:   chebyshevP->tform[3] = 1.1;
654:   chebyshevP->eststeps = 10;
655:   chebyshevP->usenoisy = PETSC_TRUE;
656:   ksp->setupnewmatrix  = PETSC_TRUE;

658:   ksp->ops->setup          = KSPSetUp_Chebyshev;
659:   ksp->ops->solve          = KSPSolve_Chebyshev;
660:   ksp->ops->destroy        = KSPDestroy_Chebyshev;
661:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
662:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
663:   ksp->ops->setfromoptions = KSPSetFromOptions_Chebyshev;
664:   ksp->ops->view           = KSPView_Chebyshev;
665:   ksp->ops->reset          = KSPReset_Chebyshev;

667:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevSetEigenvalues_C", KSPChebyshevSetEigenvalues_Chebyshev));
668:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSet_C", KSPChebyshevEstEigSet_Chebyshev));
669:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSetUseNoisy_C", KSPChebyshevEstEigSetUseNoisy_Chebyshev));
670:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigGetKSP_C", KSPChebyshevEstEigGetKSP_Chebyshev));
671:   PetscFunctionReturn(PETSC_SUCCESS);
672: }