Internal Abstract Types
Solvers
NonlinearSolveBase.AbstractNonlinearSolveAlgorithm
— TypeAbstractNonlinearSolveAlgorithm <: AbstractNonlinearAlgorithm
Abstract Type for all NonlinearSolveBase Algorithms.
Interface Functions
concrete_jac(alg)
: whether or not the algorithm uses a concrete Jacobian. Defaults tonothing
.
NonlinearSolveBase.AbstractNonlinearSolveCache
— TypeAbstractNonlinearSolveCache
Abstract Type for all NonlinearSolveBase Caches.
Interface Functions
get_fu(cache)
: get the residual.get_u(cache)
: get the current state.set_fu!(cache, fu)
: set the residual.has_time_limit(cache)
: whether or not the solver has a maximum time limit.not_terminated(cache)
: whether or not the solver has terminated.SciMLBase.set_u!(cache, u)
: set the current state.SciMLBase.reinit!(cache, u0; kwargs...)
: reinitialize the cache with the initial stateu0
and any additional keyword arguments.SciMLBase.isinplace(cache)
: whether or not the solver is inplace.CommonSolve.step!(cache; kwargs...)
: SeeCommonSolve.step!
for more details.get_abstol(cache)
: get theabstol
provided to the cache.get_reltol(cache)
: get thereltol
provided to the cache.
Additionally implements SymbolicIndexingInterface
interface Functions.
Expected Fields in Sub-Types
For the default interface implementations we expect the following fields to be present in the cache:
fu
: the residual.u
: the current state.maxiters
: the maximum number of iterations.nsteps
: the number of steps taken.force_stop
: whether or not the solver has been forced to stop.retcode
: the return code.stats
:NLStats
object.alg
: the algorithm.maxtime
: the maximum time limit for the solver. (Optional)timer
: the timer for the solver. (Optional)total_time
: the total time taken by the solver. (Optional)
Descent Directions
NonlinearSolveBase.AbstractDescentDirection
— TypeAbstractDescentDirection
Abstract Type for all Descent Directions used in NonlinearSolveBase. Given the Jacobian J
and the residual fu
, these algorithms compute 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.
InternalAPI.init
specification
InternalAPI.init(
prob::AbstractNonlinearProblem, alg::AbstractDescentDirection, J, fu, u;
pre_inverted::Val = Val(false), linsolve_kwargs = (;),
abstol = nothing, reltol = nothing, alias_J::Bool = true,
shared::Val = Val(1), kwargs...
)::AbstractDescentCache
pre_inverted
: whether or not the Jacobian has been pre_inverted.linsolve_kwargs
: keyword arguments to pass to the linear solver.abstol
: absolute tolerance for the linear solver.reltol
: relative tolerance for the linear solver.alias_J
: whether or not to alias the Jacobian.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
.
NonlinearSolveBase.AbstractDescentCache
— TypeAbstractDescentCache
Abstract Type for all Descent Caches.
InternalAPI.solve!
specification
InternalAPI.solve!(
cache::AbstractDescentCache, J, fu, u, idx::Val;
skip_solve::Bool = false, new_jacobian::Bool = true, kwargs...
)::DescentResult
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.new_jacobian
: Whether the Jacobian has been updated. Defaults totrue
.kwargs
: keyword arguments to pass to the linear solver if there is one.
Returned values
descent_result
: Result in aDescentResult
.
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
.preinverted_jacobian(cache)
: whether or not the Jacobian has been preinverted.normal_form(cache)
: whether or not the linear solver uses normal form.
Descent Results
NonlinearSolveBase.DescentResult
— TypeDescentResult(;
δu = missing, u = missing, success::Bool = true, linsolve_success::Bool = true,
extras = (;)
)
Construct a DescentResult
object.
Keyword Arguments
δu
: The descent direction.u
: The new iterate. This is provided only for multi-step methods currently.success
: Certain Descent Algorithms can reject a descent direction for exampleGeodesicAcceleration
.linsolve_success
: Whether the line search was successful.extras
: A named tuple containing intermediates computed during the solve. For example,GeodesicAcceleration
returnsNamedTuple{(:v, :a)}
containing the "velocity" and "acceleration" terms.
Approximate Jacobian
NonlinearSolveBase.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
.
NonlinearSolveBase.AbstractJacobianInitialization
— TypeAbstractJacobianInitialization
Abstract Type for all Jacobian Initialization Algorithms used in NonlinearSolveBase.
Interface Functions
jacobian_initialized_preinverted(alg)
: whether or not the Jacobian is initialized preinverted. Defaults tofalse
.
InternalAPI.init
specification
InternalAPI.init(
prob::AbstractNonlinearProblem, alg::AbstractJacobianInitialization, solver,
f, fu, u, p;
linsolve = missing, internalnorm::IN = L2_NORM, kwargs...
)::AbstractJacobianCache
All subtypes need to define (cache::AbstractJacobianCache)(alg::NewSubType, fu, u)
which reinitializes the Jacobian in cache.J
.
NonlinearSolveBase.AbstractApproximateJacobianUpdateRule
— TypeAbstractApproximateJacobianUpdateRule
Abstract Type for all Approximate Jacobian Update Rules used in NonlinearSolveBase.
Interface Functions
store_inverse_jacobian(alg)
: Returnalg.store_inverse_jacobian
InternalAPI.init
specification
InternalAPI.init(
prob::AbstractNonlinearProblem, alg::AbstractApproximateJacobianUpdateRule, J, fu, u,
du, args...; internalnorm = L2_NORM, kwargs...
)::AbstractApproximateJacobianUpdateRuleCache
NonlinearSolveBase.AbstractApproximateJacobianUpdateRuleCache
— TypeAbstractApproximateJacobianUpdateRuleCache
Abstract Type for all Approximate Jacobian Update Rule Caches used in NonlinearSolveBase.
Interface Functions
store_inverse_jacobian(cache)
: Returnstore_inverse_jacobian(cache.rule)
InternalAPI.solve!
specification
InternalAPI.solve!(
cache::AbstractApproximateJacobianUpdateRuleCache, J, fu, u, du; kwargs...
) --> J / J⁻¹
NonlinearSolveBase.AbstractResetCondition
— TypeAbstractResetCondition
Condition for resetting the Jacobian in Quasi-Newton's methods.
InternalAPI.init
specification
InternalAPI.init(
alg::AbstractResetCondition, J, fu, u, du, args...; kwargs...
)::AbstractResetConditionCache
Damping Algorithms
NonlinearSolveBase.AbstractDampingFunction
— TypeAbstractDampingFunction
Abstract Type for Damping Functions in DampedNewton.
InternalAPI.init
specification
InternalAPI.init(
prob::AbstractNonlinearProblem, f::AbstractDampingFunction, initial_damping,
J, fu, u, args...;
internalnorm = L2_NORM, kwargs...
)::AbstractDampingFunctionCache
NonlinearSolveBase.AbstractDampingFunctionCache
— TypeAbstractDampingFunctionCache
Abstract Type for the Caches created by AbstractDampingFunctions
Interface Functions
requires_normal_form_jacobian(alg)
: whether or not the Jacobian is needed in normal form. No default.requires_normal_form_rhs(alg)
: whether or not the residual is needed in normal form. No default.returns_norm_form_damping(alg)
: whether or not the damping function returns the damping factor in normal form. Defaults torequires_normal_form_jacobian(alg) || requires_normal_form_rhs(alg)
.(cache::AbstractDampingFunctionCache)(::Nothing)
: returns the damping factor. The type of the damping factor returned fromsolve!
is guaranteed to be the same as this.
InternalAPI.solve!
specification
InternalAPI.solve!(
cache::AbstractDampingFunctionCache, J, fu, u, δu, descent_stats
)
Returns the damping factor.
Trust Region
NonlinearSolveBase.AbstractTrustRegionMethod
— TypeAbstractTrustRegionMethod
Abstract Type for all Trust Region Methods used in NonlinearSolveBase.
InternalAPI.init
specification
InternalAPI.init(
prob::AbstractNonlinearProblem, alg::AbstractTrustRegionMethod, f, fu, u, p, args...;
internalnorm = L2_NORM, kwargs...
)::AbstractTrustRegionMethodCache
NonlinearSolveBase.AbstractTrustRegionMethodCache
— TypeAbstractTrustRegionMethodCache
Abstract Type for all Trust Region Method Caches used in NonlinearSolveBase.
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.
InternalAPI.solve!
specification
InternalAPI.solve!(
cache::AbstractTrustRegionMethodCache, J, fu, u, δu, descent_stats; kwargs...
)
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.