Internal Algorithm Helpers
Pseudo Transient Method
NonlinearSolve.SwitchedEvolutionRelaxation
— TypeSwitchedEvolutionRelaxation()
Method for updating the damping parameter in the PseudoTransient
method based on "switched evolution relaxation" [8] SER method.
NonlinearSolve.SwitchedEvolutionRelaxationCache
— TypeSwitchedEvolutionRelaxationCache <: AbstractDampingFunctionCache
Cache for the SwitchedEvolutionRelaxation
method.
Approximate Jacobian Methods
Initialization
NonlinearSolve.IdentityInitialization
— TypeIdentityInitialization(alpha, structure)
Initialize the Jacobian to be an Identity Matrix scaled by alpha
and maintain the structure as specified by structure
.
NonlinearSolve.TrueJacobianInitialization
— TypeTrueJacobianInitialization(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.BroydenLowRankInitialization
— TypeBroydenLowRankInitialization{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
NonlinearSolve.FullStructure
— TypeFullStructure()
Stores the full matrix.
NonlinearSolve.DiagonalStructure
— TypeDiagonalStructure()
Preserves only the Diagonal of the Matrix.
Jacobian Caches
NonlinearSolve.InitializedApproximateJacobianCache
— TypeInitializedApproximateJacobianCache(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 cacheNonlinearSolve.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.NoChangeInStateReset
— TypeNoChangeInStateReset(; 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 to3
.reset_tolerance
: the tolerance for the reset check. Defaults tosqrt(eps(real(eltype(u))))
.check_du
: whether to check the state. Defaults totrue
.check_dfu
: whether to check the function value. Defaults totrue
.
NonlinearSolve.IllConditionedJacobianReset
— TypeIllConditionedJacobianReset()
Recommend resetting the Jacobian if the current jacobian is ill-conditioned. This is used in Klement
.
Update Rules
NonlinearSolve.GoodBroydenUpdateRule
— TypeGoodBroydenUpdateRule()
Broyden Update Rule corresponding to "good broyden's method" [3].
NonlinearSolve.BadBroydenUpdateRule
— TypeBadBroydenUpdateRule()
Broyden Update Rule corresponding to "bad broyden's method" [3].
NonlinearSolve.KlementUpdateRule
— TypeKlementUpdateRule()
Update rule for Klement
.
Levenberg Marquardt Method
NonlinearSolve.LevenbergMarquardtTrustRegion
— TypeLevenbergMarquardtTrustRegion(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 forb_uphill
are1.0
or2.0
, withb_uphill = 2.0
allowing higher uphill moves thanb_uphill = 1.0
. Whenb_uphill = 0.0
, no uphill moves will be accepted. Defaults to1.0
. See Section 4 of Transtrum and Sethna [1].
Trust Region Method
NonlinearSolve.GenericTrustRegionScheme
— TypeGenericTrustRegionScheme(; 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 toRadiusUpdateSchemes.Simple
which follows the conventional approach. Other available schemes are documented inRadiusUpdateSchemes
,. 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 tomax(norm(fu), maximum(u) - minimum(u))
, except forRadiusUpdateSchemes.NLsolve
where it defaults toInf
.initial_trust_radius
: the initial trust region radius. Defaults tomax_trust_radius / 11
, except forRadiusUpdateSchemes.NLsolve
where it defaults tou0_norm > 0 ? u0_norm : 1
.step_threshold
: the threshold for taking a step. In every iteration, the threshold is compared with a valuer
, which is the actual reduction in the objective function divided by the predicted reduction. Ifstep_threshold > r
the model is not a good approximation, and the step is rejected. Defaults tonothing
.shrink_threshold
: the threshold for shrinking the trust region radius. In every iteration, the threshold is compared with a valuer
which is the actual reduction in the objective function divided by the predicted reduction. Ifshrink_threshold > r
the trust region radius is shrunk byshrink_factor
. Defaults tonothing
.expand_threshold
: the threshold for expanding the trust region radius. If a step is taken, i.estep_threshold < r
(withr
defined inshrink_threshold
), a check is also made to see ifexpand_threshold < r
. If that is true, the trust region radius is expanded byexpand_factor
. Defaults tonothing
.shrink_factor
: the factor to shrink the trust region radius with ifshrink_threshold > r
(withr
defined inshrink_threshold
). Defaults to0.25
.expand_factor
: the factor to expand the trust region radius with ifexpand_threshold < r
(withr
defined inshrink_threshold
). Defaults to2.0
.
Miscellaneous
SimpleNonlinearSolve.__nextfloat_tdir
— Function__nextfloat_tdir(x, x0, x1)
Move x
one floating point towards x1.
SimpleNonlinearSolve.__prevfloat_tdir
— Function__prevfloat_tdir(x, x0, x1)
Move x
one floating point towards x0.
SimpleNonlinearSolve.__max_tdir
— Function__max_tdir(a, b, x0, x1)
Return the maximum of a
and b
if x1 > x0
, otherwise return the minimum.
NonlinearSolve.callback_into_cache!
— Functioncallback_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_jac
— Functionconcrete_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.