Internal Abstract Types

Solvers

NonlinearSolve.AbstractNonlinearSolveAlgorithmType
AbstractNonlinearSolveAlgorithm{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 to nothing.
  • get_name(alg): get the name of the algorithm.
NonlinearSolve.AbstractNonlinearSolveCacheType
AbstractNonlinearSolveCache{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 state u0 and any additional keyword arguments.
  • step!(cache; kwargs...): See SciMLBase.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.AbstractDescentAlgorithmType
AbstractDescentAlgorithm

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 to False. Note that for most algorithms except NewtonDescent setting it to Val(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 to nothing.
  • reltol: relative tolerance for the linear solver. Defaults to nothing.
  • alias_J: whether or not to alias the Jacobian. Defaults to true.
  • 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 to false.
  • supports_line_search(alg): whether or not the algorithm supports line search methods. Defaults to false.

See also NewtonDescent, Dogleg, SteepestDescent, DampedNewtonDescent.

NonlinearSolve.AbstractDescentCacheType
AbstractDescentCache

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 (if pre_inverted = Val(true)).
  • fu: residual.
  • u: current state.
  • idx: index of the descent problem to solve and return. Defaults to Val(1).
  • skip_solve: Skip the direction computation and return the previous direction. Defaults to false. 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 example GeodesicAcceleration.
  • intermediates: A named tuple containing intermediates computed during the solve. For example, GeodesicAcceleration returns NamedTuple{(:v, :a)} containing the "velocity" and "acceleration" terms.

Interface Functions

  • get_du(cache): get the descent direction.
  • get_du(cache, ::Val{N}): get the Nth descent direction.
  • set_du!(cache, δu): set the descent direction.
  • set_du!(cache, δu, ::Val{N}): set the Nth descent direction.
  • last_step_accepted(cache): whether or not the last step was accepted. Checks if the cache has a last_step_accepted field and returns it if it does, else returns true.

Approximate Jacobian

NonlinearSolve.AbstractApproximateJacobianStructureType
AbstractApproximateJacobianStructure

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 to false.
  • get_full_jacobian(cache, alg, J): get the full Jacobian. Defaults to throwing an error if stores_full_jacobian(alg) is false.
NonlinearSolve.AbstractJacobianInitializationType
AbstractJacobianInitialization

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 to false.

__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.AbstractApproximateJacobianUpdateRuleType
AbstractApproximateJacobianUpdateRule{INV}

Abstract Type for all Approximate Jacobian Update Rules used in NonlinearSolve.jl.

Interface Functions

  • store_inverse_jacobian(alg): Return INV

__internal_init specification

__internal_init(prob::AbstractNonlinearProblem,
    alg::AbstractApproximateJacobianUpdateRule, J, fu, u, du, args...;
    internalnorm::F = DEFAULT_NORM,
    kwargs...) where {F} --> AbstractApproximateJacobianUpdateRuleCache{INV}
NonlinearSolve.AbstractApproximateJacobianUpdateRuleCacheType
AbstractApproximateJacobianUpdateRuleCache{INV}

Abstract Type for all Approximate Jacobian Update Rule Caches used in NonlinearSolve.jl.

Interface Functions

  • store_inverse_jacobian(alg): Return INV

__internal_solve! specification

__internal_solve!(cache::AbstractApproximateJacobianUpdateRuleCache, J, fu, u, du;
    kwargs...) --> J / J⁻¹
NonlinearSolve.AbstractResetConditionType
AbstractResetCondition

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.AbstractDampingFunctionType
AbstractDampingFunction

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.AbstractDampingFunctionCacheType
AbstractDampingFunctionCache

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 to requires_normal_form_jacobian(f) || requires_normal_form_rhs(f).
  • (cache::AbstractDampingFunctionCache)(::Nothing): returns the damping factor. The type of the damping factor returned from solve! is guaranteed to be the same as this.

__internal_solve! specification

__internal_solve!(cache::AbstractDampingFunctionCache, J, fu, args...; kwargs...)

Returns the damping factor.

NonlinearSolve.AbstractNonlinearSolveLineSearchAlgorithmType
AbstractNonlinearSolveLineSearchAlgorithm

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.AbstractNonlinearSolveLineSearchCacheType
AbstractNonlinearSolveLineSearchCache

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: If true it means that the Line Search Failed.
  • alpha: The step size.

Trust Region

NonlinearSolve.AbstractTrustRegionMethodType
AbstractTrustRegionMethod

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.AbstractTrustRegionMethodCacheType
AbstractTrustRegionMethodCache

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 to cache.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.AbstractNonlinearSolveTraceLevelType
AbstractNonlinearSolveTraceLevel

Common Arguments

  • freq: Sets both print_frequency and store_frequency to freq.

Common Keyword Arguments

  • print_frequency: Print the trace every print_frequency iterations if show_trace == Val(true).
  • store_frequency: Store the trace every store_frequency iterations if store_trace == Val(true).