Internal Algorithm Helpers

Pseudo Transient Method

Approximate Jacobian Methods

Initialization

NonlinearSolve.IdentityInitializationType
IdentityInitialization(alpha, structure)

Initialize the Jacobian to be an Identity Matrix scaled by alpha and maintain the structure as specified by structure.

NonlinearSolve.TrueJacobianInitializationType
TrueJacobianInitialization(structure, autodiff)

Initialize the Jacobian to be the true Jacobian and maintain the structure as specified by structure. autodiff is used to compute the true Jacobian and if not specified we make a selection automatically.

NonlinearSolve.BroydenLowRankInitializationType
BroydenLowRankInitialization{T}(alpha, threshold::Val{T})

An initialization for LimitedMemoryBroyden that uses a low rank approximation of the Jacobian. The low rank updates to the Jacobian matrix corresponds to what SciPy calls "simple".

Jacobian Structure

Jacobian Caches

NonlinearSolve.InitializedApproximateJacobianCacheType
InitializedApproximateJacobianCache(J, structure, alg, cache, initialized::Bool,
    internalnorm)

A cache for Approximate Jacobian.

Arguments

  • J: The current Jacobian.
  • structure: The structure of the Jacobian.
  • alg: The initialization algorithm.
  • cache: The Jacobian cache NonlinearSolve.JacobianCache (if needed).
  • initialized: A boolean indicating whether the Jacobian has been initialized.
  • internalnorm: The norm to be used.

Interface

(cache::InitializedApproximateJacobianCache)(::Nothing)

Returns the current Jacobian cache.J with the proper structure.

__internal_solve!(cache::InitializedApproximateJacobianCache, fu, u, ::Val{reinit})

Solves for the Jacobian cache.J and returns it. If reinit is true, then the Jacobian is reinitialized.

Reset Methods

NonlinearSolve.NoChangeInStateResetType
NoChangeInStateReset(; nsteps::Int = 3, reset_tolerance = nothing,
    check_du::Bool = true, check_dfu::Bool = true)

Recommends a reset if the state or the function value has not changed significantly in nsteps steps. This is used in Broyden.

Keyword Arguments

  • nsteps: the number of steps to check for no change. Defaults to 3.
  • reset_tolerance: the tolerance for the reset check. Defaults to sqrt(eps(real(eltype(u)))).
  • check_du: whether to check the state. Defaults to true.
  • check_dfu: whether to check the function value. Defaults to true.

Update Rules

Levenberg Marquardt Method

NonlinearSolve.LevenbergMarquardtTrustRegionType
LevenbergMarquardtTrustRegion(b_uphill)

Trust Region method for LevenbergMarquardt. This method is tightly coupled with the Levenberg-Marquardt method and works by directly updating the damping parameter instead of specifying a trust region radius.

Arguments

  • b_uphill: a factor that determines if a step is accepted or rejected. The standard choice in the Levenberg-Marquardt method is to accept all steps that decrease the cost and reject all steps that increase the cost. Although this is a natural and safe choice, it is often not the most efficient. Therefore downhill moves are always accepted, but uphill moves are only conditionally accepted. To decide whether an uphill move will be accepted at each iteration $i$, we compute $\beta_i = \cos(v_{\text{new}}, v_{\text{old}})$, which denotes the cosine angle between the proposed velocity $v_{\text{new}}$ and the velocity of the last accepted step $v_{\text{old}}$. The idea is to accept uphill moves if the angle is small. To specify, uphill moves are accepted if $(1-\beta_i)^{b_{\text{uphill}}} C_{i+1} \le C_i$, where $C_i$ is the cost at iteration $i$. Reasonable choices for b_uphill are 1.0 or 2.0, with b_uphill = 2.0 allowing higher uphill moves than b_uphill = 1.0. When b_uphill = 0.0, no uphill moves will be accepted. Defaults to 1.0. See Section 4 of Transtrum and Sethna [1].

Trust Region Method

NonlinearSolve.GenericTrustRegionSchemeType
GenericTrustRegionScheme(; method = RadiusUpdateSchemes.Simple,
    max_trust_radius = nothing, initial_trust_radius = nothing,
    step_threshold = nothing, shrink_threshold = nothing, expand_threshold = nothing,
    shrink_factor = nothing, expand_factor = nothing, forward_ad = nothing,
    reverse_ad = nothing)

Trust Region Method that updates and stores the current trust region radius in trust_region. For any of the keyword arguments, if the value is nothing, then we use the value used in the respective paper.

Keyword Arguments

  • radius_update_scheme: the choice of radius update scheme to be used. Defaults to RadiusUpdateSchemes.Simple which follows the conventional approach. Other available schemes are documented in RadiusUpdateSchemes,. These schemes have the trust region radius converging to zero that is seen to improve convergence. For more details, see [1].
  • max_trust_radius: the maximal trust region radius. Defaults to max(norm(fu), maximum(u) - minimum(u)), except for RadiusUpdateSchemes.NLsolve where it defaults to Inf.
  • initial_trust_radius: the initial trust region radius. Defaults to max_trust_radius / 11, except for RadiusUpdateSchemes.NLsolve where it defaults to u0_norm > 0 ? u0_norm : 1.
  • step_threshold: the threshold for taking a step. In every iteration, the threshold is compared with a value r, which is the actual reduction in the objective function divided by the predicted reduction. If step_threshold > r the model is not a good approximation, and the step is rejected. Defaults to nothing.
  • shrink_threshold: the threshold for shrinking the trust region radius. In every iteration, the threshold is compared with a value r which is the actual reduction in the objective function divided by the predicted reduction. If shrink_threshold > r the trust region radius is shrunk by shrink_factor. Defaults to nothing.
  • expand_threshold: the threshold for expanding the trust region radius. If a step is taken, i.e step_threshold < r (with r defined in shrink_threshold), a check is also made to see if expand_threshold < r. If that is true, the trust region radius is expanded by expand_factor. Defaults to nothing.
  • shrink_factor: the factor to shrink the trust region radius with if shrink_threshold > r (with r defined in shrink_threshold). Defaults to 0.25.
  • expand_factor: the factor to expand the trust region radius with if expand_threshold < r (with r defined in shrink_threshold). Defaults to 2.0.

Miscellaneous

NonlinearSolve.callback_into_cache!Function
callback_into_cache!(cache, internalcache, args...)

Define custom operations on internalcache tightly coupled with the calling cache. args... contain the sequence of caches calling into internalcache.

This unfortunately makes code very tightly coupled and not modular. It is recommended to not use this functionality unless it can't be avoided (like in LevenbergMarquardt).

NonlinearSolve.concrete_jacFunction
concrete_jac(alg::AbstractNonlinearSolveAlgorithm)

Whether the algorithm uses a concrete Jacobian. Defaults to nothing if it is unknown or not applicable. Else a boolean value is returned.