Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Scalar nonlinear solve AD #11

Merged
merged 5 commits into from
Nov 22, 2020
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
1 change: 1 addition & 0 deletions Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -20,6 +20,7 @@ Reexport = "0.2"
Setfield = "0.7"
StaticArrays = "0.11, 0.12"
UnPack = "0.1, 1.0"
julia = "1"

[extras]
BenchmarkTools = "6e4b80f9-dd63-53aa-95a3-0cdb28fa8baf"
Expand Down
1 change: 1 addition & 0 deletions src/NonlinearSolve.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,7 @@ module NonlinearSolve
using Reexport
using UnPack: @unpack
using FiniteDiff, ForwardDiff
using ForwardDiff: Dual
using Setfield
using StaticArrays
using RecursiveArrayTools
Expand Down
45 changes: 45 additions & 0 deletions src/scalar.jl
Original file line number Diff line number Diff line change
Expand Up @@ -19,6 +19,51 @@ function solve(prob::NonlinearProblem{<:Number}, ::NewtonRaphson, args...; xatol
return NewtonSolution(x, MAXITERS_EXCEED)
end

function scalar_nlsolve_ad(prob, alg, args...; kwargs...)
f = prob.f
p = value(prob.p)
u0 = value(prob.u0)

newprob = NonlinearProblem(f, u0, p; prob.kwargs...)
sol = solve(newprob, alg, args...; kwargs...)

uu = getsolution(sol)
if p isa Number
f_p = ForwardDiff.derivative(Base.Fix1(f, uu), p)
else
f_p = ForwardDiff.gradient(Base.Fix1(f, uu), p)
end

f_x = ForwardDiff.derivative(Base.Fix2(f, p), uu)
pp = prob.p
sumfun = let f_x′ = -f_x
((fp, p),) -> (fp / f_x′) * ForwardDiff.partials(p)
end
partials = sum(sumfun, zip(f_p, pp))
return sol, partials
end

function solve(prob::NonlinearProblem{<:Number, iip, <:Dual{T,V,P}}, alg::NewtonRaphson, args...; kwargs...) where {iip, T, V, P}
sol, partials = scalar_nlsolve_ad(prob, alg, args...; kwargs...)
return NewtonSolution(Dual{T,V,P}(sol.u, partials), sol.retcode)
end
function solve(prob::NonlinearProblem{<:Number, iip, <:AbstractArray{<:Dual{T,V,P}}}, alg::NewtonRaphson, args...; kwargs...) where {iip, T, V, P}
sol, partials = scalar_nlsolve_ad(prob, alg, args...; kwargs...)
return NewtonSolution(Dual{T,V,P}(sol.u, partials), sol.retcode)
end

# avoid ambiguities
for Alg in [Bisection, Falsi]
@eval function solve(prob::NonlinearProblem{uType, iip, <:Dual{T,V,P}}, alg::$Alg, args...; kwargs...) where {uType, iip, T, V, P}
sol, partials = scalar_nlsolve_ad(prob, alg, args...; kwargs...)
return BracketingSolution(Dual{T,V,P}(sol.left, partials), Dual{T,V,P}(sol.right, partials), sol.retcode)
end
@eval function solve(prob::NonlinearProblem{uType, iip, <:AbstractArray{<:Dual{T,V,P}}}, alg::$Alg, args...; kwargs...) where {uType, iip, T, V, P}
sol, partials = scalar_nlsolve_ad(prob, alg, args...; kwargs...)
return BracketingSolution(Dual{T,V,P}(sol.left, partials), Dual{T,V,P}(sol.right, partials), sol.retcode)
end
end

function solve(prob::NonlinearProblem, ::Bisection, args...; maxiters = 1000, kwargs...)
f = Base.Fix2(prob.f, prob.p)
left, right = prob.u0
Expand Down
3 changes: 3 additions & 0 deletions src/types.jl
Original file line number Diff line number Diff line change
Expand Up @@ -78,3 +78,6 @@ function sync_residuals!(solver::BracketingImmutableSolver)
@set! solver.fr = solver.f(solver.right, solver.p)
solver
end

getsolution(sol::NewtonSolution) = sol.u
getsolution(sol::BracketingSolution) = sol.left
6 changes: 5 additions & 1 deletion src/utils.jl
Original file line number Diff line number Diff line change
Expand Up @@ -216,7 +216,7 @@ Move `x` one floating point towards x0.
function prevfloat_tdir(x::T, x0::T, x1::T)::T where {T}
x1 > x0 ? prevfloat(x) : nextfloat(x)
end

function nextfloat_tdir(x::T, x0::T, x1::T)::T where {T}
x1 > x0 ? nextfloat(x) : prevfloat(x)
end
Expand All @@ -234,3 +234,7 @@ function value_derivative(f::F, x::R) where {F,R}
out = f(ForwardDiff.Dual{T}(x, one(x)))
ForwardDiff.value(out), ForwardDiff.extract_derivative(T, out)
end

value(x) = x
value(x::Dual) = ForwardDiff.value(x)
value(x::AbstractArray{<:Dual}) = map(ForwardDiff.value, x)
27 changes: 25 additions & 2 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -57,18 +57,41 @@ end
f, u0 = (u, p) -> u * u - p, 1.0

g = function (p)
probN = NonlinearProblem{false}(f, u0, p)
probN = NonlinearProblem{false}(f, oftype(p, u0), p)
sol = solve(probN, NewtonRaphson())
return sol.u
end

@test_broken ForwardDiff.derivative(g, 1.0) ≈ 0.5
@test ForwardDiff.derivative(g, 1.0) ≈ 0.5

for p in 1.1:0.1:100.0
@test g(p) ≈ sqrt(p)
@test ForwardDiff.derivative(g, p) ≈ 1/(2*sqrt(p))
end

f, u0 = (u, p) -> p[1] * u * u - p[2], (1.0, 100.0)
t = (p) -> [sqrt(p[2] / p[1])]
p = [0.9, 50.0]
for alg in [Bisection(), Falsi()]
global g, p
g = function (p)
probN = NonlinearProblem{false}(f, u0, p)
sol = solve(probN, Bisection())
return [sol.left]
end

@test g(p) ≈ [sqrt(p[2] / p[1])]
@test ForwardDiff.jacobian(g, p) ≈ ForwardDiff.jacobian(t, p)
end

gnewton = function (p)
probN = NonlinearProblem{false}(f, 0.5, p)
sol = solve(probN, NewtonRaphson())
return [sol.u]
end
@test gnewton(p) ≈ [sqrt(p[2] / p[1])]
@test ForwardDiff.jacobian(gnewton, p) ≈ ForwardDiff.jacobian(t, p)

# Error Checks

f, u0 = (u, p) -> u .* u .- 2.0, @SVector[1.0, 1.0]
Expand Down