Skip to content
This repository was archived by the owner on May 15, 2025. It is now read-only.

Commit bc5244d

Browse files
Merge pull request #49 from yash2798/ys/halley
Updated halley's method for multivariate functions
2 parents 77fcd76 + 246be18 commit bc5244d

File tree

2 files changed

+47
-26
lines changed

2 files changed

+47
-26
lines changed

src/halley.jl

Lines changed: 37 additions & 15 deletions
Original file line numberDiff line numberDiff line change
@@ -43,9 +43,8 @@ function SciMLBase.__solve(prob::NonlinearProblem,
4343
f = Base.Fix2(prob.f, prob.p)
4444
x = float(prob.u0)
4545
fx = f(x)
46-
# fx = float(prob.u0)
47-
if !isa(fx, Number) || !isa(x, Number)
48-
error("Halley currently only supports scalar-valued single-variable functions")
46+
if isa(x, AbstractArray)
47+
n = length(x)
4948
end
5049
T = typeof(x)
5150

@@ -65,22 +64,45 @@ function SciMLBase.__solve(prob::NonlinearProblem,
6564

6665
for i in 1:maxiters
6766
if alg_autodiff(alg)
68-
fx = f(x)
69-
dfdx(x) = ForwardDiff.derivative(f, x)
70-
dfx = dfdx(x)
71-
d2fx = ForwardDiff.derivative(dfdx, x)
67+
if isa(x, Number)
68+
fx = f(x)
69+
dfx = ForwardDiff.derivative(f, x)
70+
d2fx = ForwardDiff.derivative(x -> ForwardDiff.derivative(f, x), x)
71+
else
72+
fx = f(x)
73+
dfx = ForwardDiff.jacobian(f, x)
74+
d2fx = ForwardDiff.jacobian(x -> ForwardDiff.jacobian(f, x), x)
75+
ai = -(dfx \ fx)
76+
A = reshape(d2fx * ai, (n, n))
77+
bi = (dfx) \ (A * ai)
78+
ci = (ai .* ai) ./ (ai .+ (0.5 .* bi))
79+
end
7280
else
73-
fx = f(x)
74-
dfx = FiniteDiff.finite_difference_derivative(f, x, diff_type(alg), eltype(x),
75-
fx)
76-
d2fx = FiniteDiff.finite_difference_derivative(x -> FiniteDiff.finite_difference_derivative(f,
77-
x),
78-
x, diff_type(alg), eltype(x), fx)
81+
if isa(x, Number)
82+
fx = f(x)
83+
dfx = FiniteDiff.finite_difference_derivative(f, x, diff_type(alg), eltype(x))
84+
d2fx = FiniteDiff.finite_difference_derivative(x -> FiniteDiff.finite_difference_derivative(f, x), x,
85+
diff_type(alg), eltype(x))
86+
else
87+
fx = f(x)
88+
dfx = FiniteDiff.finite_difference_jacobian(f, x, diff_type(alg), eltype(x))
89+
d2fx = FiniteDiff.finite_difference_jacobian(x -> FiniteDiff.finite_difference_jacobian(f, x), x,
90+
diff_type(alg), eltype(x))
91+
ai = -(dfx \ fx)
92+
A = reshape(d2fx * ai, (n, n))
93+
bi = (dfx) \ (A * ai)
94+
ci = (ai .* ai) ./ (ai .+ (0.5 .* bi))
95+
end
7996
end
8097
iszero(fx) &&
8198
return SciMLBase.build_solution(prob, alg, x, fx; retcode = ReturnCode.Success)
82-
Δx = (2 * dfx^2 - fx * d2fx) \ (2fx * dfx)
83-
x -= Δx
99+
if isa(x, Number)
100+
Δx = (2 * dfx^2 - fx * d2fx) \ (2fx * dfx)
101+
x -= Δx
102+
else
103+
Δx = ci
104+
x += Δx
105+
end
84106
if isapprox(x, xo, atol = atol, rtol = rtol)
85107
return SciMLBase.build_solution(prob, alg, x, fx; retcode = ReturnCode.Success)
86108
end

test/basictests.jl

Lines changed: 10 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -49,10 +49,10 @@ function benchmark_scalar(f, u0)
4949
sol = (solve(probN, Halley()))
5050
end
5151

52-
# function ff(u, p)
53-
# u .* u .- 2
54-
# end
55-
# const cu0 = @SVector[1.0, 1.0]
52+
function ff(u, p)
53+
u .* u .- 2
54+
end
55+
const cu0 = @SVector[1.0, 1.0]
5656
function sf(u, p)
5757
u * u - 2
5858
end
@@ -62,6 +62,10 @@ sol = benchmark_scalar(sf, csu0)
6262
@test sol.retcode === ReturnCode.Success
6363
@test sol.u * sol.u - 2 < 1e-9
6464

65+
sol = benchmark_scalar(ff, cu0)
66+
@test sol.retcode === ReturnCode.Success
67+
@test sol.u .* sol.u .- 2 < [1e-9, 1e-9]
68+
6569
if VERSION >= v"1.7"
6670
@test (@ballocated benchmark_scalar(sf, csu0)) == 0
6771
end
@@ -122,7 +126,7 @@ using ForwardDiff
122126
f, u0 = (u, p) -> u .* u .- p, @SVector[1.0, 1.0]
123127

124128
for alg in (SimpleNewtonRaphson(), LBroyden(), Klement(), SimpleTrustRegion(),
125-
SimpleDFSane(), BROYDEN_SOLVERS...)
129+
SimpleDFSane(), Halley(), BROYDEN_SOLVERS...)
126130
g = function (p)
127131
probN = NonlinearProblem{false}(f, csu0, p)
128132
sol = solve(probN, alg, abstol = 1e-9)
@@ -221,19 +225,14 @@ probN = NonlinearProblem(f, u0)
221225

222226
for alg in (SimpleNewtonRaphson(), SimpleNewtonRaphson(; autodiff = false),
223227
SimpleTrustRegion(),
224-
SimpleTrustRegion(; autodiff = false), LBroyden(), Klement(), SimpleDFSane(),
228+
SimpleTrustRegion(; autodiff = false), Halley(), Halley(; autodiff = false), LBroyden(), Klement(), SimpleDFSane(),
225229
BROYDEN_SOLVERS...)
226230
sol = solve(probN, alg)
227231

228232
@test sol.retcode == ReturnCode.Success
229233
@test sol.u[end] sqrt(2.0)
230234
end
231235

232-
# Separate Error check for Halley; will be included in above error checks for the improved Halley
233-
f, u0 = (u, p) -> u * u - 2.0, 1.0
234-
probN = NonlinearProblem(f, u0)
235-
236-
@test solve(probN, Halley()).u sqrt(2.0)
237236

238237
for u0 in [1.0, [1, 1.0]]
239238
local f, probN, sol

0 commit comments

Comments
 (0)