Skip to content
This repository was archived by the owner on May 15, 2025. It is now read-only.
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
10 changes: 5 additions & 5 deletions Project.toml
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
name = "SimpleNonlinearSolve"
uuid = "727e6d20-b764-4bd8-a329-72de5adea6c7"
authors = ["SciML"]
version = "1.6.0"
version = "1.7.0"

[deps]
ADTypes = "47edcb42-4c32-4615-8424-f2b9edc5f35b"
Expand Down Expand Up @@ -39,13 +39,13 @@ SimpleNonlinearSolveZygoteExt = "Zygote"
ADTypes = "0.2.6"
AllocCheck = "0.1.1"
Aqua = "0.8"
ArrayInterface = "7.7"
ArrayInterface = "7.8"
CUDA = "5.2"
ChainRulesCore = "1.22"
ConcreteStructs = "0.2.3"
DiffEqBase = "6.146"
DiffEqBase = "6.149"
DiffResults = "1.1"
FastClosures = "0.3"
FastClosures = "0.3.2"
FiniteDiff = "2.22"
ForwardDiff = "0.10.36"
LinearAlgebra = "1.10"
Expand All @@ -59,7 +59,7 @@ Random = "1.10"
ReTestItems = "1.23"
Reexport = "1.2"
ReverseDiff = "1.15"
SciMLBase = "2.26.3"
SciMLBase = "2.28.0"
SciMLSensitivity = "7.56"
StaticArrays = "1.9"
StaticArraysCore = "1.4.2"
Expand Down
2 changes: 1 addition & 1 deletion src/nlsolve/broyden.jl
Original file line number Diff line number Diff line change
Expand Up @@ -48,7 +48,7 @@ function SciMLBase.__solve(prob::NonlinearProblem, alg::SimpleBroyden, args...;
@bb δJ⁻¹n = copy(x)
@bb δJ⁻¹ = copy(J⁻¹)

abstol, reltol, tc_cache = init_termination_cache(abstol, reltol, fx, x,
abstol, reltol, tc_cache = init_termination_cache(prob, abstol, reltol, fx, x,
termination_condition)

ls_cache = __get_linesearch(alg) === Val(true) ?
Expand Down
2 changes: 1 addition & 1 deletion src/nlsolve/dfsane.jl
Original file line number Diff line number Diff line change
Expand Up @@ -70,7 +70,7 @@ function SciMLBase.__solve(prob::NonlinearProblem, alg::SimpleDFSane{M}, args...
τ_min = T(alg.τ_min)
τ_max = T(alg.τ_max)

abstol, reltol, tc_cache = init_termination_cache(abstol, reltol, fx, x,
abstol, reltol, tc_cache = init_termination_cache(prob, abstol, reltol, fx, x,
termination_condition)

fx_norm = NONLINEARSOLVE_DEFAULT_NORM(fx)^nexp
Expand Down
2 changes: 1 addition & 1 deletion src/nlsolve/halley.jl
Original file line number Diff line number Diff line change
Expand Up @@ -34,7 +34,7 @@ function SciMLBase.__solve(prob::NonlinearProblem, alg::SimpleHalley, args...;
T = eltype(x)

autodiff = __get_concrete_autodiff(prob, alg.autodiff)
abstol, reltol, tc_cache = init_termination_cache(abstol, reltol, fx, x,
abstol, reltol, tc_cache = init_termination_cache(prob, abstol, reltol, fx, x,
termination_condition)

@bb xo = copy(x)
Expand Down
2 changes: 1 addition & 1 deletion src/nlsolve/klement.jl
Original file line number Diff line number Diff line change
Expand Up @@ -13,7 +13,7 @@ function SciMLBase.__solve(prob::NonlinearProblem, alg::SimpleKlement, args...;
T = eltype(x)
fx = _get_fx(prob, x)

abstol, reltol, tc_cache = init_termination_cache(abstol, reltol, fx, x,
abstol, reltol, tc_cache = init_termination_cache(prob, abstol, reltol, fx, x,
termination_condition)

@bb δx = copy(x)
Expand Down
2 changes: 1 addition & 1 deletion src/nlsolve/lbroyden.jl
Original file line number Diff line number Diff line change
Expand Up @@ -61,7 +61,7 @@ end

U, Vᵀ = __init_low_rank_jacobian(x, fx, x isa StaticArray ? threshold : Val(η))

abstol, reltol, tc_cache = init_termination_cache(abstol, reltol, fx, x,
abstol, reltol, tc_cache = init_termination_cache(prob, abstol, reltol, fx, x,
termination_condition)

@bb xo = copy(x)
Expand Down
2 changes: 1 addition & 1 deletion src/nlsolve/raphson.jl
Original file line number Diff line number Diff line change
Expand Up @@ -32,7 +32,7 @@ function SciMLBase.__solve(prob::Union{NonlinearProblem, NonlinearLeastSquaresPr
@bb xo = copy(x)
J, jac_cache = jacobian_cache(autodiff, prob.f, fx, x, prob.p)

abstol, reltol, tc_cache = init_termination_cache(abstol, reltol, fx, x,
abstol, reltol, tc_cache = init_termination_cache(prob, abstol, reltol, fx, x,
termination_condition)

for i in 1:maxiters
Expand Down
2 changes: 1 addition & 1 deletion src/nlsolve/trustRegion.jl
Original file line number Diff line number Diff line change
Expand Up @@ -88,7 +88,7 @@ function SciMLBase.__solve(prob::NonlinearProblem, alg::SimpleTrustRegion, args.
J, jac_cache = jacobian_cache(autodiff, prob.f, fx, x, prob.p)
fx, ∇f = value_and_jacobian(autodiff, prob.f, fx, x, prob.p, jac_cache; J)

abstol, reltol, tc_cache = init_termination_cache(abstol, reltol, fx, x,
abstol, reltol, tc_cache = init_termination_cache(prob, abstol, reltol, fx, x,
termination_condition)

# Set default trust region radius if not specified by user.
Expand Down
54 changes: 25 additions & 29 deletions src/utils.jl
Original file line number Diff line number Diff line change
Expand Up @@ -288,14 +288,30 @@ end
# different. NonlinearSolve is more for robust / cached solvers while SimpleNonlinearSolve
# is meant for low overhead solvers, users can opt into the other termination modes but the
# default is to use the least overhead version.
function init_termination_cache(abstol, reltol, du, u, ::Nothing)
return init_termination_cache(abstol, reltol, du, u, AbsNormTerminationMode())
function init_termination_cache(prob::NonlinearProblem, abstol, reltol, du, u, ::Nothing)
return init_termination_cache(prob, abstol, reltol, du, u,
AbsNormTerminationMode(Base.Fix1(maximum, abs)))
end
function init_termination_cache(abstol, reltol, du, u, tc::AbstractNonlinearTerminationMode)
function init_termination_cache(
prob::NonlinearLeastSquaresProblem, abstol, reltol, du, u, ::Nothing)
return init_termination_cache(prob, abstol, reltol, du, u,
AbsNormTerminationMode(Base.Fix2(norm, 2)))
end

function init_termination_cache(
prob::Union{NonlinearProblem, NonlinearLeastSquaresProblem},
abstol, reltol, du, u, tc::AbstractNonlinearTerminationMode)
T = promote_type(eltype(du), eltype(u))
abstol = __get_tolerance(u, abstol, T)
reltol = __get_tolerance(u, reltol, T)
tc_cache = init(du, u, tc; abstol, reltol)
tc_ = if hasfield(typeof(tc), :internalnorm) && tc.internalnorm === nothing
internalnorm = ifelse(
prob isa NonlinearProblem, Base.Fix1(maximum, abs), Base.Fix2(norm, 2))
DiffEqBase.set_termination_mode_internalnorm(tc, internalnorm)
else
tc
end
tc_cache = init(du, u, tc_; abstol, reltol, use_deprecated_retcodes = Val(false))
return DiffEqBase.get_abstol(tc_cache), DiffEqBase.get_reltol(tc_cache), tc_cache
end

Expand All @@ -305,45 +321,25 @@ function check_termination(tc_cache, fx, x, xo, prob, alg)
end
function check_termination(tc_cache, fx, x, xo, prob, alg,
::AbstractNonlinearTerminationMode)
if Bool(tc_cache(fx, x, xo))
tc_cache(fx, x, xo) &&
return build_solution(prob, alg, x, fx; retcode = ReturnCode.Success)
end
return nothing
end
function check_termination(tc_cache, fx, x, xo, prob, alg,
::AbstractSafeNonlinearTerminationMode)
if Bool(tc_cache(fx, x, xo))
if tc_cache.retcode == NonlinearSafeTerminationReturnCode.Success
retcode = ReturnCode.Success
elseif tc_cache.retcode == NonlinearSafeTerminationReturnCode.PatienceTermination
retcode = ReturnCode.ConvergenceFailure
elseif tc_cache.retcode == NonlinearSafeTerminationReturnCode.ProtectiveTermination
retcode = ReturnCode.Unstable
else
error("Unknown termination code: $(tc_cache.retcode)")
end
return build_solution(prob, alg, x, fx; retcode)
end
tc_cache(fx, x, xo) &&
return build_solution(prob, alg, x, fx; retcode = tc_cache.retcode)
return nothing
end
function check_termination(tc_cache, fx, x, xo, prob, alg,
::AbstractSafeBestNonlinearTerminationMode)
if Bool(tc_cache(fx, x, xo))
if tc_cache.retcode == NonlinearSafeTerminationReturnCode.Success
retcode = ReturnCode.Success
elseif tc_cache.retcode == NonlinearSafeTerminationReturnCode.PatienceTermination
retcode = ReturnCode.ConvergenceFailure
elseif tc_cache.retcode == NonlinearSafeTerminationReturnCode.ProtectiveTermination
retcode = ReturnCode.Unstable
else
error("Unknown termination code: $(tc_cache.retcode)")
end
if tc_cache(fx, x, xo)
if isinplace(prob)
prob.f(fx, x, prob.p)
else
fx = prob.f(x, prob.p)
end
return build_solution(prob, alg, tc_cache.u, fx; retcode)
return build_solution(prob, alg, tc_cache.u, fx; retcode = tc_cache.retcode)
end
return nothing
end
Expand Down