NonlinearSolve.jl Solvers
These are the native solvers of NonlinearSolve.jl.
NonlinearSolveBase.NonlinearSolvePolyAlgorithm
NonlinearSolveFirstOrder.GeneralizedFirstOrderAlgorithm
NonlinearSolveQuasiNewton.QuasiNewtonAlgorithm
NonlinearSolveSpectralMethods.GeneralizedDFSane
NonlinearSolve.FastShortcutNLLSPolyalg
NonlinearSolve.FastShortcutNonlinearPolyalg
NonlinearSolveFirstOrder.GaussNewton
NonlinearSolveFirstOrder.LevenbergMarquardt
NonlinearSolveFirstOrder.NewtonRaphson
NonlinearSolveFirstOrder.PseudoTransient
NonlinearSolveFirstOrder.RobustMultiNewton
NonlinearSolveFirstOrder.TrustRegion
NonlinearSolveQuasiNewton.Broyden
NonlinearSolveQuasiNewton.Klement
NonlinearSolveQuasiNewton.LimitedMemoryBroyden
NonlinearSolveSpectralMethods.DFSane
General Keyword Arguments
Several Algorithms share the same specification for common keyword arguments. Those are documented in this section to avoid repetition. Certain algorithms might have additional considerations for these keyword arguments, which are documented in the algorithm's documentation.
linsolve
: the LinearSolve.jl solvers used for the linear solves within the Newton method. Defaults tonothing
, which means it uses the LinearSolve.jl default algorithm choice. For more information on available algorithm choices, see the LinearSolve.jl documentation.linesearch
: the line search algorithm to use. Defaults toNoLineSearch()
, which means that no line search is performed.autodiff
: determines the backend used for the Jacobian. Note that this argument is ignored if an analytical Jacobian is passed, as that will be used instead. Defaults tonothing
which means that a default is selected according to the problem specification! Valid choices are types from ADTypes.jl.vjp_autodiff
: similar toautodiff
, but is used to compute Jacobian Vector Products. Ignored if the NonlinearFunction contains thejvp
function.vjp_autodiff
: similar toautodiff
, but is used to compute Vector Jacobian Products. Ignored if the NonlinearFunction contains thevjp
function.concrete_jac
: whether to build a concrete Jacobian. If a Krylov-subspace method is used, then the Jacobian will not be constructed and instead direct Jacobian-Vector productsJ*v
are computed using forward-mode automatic differentiation or finite differencing tricks (without ever constructing the Jacobian). However, if the Jacobian is still needed, for example for a preconditioner,concrete_jac = true
can be passed in order to force the construction of the Jacobian.
Nonlinear Solvers
NonlinearSolveFirstOrder.NewtonRaphson
— FunctionNewtonRaphson(;
concrete_jac = nothing, linsolve = nothing, linesearch = missing,
autodiff = nothing, vjp_autodiff = nothing, jvp_autodiff = nothing
)
An advanced NewtonRaphson implementation with support for efficient handling of sparse matrices via colored automatic differentiation and preconditioned linear solvers. Designed for large-scale and numerically-difficult nonlinear systems.
NonlinearSolveSpectralMethods.DFSane
— FunctionDFSane(;
sigma_min = 1 // 10^10, sigma_max = 1e10, sigma_1 = 1, M::Int = 10,
gamma = 1 // 10^4, tau_min = 1 // 10, tau_max = 1 // 2, n_exp::Int = 2,
max_inner_iterations::Int = 100, eta_strategy = (fn_1, n, x_n, f_n) -> fn_1 / n^2
)
A low-overhead and allocation-free implementation of the df-sane method for solving large-scale nonlinear systems of equations. For in depth information about all the parameters and the algorithm, see La Cruz et al. [2].
Keyword Arguments
sigma_min
: the minimum value of the spectral coefficientσ
which is related to the step size in the algorithm. Defaults to1e-10
.sigma_max
: the maximum value of the spectral coefficientσₙ
which is related to the step size in the algorithm. Defaults to1e10
.
For other keyword arguments, see RobustNonMonotoneLineSearch
.
NonlinearSolveQuasiNewton.Broyden
— FunctionBroyden(;
max_resets::Int = 100, linesearch = nothing, reset_tolerance = nothing,
init_jacobian::Val = Val(:identity), autodiff = nothing, alpha = nothing,
update_rule = Val(:good_broyden)
)
An implementation of Broyden
's Method [3] with resetting and line search.
Keyword Arguments
max_resets
: the maximum number of resets to perform. Defaults to100
.reset_tolerance
: the tolerance for the reset check. Defaults tosqrt(eps(real(eltype(u))))
.alpha
: Ifinit_jacobian
is set toVal(:identity)
, then the initial Jacobian inverse is set to be(αI)⁻¹
. Defaults tonothing
which impliesα = max(norm(u), 1) / (2 * norm(fu))
.init_jacobian
: the method to use for initializing the jacobian. Defaults toVal(:identity)
. Choices include:Val(:identity)
: Identity Matrix.Val(:true_jacobian)
: True Jacobian. This is a good choice for differentiable problems.
update_rule
: Update Rule for the Jacobian. Choices are:Val(:good_broyden)
: Good Broyden's Update RuleVal(:bad_broyden)
: Bad Broyden's Update RuleVal(:diagonal)
: Only update the diagonal of the Jacobian. This algorithm may be useful for specific problems, but whether it will work may depend strongly on the problem
NonlinearSolveQuasiNewton.Klement
— FunctionKlement(;
max_resets = 100, linsolve = nothing, linesearch = nothing,
alpha = nothing, init_jacobian::Val = Val(:identity),
autodiff = nothing
)
An implementation of Klement
[4] with line search, preconditioning and customizable linear solves. It is recommended to use Broyden
for most problems over this.
Keyword Arguments
max_resets
: the maximum number of resets to perform. Defaults to100
.alpha
: Ifinit_jacobian
is set toVal(:identity)
, then the initial Jacobian inverse is set to beαI
. Defaults to1
. Can be set tonothing
which impliesα = max(norm(u), 1) / (2 * norm(fu))
.init_jacobian
: the method to use for initializing the jacobian. Defaults toVal(:identity)
. Choices include:Val(:identity)
: Identity Matrix.Val(:true_jacobian)
: True Jacobian. Our tests suggest that this is not very stable. Instead usingBroyden
withVal(:true_jacobian)
gives faster and more reliable convergence.Val(:true_jacobian_diagonal)
: Diagonal of True Jacobian. This is a good choice for differentiable problems.
NonlinearSolveQuasiNewton.LimitedMemoryBroyden
— FunctionLimitedMemoryBroyden(;
max_resets::Int = 3, linesearch = nothing, threshold::Val = Val(10),
reset_tolerance = nothing, alpha = nothing
)
An implementation of LimitedMemoryBroyden
[5] with resetting and line search.
Keyword Arguments
max_resets
: the maximum number of resets to perform. Defaults to3
.reset_tolerance
: the tolerance for the reset check. Defaults tosqrt(eps(real(eltype(u))))
.threshold
: the number of vectors to store in the low rank approximation. Defaults toVal(10)
.alpha
: The initial Jacobian inverse is set to be(αI)⁻¹
. Defaults tonothing
which impliesα = max(norm(u), 1) / (2 * norm(fu))
.
Nonlinear Least Squares Solvers
NonlinearSolveFirstOrder.GaussNewton
— FunctionGaussNewton(;
concrete_jac = nothing, linsolve = nothing, linesearch = missing,
autodiff = nothing, vjp_autodiff = nothing, jvp_autodiff = nothing
)
An advanced GaussNewton implementation with support for efficient handling of sparse matrices via colored automatic differentiation and preconditioned linear solvers. Designed for large-scale and numerically-difficult nonlinear systems.
Both Nonlinear & Nonlinear Least Squares Solvers
These solvers can be used for both nonlinear and nonlinear least squares problems.
NonlinearSolveFirstOrder.TrustRegion
— FunctionTrustRegion(;
concrete_jac = nothing, linsolve = nothing,
radius_update_scheme = RadiusUpdateSchemes.Simple, max_trust_radius::Real = 0 // 1,
initial_trust_radius::Real = 0 // 1, step_threshold::Real = 1 // 10000,
shrink_threshold::Real = 1 // 4, expand_threshold::Real = 3 // 4,
shrink_factor::Real = 1 // 4, expand_factor::Real = 2 // 1,
max_shrink_times::Int = 32,
vjp_autodiff = nothing, autodiff = nothing, jvp_autodiff = nothing
)
An advanced TrustRegion implementation with support for efficient handling of sparse matrices via colored automatic differentiation and preconditioned linear solvers. Designed for large-scale and numerically-difficult nonlinear systems.
Keyword Arguments
radius_update_scheme
: the scheme used to update the trust region radius. Defaults toRadiusUpdateSchemes.Simple
. SeeRadiusUpdateSchemes
for more details. For a review on trust region radius update schemes, see Yuan [6].
For the remaining arguments, see NonlinearSolveFirstOrder.GenericTrustRegionScheme
documentation.
NonlinearSolveFirstOrder.LevenbergMarquardt
— FunctionLevenbergMarquardt(;
linsolve = nothing,
damping_initial::Real = 1.0, α_geodesic::Real = 0.75, disable_geodesic = Val(false),
damping_increase_factor::Real = 2.0, damping_decrease_factor::Real = 3.0,
finite_diff_step_geodesic = 0.1, b_uphill::Real = 1.0, min_damping_D::Real = 1e-8,
autodiff = nothing, vjp_autodiff = nothing, jvp_autodiff = nothing
)
An advanced Levenberg-Marquardt implementation with the improvements suggested in Transtrum and Sethna [1]. Designed for large-scale and numerically-difficult nonlinear systems.
Keyword Arguments
damping_initial
: the starting value for the damping factor. The damping factor is inversely proportional to the step size. The damping factor is adjusted during each iteration. Defaults to1.0
. See Section 2.1 of Transtrum and Sethna [1].damping_increase_factor
: the factor by which the damping is increased if a step is rejected. Defaults to2.0
. See Section 2.1 of Transtrum and Sethna [1].damping_decrease_factor
: the factor by which the damping is decreased if a step is accepted. Defaults to3.0
. See Section 2.1 of Transtrum and Sethna [1].min_damping_D
: the minimum value of the damping terms in the diagonal damping matrixDᵀD
, whereDᵀD
is given by the largest diagonal entries ofJᵀJ
yet encountered, whereJ
is the Jacobian. It is suggested by Transtrum and Sethna [1] to use a minimum value of the elements inDᵀD
to prevent the damping from being too small. Defaults to1e-8
.disable_geodesic
: Disables Geodesic Acceleration if set toVal(true)
. It provides a way to trade-off robustness for speed, though in most situations Geodesic Acceleration should not be disabled.
For the remaining arguments, see GeodesicAcceleration
and NonlinearSolveFirstOrder.LevenbergMarquardtTrustRegion
documentations.
NonlinearSolveFirstOrder.PseudoTransient
— FunctionPseudoTransient(;
concrete_jac = nothing, linesearch = missing, alpha_initial = 1e-3,
linsolve = nothing,
autodiff = nothing, jvp_autodiff = nothing, vjp_autodiff = nothing
)
An implementation of PseudoTransient Method [7] that is used to solve steady state problems in an accelerated manner. It uses an adaptive time-stepping to integrate an initial value of nonlinear problem until sufficient accuracy in the desired steady-state is achieved to switch over to Newton's method and gain a rapid convergence. This implementation specifically uses "switched evolution relaxation" [8] SER method.
Keyword Arguments
alpha_initial
: the initial pseudo time step. It defaults to1e-3
. If it is small, you are going to need more iterations to converge but it can be more stable.
Polyalgorithms
NonlinearSolveBase.NonlinearSolvePolyAlgorithm
— TypeNonlinearSolvePolyAlgorithm(algs; start_index::Int = 1)
A general way to define PolyAlgorithms for NonlinearProblem
and NonlinearLeastSquaresProblem
. This is a container for a tuple of algorithms that will be tried in order until one succeeds. If none succeed, then the algorithm with the lowest residual is returned.
Arguments
algs
: a tuple of algorithms to try in-order! (If this is not a Tuple, then the returned algorithm is not type-stable).
Keyword Arguments
start_index
: the index to start at. Defaults to1
.
Example
using NonlinearSolve
alg = NonlinearSolvePolyAlgorithm((NewtonRaphson(), Broyden()))
NonlinearSolve.FastShortcutNonlinearPolyalg
— FunctionFastShortcutNonlinearPolyalg(
::Type{T} = Float64;
concrete_jac = nothing,
linsolve = nothing,
must_use_jacobian::Val = Val(false),
prefer_simplenonlinearsolve::Val = Val(false),
autodiff = nothing, vjp_autodiff = nothing, jvp_autodiff = nothing,
u0_len::Union{Int, Nothing} = nothing
) where {T}
A polyalgorithm focused on balancing speed and robustness. It first tries less robust methods for more performance and then tries more robust techniques if the faster ones fail.
Arguments
T
: The eltype of the initial guess. It is only used to check if some of the algorithms are compatible with the problem type. Defaults toFloat64
.
Keyword Arguments
u0_len
: The length of the initial guess. If this isnothing
, then the length of the initial guess is not checked. If this is an integer and it is less than25
, we use jacobian based methods.
NonlinearSolve.FastShortcutNLLSPolyalg
— FunctionFastShortcutNLLSPolyalg(
::Type{T} = Float64;
concrete_jac = nothing,
linsolve = nothing,
autodiff = nothing, vjp_autodiff = nothing, jvp_autodiff = nothing
)
A polyalgorithm focused on balancing speed and robustness. It first tries less robust methods for more performance and then tries more robust techniques if the faster ones fail.
Arguments
T
: The eltype of the initial guess. It is only used to check if some of the algorithms are compatible with the problem type. Defaults toFloat64
.
NonlinearSolveFirstOrder.RobustMultiNewton
— FunctionRobustMultiNewton(
::Type{T} = Float64;
concrete_jac = nothing,
linsolve = nothing,
autodiff = nothing, vjp_autodiff = nothing, jvp_autodiff = nothing
)
A polyalgorithm focused on robustness. It uses a mixture of Newton methods with different globalizing techniques (trust region updates, line searches, etc.) in order to find a method that is able to adequately solve the minimization problem.
Basically, if this algorithm fails, then "most" good ways of solving your problem fail and you may need to think about reformulating the model (either there is an issue with the model, or more precision / more stable linear solver choice is required).
Arguments
T
: The eltype of the initial guess. It is only used to check if some of the algorithms are compatible with the problem type. Defaults toFloat64
.
Advanced Solvers
All of the previously mentioned solvers are wrappers around the following solvers. These are meant for advanced users and allow building custom solvers.
NonlinearSolveQuasiNewton.QuasiNewtonAlgorithm
— TypeQuasiNewtonAlgorithm(;
linesearch = missing, trustregion = missing, descent, update_rule, reinit_rule,
initialization, max_resets::Int = typemax(Int), name::Symbol = :unknown,
max_shrink_times::Int = typemax(Int), concrete_jac = Val(false)
)
Nonlinear Solve Algorithms using an Iterative Approximation of the Jacobian. Most common examples include Broyden
's Method.
Keyword Arguments
trustregion
: Globalization using a Trust Region Method. This needs to follow theNonlinearSolveBase.AbstractTrustRegionMethod
interface.descent
: The descent method to use to compute the step. This needs to follow theNonlinearSolveBase.AbstractDescentDirection
interface.max_shrink_times
: The maximum number of times the trust region radius can be shrunk before the algorithm terminates.update_rule
: The update rule to use to update the Jacobian. This needs to follow theNonlinearSolveBase.AbstractApproximateJacobianUpdateRule
interface.reinit_rule
: The reinitialization rule to use to reinitialize the Jacobian. This needs to follow theNonlinearSolveBase.AbstractResetCondition
interface.initialization
: The initialization method to use to initialize the Jacobian. This needs to follow theNonlinearSolveBase.AbstractJacobianInitialization
interface.
NonlinearSolveFirstOrder.GeneralizedFirstOrderAlgorithm
— TypeGeneralizedFirstOrderAlgorithm(;
descent, linesearch = missing,
trustregion = missing, autodiff = nothing, vjp_autodiff = nothing,
jvp_autodiff = nothing, max_shrink_times::Int = typemax(Int),
concrete_jac = Val(false), name::Symbol = :unknown
)
This is a Generalization of First-Order (uses Jacobian) Nonlinear Solve Algorithms. The most common example of this is Newton-Raphson Method.
First Order here refers to the order of differentiation, and should not be confused with the order of convergence.
Keyword Arguments
trustregion
: Globalization using a Trust Region Method. This needs to follow theNonlinearSolveBase.AbstractTrustRegionMethod
interface.descent
: The descent method to use to compute the step. This needs to follow theNonlinearSolveBase.AbstractDescentDirection
interface.max_shrink_times
: The maximum number of times the trust region radius can be shrunk before the algorithm terminates.
NonlinearSolveSpectralMethods.GeneralizedDFSane
— TypeGeneralizedDFSane(; linesearch, sigma_min, sigma_max, sigma_1, name::Symbol = :unknown)
A generalized version of the DF-SANE algorithm. This algorithm is a Jacobian-Free Spectral Method.
Arguments
linesearch
: Globalization using a Line Search Method. This is not optional currently, but that restriction might be lifted in the future.sigma_min
: The minimum spectral parameter allowed. This is used to ensure that the spectral parameter is not too small.sigma_max
: The maximum spectral parameter allowed. This is used to ensure that the spectral parameter is not too large.sigma_1
: The initial spectral parameter. If this is not provided, then the algorithm initializes it assigma_1 = <u, u> / <u, f(u)>
.