KSPGMRESSetHapTol#

Sets the tolerance for detecting a happy ending in GMRES (KSPGMRES, KSPFGMRES and KSPLGMRES and others)

Synopsis#

#include "petscksp.h"  
PetscErrorCode KSPGMRESSetHapTol(KSP ksp, PetscReal tol)

Logically Collective

Input Parameters#

  • ksp - the Krylov space solver context

  • tol - the tolerance for detecting a happy ending

Options Database Key#

  • -ksp_gmres_haptol - set tolerance for determining happy breakdown

Note#

Happy ending is the rare case in KSPGMRES where a very near zero matrix entry is generated in the upper Hessenberg matrix indicating an ‘exact’ solution has been obtained. If you attempt more iterations after this point with GMRES unstable things can happen.

The default tolerance value for detecting a happy ending with GMRES in PETSc is 1.0e-30.

See Also#

KSP: Linear System Solvers, KSPGMRES, KSPSetTolerances()

Level#

intermediate

Location#

src/ksp/ksp/impls/gmres/gmres.c

Implementations#

KSPGMRESSetHapTol_GMRES() in src/ksp/ksp/impls/gmres/gmres.c


Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages