petsc-3.7.7 2017-09-25
KSPGMRESSetHapTol
Sets tolerance for determining happy breakdown in GMRES, FGMRES and LGMRES.
Synopsis
#include "petscksp.h"
PetscErrorCode KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
Logically Collective on KSP
Input Parameters
| ksp | - the Krylov space context
|
| tol | - the tolerance
|
Options Database
-ksp_gmres_haptol <positive real value> -
Note: Happy breakdown is the rare case in GMRES where an 'exact' solution is obtained after
a certain number of iterations. If you attempt more iterations after this point unstable
things can happen hence very occasionally you may need to set this value to detect this condition
Keywords
KSP, GMRES, tolerance
See Also
KSPSetTolerances()
Level:intermediate
Location:src/ksp/ksp/impls/gmres/gmres.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages