SNESConvergedDefault#
Default onvergence test of the solvers for systems of nonlinear equations.
Synopsis#
#include "petsc/private/snesimpl.h"
PetscErrorCode SNESConvergedDefault(SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
Collective
Input Parameters#
snes - the
SNES
contextit - the iteration (0 indicates before any Newton steps)
xnorm - 2-norm of current iterate
snorm - 2-norm of current step
fnorm - 2-norm of function at current iterate
dummy - unused context
Output Parameter#
reason - one of
SNES_CONVERGED_FNORM_ABS - (fnorm < abstol),
SNES_CONVERGED_SNORM_RELATIVE - (snorm < stol*xnorm),
SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0),
SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf),
SNES_DIVERGED_FNORM_NAN - (fnorm == NaN),
SNES_CONVERGED_ITERATING - (otherwise),
SNES_DIVERGED_DTOL - (fnorm > divtol*snes->fnorm0)
where
maxf - maximum number of function evaluations, set with
SNESSetTolerances()
nfct - number of function evaluations,
abstol - absolute function norm tolerance, set with
SNESSetTolerances()
rtol - relative function norm tolerance, set with
SNESSetTolerances()
divtol - divergence tolerance, set with
SNESSetDivergenceTolerance()
fnorm0 - 2-norm of the function at the initial solution (initial guess; zeroth iteration)
Options Database Keys#
-snes_convergence_test default - see
SNESSetFromOptions()
-snes_stol - convergence tolerance in terms of the norm of the change in the solution between steps
-snes_atol
- absolute tolerance of residual norm-snes_rtol
- relative decrease in tolerance norm from the initial 2-norm of the solution-snes_divergence_tolerance
- if the residual goes above divtol*rnorm0, exit with divergence-snes_max_funcs <max_funcs> - maximum number of function evaluations
-snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
-snes_max_linear_solve_fail - number of linear solver failures before
SNESSolve()
stops
See Also#
SNES: Nonlinear Solvers, SNES
, SNESSolve()
, SNESSetConvergenceTest()
, SNESConvergedSkip()
, SNESSetTolerances()
, SNESSetDivergenceTolerance()
,
SNESConvergedReason
Level#
intermediate
Location#
Examples#
Implementations#
SNESConvergedDefault_VI in src/snes/impls/vi/vi.c
Index of all SNES routines
Table of Contents for all manual pages
Index of all manual pages