Internal Abstract Types
Solvers
NonlinearSolve.AbstractNonlinearSolveAlgorithm
— TypeAbstractNonlinearSolveAlgorithm{name} <: AbstractNonlinearAlgorithm
Abstract Type for all NonlinearSolve.jl Algorithms. name
can be used to define custom dispatches by wrapped solvers.
Interface Functions
concrete_jac(alg)
: whether or not the algorithm uses a concrete Jacobian. Defaults tonothing
.get_name(alg)
: get the name of the algorithm.
NonlinearSolve.AbstractNonlinearSolveExtensionAlgorithm
— TypeAbstractNonlinearSolveExtensionAlgorithm <: AbstractNonlinearSolveAlgorithm{:Extension}
Abstract Type for all NonlinearSolve.jl Extension Algorithms, i.e. wrappers over 3rd party solvers.
NonlinearSolve.AbstractNonlinearSolveCache
— TypeAbstractNonlinearSolveCache{iip, timeit}
Abstract Type for all NonlinearSolve.jl Caches.
Interface Functions
get_fu(cache)
: get the residual.get_u(cache)
: get the current state.set_fu!(cache, fu)
: set the residual.set_u!(cache, u)
: set the current state.reinit!(cache, u0; kwargs...)
: reinitialize the cache with the initial stateu0
and any additional keyword arguments.step!(cache; kwargs...)
: SeeSciMLBase.step!
for more details.not_terminated(cache)
: whether or not the solver has terminated.isinplace(cache)
: whether or not the solver is inplace.
Descent Algorithms
NonlinearSolve.AbstractDescentAlgorithm
— TypeAbstractDescentAlgorithm
Given the Jacobian J
and the residual fu
, this type of algorithm computes the descent direction δu
.
For non-square Jacobian problems, if we need to solve a linear solve problem, we use a least squares solver by default, unless the provided linsolve
can't handle non-square matrices, in which case we use the normal form equations $JᵀJ δu = Jᵀ fu$. Note that this factorization is often the faster choice, but it is not as numerically stable as the least squares solver.
__internal_init
specification
__internal_init(prob::NonlinearProblem{uType, iip}, alg::AbstractDescentAlgorithm, J, fu, u;
pre_inverted::Val{INV} = Val(false), linsolve_kwargs = (;), abstol = nothing,
reltol = nothing, alias_J::Bool = true, shared::Val{N} = Val(1),
kwargs...) where {INV, N, uType, iip} --> AbstractDescentCache
__internal_init(prob::NonlinearLeastSquaresProblem{uType, iip},
alg::AbstractDescentAlgorithm, J, fu, u; pre_inverted::Val{INV} = Val(false),
linsolve_kwargs = (;), abstol = nothing, reltol = nothing, alias_J::Bool = true,
shared::Val{N} = Val(1), kwargs...) where {INV, N, uType, iip} --> AbstractDescentCache
pre_inverted
: whether or not the Jacobian has been pre_inverted. Defaults toFalse
. Note that for most algorithms exceptNewtonDescent
setting it toVal(true)
is generally a bad idea.linsolve_kwargs
: keyword arguments to pass to the linear solver. Defaults to(;)
.abstol
: absolute tolerance for the linear solver. Defaults tonothing
.reltol
: relative tolerance for the linear solver. Defaults tonothing
.alias_J
: whether or not to alias the Jacobian. Defaults totrue
.shared
: Store multiple descent directions in the cache. Allows efficient and correct reuse of factorizations if needed,
Some of the algorithms also allow additional keyword arguments. See the documentation for the specific algorithm for more information.
Interface Functions
supports_trust_region(alg)
: whether or not the algorithm supports trust region methods. Defaults tofalse
.supports_line_search(alg)
: whether or not the algorithm supports line search methods. Defaults tofalse
.
See also NewtonDescent
, Dogleg
, SteepestDescent
, DampedNewtonDescent
.
NonlinearSolve.AbstractDescentCache
— TypeAbstractDescentCache
Abstract Type for all Descent Caches.
__internal_solve!
specification
δu, success, intermediates = __internal_solve!(cache::AbstractDescentCache, J, fu, u,
idx::Val; skip_solve::Bool = false, kwargs...)
J
: Jacobian or Inverse Jacobian (ifpre_inverted = Val(true)
).fu
: residual.u
: current state.idx
: index of the descent problem to solve and return. Defaults toVal(1)
.skip_solve
: Skip the direction computation and return the previous direction. Defaults tofalse
. This is useful for Trust Region Methods where the previous direction was rejected and we want to try with a modified trust region.kwargs
: keyword arguments to pass to the linear solver if there is one.
Returned values
δu
: the descent direction.success
: Certain Descent Algorithms can reject a descent direction for exampleGeodesicAcceleration
.intermediates
: A named tuple containing intermediates computed during the solve. For example,GeodesicAcceleration
returnsNamedTuple{(:v, :a)}
containing the "velocity" and "acceleration" terms.
Interface Functions
get_du(cache)
: get the descent direction.get_du(cache, ::Val{N})
: get theN
th descent direction.set_du!(cache, δu)
: set the descent direction.set_du!(cache, δu, ::Val{N})
: set theN
th descent direction.last_step_accepted(cache)
: whether or not the last step was accepted. Checks if the cache has alast_step_accepted
field and returns it if it does, else returnstrue
.
Approximate Jacobian
NonlinearSolve.AbstractApproximateJacobianStructure
— TypeAbstractApproximateJacobianStructure
Abstract Type for all Approximate Jacobian Structures used in NonlinearSolve.jl.
Interface Functions
stores_full_jacobian(alg)
: whether or not the algorithm stores the full Jacobian. Defaults tofalse
.get_full_jacobian(cache, alg, J)
: get the full Jacobian. Defaults to throwing an error ifstores_full_jacobian(alg)
isfalse
.
NonlinearSolve.AbstractJacobianInitialization
— TypeAbstractJacobianInitialization
Abstract Type for all Jacobian Initialization Algorithms used in NonlinearSolve.jl.
Interface Functions
jacobian_initialized_preinverted(alg)
: whether or not the Jacobian is initialized preinverted. Defaults tofalse
.
__internal_init
specification
__internal_init(prob::AbstractNonlinearProblem, alg::AbstractJacobianInitialization,
solver, f::F, fu, u, p; linsolve = missing, internalnorm::IN = DEFAULT_NORM,
kwargs...)
Returns a NonlinearSolve.InitializedApproximateJacobianCache
.
All subtypes need to define (cache::InitializedApproximateJacobianCache)(alg::NewSubType, fu, u)
which reinitializes the Jacobian in cache.J
.
NonlinearSolve.AbstractApproximateJacobianUpdateRule
— TypeAbstractApproximateJacobianUpdateRule{INV}
Abstract Type for all Approximate Jacobian Update Rules used in NonlinearSolve.jl.
Interface Functions
store_inverse_jacobian(alg)
: ReturnINV
__internal_init
specification
__internal_init(prob::AbstractNonlinearProblem,
alg::AbstractApproximateJacobianUpdateRule, J, fu, u, du, args...;
internalnorm::F = DEFAULT_NORM,
kwargs...) where {F} --> AbstractApproximateJacobianUpdateRuleCache{INV}
NonlinearSolve.AbstractApproximateJacobianUpdateRuleCache
— TypeAbstractApproximateJacobianUpdateRuleCache{INV}
Abstract Type for all Approximate Jacobian Update Rule Caches used in NonlinearSolve.jl.
Interface Functions
store_inverse_jacobian(alg)
: ReturnINV
__internal_solve!
specification
__internal_solve!(cache::AbstractApproximateJacobianUpdateRuleCache, J, fu, u, du;
kwargs...) --> J / J⁻¹
NonlinearSolve.AbstractResetCondition
— TypeAbstractResetCondition
Condition for resetting the Jacobian in Quasi-Newton's methods.
__internal_init
specification
__internal_init(alg::AbstractResetCondition, J, fu, u, du, args...;
kwargs...) --> ResetCache
__internal_solve!
specification
__internal_solve!(cache::ResetCache, J, fu, u, du) --> Bool
Damping Algorithms
NonlinearSolve.AbstractDampingFunction
— TypeAbstractDampingFunction
Abstract Type for Damping Functions in DampedNewton.
__internal_init
specification
__internal_init(prob::AbstractNonlinearProblem, f::AbstractDampingFunction, initial_damping,
J, fu, u, args...; internal_norm = DEFAULT_NORM,
kwargs...) --> AbstractDampingFunctionCache
Returns a AbstractDampingFunctionCache
.
NonlinearSolve.AbstractDampingFunctionCache
— TypeAbstractDampingFunctionCache
Abstract Type for the Caches created by AbstractDampingFunctions
Interface Functions
requires_normal_form_jacobian(f)
: whether or not the Jacobian is needed in normal form. No default.requires_normal_form_rhs(f)
: whether or not the residual is needed in normal form. No default.returns_norm_form_damping(f)
: whether or not the damping function returns the damping factor in normal form. Defaults torequires_normal_form_jacobian(f) || requires_normal_form_rhs(f)
.(cache::AbstractDampingFunctionCache)(::Nothing)
: returns the damping factor. The type of the damping factor returned fromsolve!
is guaranteed to be the same as this.
__internal_solve!
specification
__internal_solve!(cache::AbstractDampingFunctionCache, J, fu, args...; kwargs...)
Returns the damping factor.
Line Search
NonlinearSolve.AbstractNonlinearSolveLineSearchAlgorithm
— TypeAbstractNonlinearSolveLineSearchAlgorithm
Abstract Type for all Line Search Algorithms used in NonlinearSolve.jl.
__internal_init
specification
__internal_init(prob::AbstractNonlinearProblem,
alg::AbstractNonlinearSolveLineSearchAlgorithm, f::F, fu, u, p, args...;
internalnorm::IN = DEFAULT_NORM,
kwargs...) where {F, IN} --> AbstractNonlinearSolveLineSearchCache
NonlinearSolve.AbstractNonlinearSolveLineSearchCache
— TypeAbstractNonlinearSolveLineSearchCache
Abstract Type for all Line Search Caches used in NonlinearSolve.jl.
__internal_solve!
specification
__internal_solve!(cache::AbstractNonlinearSolveLineSearchCache, u, du; kwargs...)
Returns 2 values:
unsuccessful
: Iftrue
it means that the Line Search Failed.alpha
: The step size.
Trust Region
NonlinearSolve.AbstractTrustRegionMethod
— TypeAbstractTrustRegionMethod
Abstract Type for all Trust Region Methods used in NonlinearSolve.jl.
__internal_init
specification
__internal_init(prob::AbstractNonlinearProblem, alg::AbstractTrustRegionMethod,
f::F, fu, u, p, args...; internalnorm::IF = DEFAULT_NORM,
kwargs...) where {F, IF} --> AbstractTrustRegionMethodCache
NonlinearSolve.AbstractTrustRegionMethodCache
— TypeAbstractTrustRegionMethodCache
Abstract Type for all Trust Region Method Caches used in NonlinearSolve.jl.
Interface Functions
last_step_accepted(cache)
: whether or not the last step was accepted. Defaults tocache.last_step_accepted
. Should if overloaded if the field is not present.
__internal_solve!
specification
__internal_solve!(cache::AbstractTrustRegionMethodCache, J, fu, u, δu, descent_stats)
Returns last_step_accepted
, updated u_cache
and fu_cache
. If the last step was accepted then these values should be copied into the toplevel cache.
Tracing
NonlinearSolve.AbstractNonlinearSolveTraceLevel
— TypeAbstractNonlinearSolveTraceLevel
Common Arguments
freq
: Sets bothprint_frequency
andstore_frequency
tofreq
.
Common Keyword Arguments
print_frequency
: Print the trace everyprint_frequency
iterations ifshow_trace == Val(true)
.store_frequency
: Store the trace everystore_frequency
iterations ifstore_trace == Val(true)
.