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

Commit 325547e

Browse files
Merge pull request #14 from CCsimon123/main
Add Broyden solver
2 parents 92b2eb2 + aca04ad commit 325547e

File tree

4 files changed

+104
-30
lines changed

4 files changed

+104
-30
lines changed

Project.toml

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -7,6 +7,7 @@ version = "0.1.4"
77
ArrayInterfaceCore = "30b0a656-2188-435a-8636-2ec0e6a096e2"
88
FiniteDiff = "6a86dc24-6348-571c-b903-95158fe2bd41"
99
ForwardDiff = "f6369f11-7733-5829-9624-2563aa707210"
10+
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
1011
Reexport = "189a3867-3050-52da-a836-e630ba90ab69"
1112
SciMLBase = "0bca4576-84f4-4d90-8ffe-ffa030f20462"
1213
SnoopPrecompile = "66db9d55-30c0-4569-8b51-7e840670fc0c"

src/SimpleNonlinearSolve.jl

Lines changed: 4 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -4,6 +4,7 @@ using Reexport
44
using FiniteDiff, ForwardDiff
55
using ForwardDiff: Dual
66
using StaticArraysCore
7+
using LinearAlgebra
78
import ArrayInterfaceCore
89

910
@reexport using SciMLBase
@@ -18,12 +19,13 @@ include("bisection.jl")
1819
include("falsi.jl")
1920
include("raphson.jl")
2021
include("ad.jl")
22+
include("broyden.jl")
2123

2224
import SnoopPrecompile
2325

2426
SnoopPrecompile.@precompile_all_calls begin for T in (Float32, Float64)
2527
prob_no_brack = NonlinearProblem{false}((u, p) -> u .* u .- p, T(0.1), T(2))
26-
for alg in (SimpleNewtonRaphson,)
28+
for alg in (SimpleNewtonRaphson, Broyden)
2729
solve(prob_no_brack, alg(), tol = T(1e-2))
2830
end
2931

@@ -44,6 +46,6 @@ SnoopPrecompile.@precompile_all_calls begin for T in (Float32, Float64)
4446
end end
4547

4648
# DiffEq styled algorithms
47-
export Bisection, Falsi, SimpleNewtonRaphson
49+
export Bisection, Broyden, Falsi, SimpleNewtonRaphson
4850

4951
end # module

src/broyden.jl

Lines changed: 52 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,52 @@
1+
"""
2+
```julia
3+
Broyden()
4+
```
5+
6+
A low-overhead implementation of Broyden. This method is non-allocating on scalar
7+
and static array problems.
8+
"""
9+
struct Broyden <: AbstractSimpleNonlinearSolveAlgorithm end
10+
11+
function SciMLBase.solve(prob::NonlinearProblem,
12+
alg::Broyden, args...; abstol = nothing,
13+
reltol = nothing,
14+
maxiters = 1000, kwargs...)
15+
f = Base.Fix2(prob.f, prob.p)
16+
x = float(prob.u0)
17+
fₙ = f(x)
18+
T = eltype(x)
19+
J⁻¹ = ArrayInterfaceCore.zeromatrix(x) + I
20+
21+
if SciMLBase.isinplace(prob)
22+
error("Broyden currently only supports out-of-place nonlinear problems")
23+
end
24+
25+
atol = abstol !== nothing ? abstol :
26+
real(oneunit(eltype(T))) * (eps(real(one(eltype(T)))))^(4 // 5)
27+
rtol = reltol !== nothing ? reltol : eps(real(one(eltype(T))))^(4 // 5)
28+
29+
xₙ = x
30+
xₙ₋₁ = x
31+
fₙ₋₁ = fₙ
32+
for _ in 1:maxiters
33+
xₙ = xₙ₋₁ - J⁻¹ * fₙ₋₁
34+
fₙ = f(xₙ)
35+
Δxₙ = xₙ - xₙ₋₁
36+
Δfₙ = fₙ - fₙ₋₁
37+
J⁻¹ += ((Δxₙ - J⁻¹ * Δfₙ) ./ (Δxₙ' * J⁻¹ * Δfₙ)) * (Δxₙ' * J⁻¹)
38+
39+
iszero(fₙ) &&
40+
return SciMLBase.build_solution(prob, alg, xₙ, fₙ;
41+
retcode = ReturnCode.Success)
42+
43+
if isapprox(xₙ, xₙ₋₁, atol = atol, rtol = rtol)
44+
return SciMLBase.build_solution(prob, alg, xₙ, fₙ;
45+
retcode = ReturnCode.Success)
46+
end
47+
xₙ₋₁ = xₙ
48+
fₙ₋₁ = fₙ
49+
end
50+
51+
return SciMLBase.build_solution(prob, alg, xₙ, fₙ; retcode = ReturnCode.MaxIters)
52+
end

test/basictests.jl

Lines changed: 47 additions & 28 deletions
Original file line numberDiff line numberDiff line change
@@ -3,6 +3,7 @@ using StaticArrays
33
using BenchmarkTools
44
using Test
55

6+
# SimpleNewtonRaphson
67
function benchmark_scalar(f, u0)
78
probN = NonlinearProblem{false}(f, u0)
89
sol = (solve(probN, SimpleNewtonRaphson()))
@@ -23,38 +24,49 @@ sol = benchmark_scalar(sf, csu0)
2324

2425
@test (@ballocated benchmark_scalar(sf, csu0)) == 0
2526

27+
# Broyden
28+
function benchmark_scalar(f, u0)
29+
probN = NonlinearProblem{false}(f, u0)
30+
sol = (solve(probN, Broyden()))
31+
end
32+
33+
sol = benchmark_scalar(sf, csu0)
34+
@test sol.retcode === ReturnCode.Success
35+
@test sol.u * sol.u - 2 < 1e-9
36+
@test (@ballocated benchmark_scalar(sf, csu0)) == 0
37+
2638
# AD Tests
2739
using ForwardDiff
2840

2941
# Immutable
3042
f, u0 = (u, p) -> u .* u .- p, @SVector[1.0, 1.0]
3143

32-
g = function (p)
33-
probN = NonlinearProblem{false}(f, csu0, p)
34-
sol = solve(probN, SimpleNewtonRaphson(), tol = 1e-9)
35-
return sol.u[end]
36-
end
44+
for alg in [SimpleNewtonRaphson(), Broyden()]
45+
g = function (p)
46+
probN = NonlinearProblem{false}(f, csu0, p)
47+
sol = solve(probN, alg, tol = 1e-9)
48+
return sol.u[end]
49+
end
3750

38-
for p in 1.0:0.1:100.0
39-
@test g(p) sqrt(p)
40-
@test ForwardDiff.derivative(g, p) 1 / (2 * sqrt(p))
51+
for p in 1.1:0.1:100.0
52+
@test g(p) sqrt(p)
53+
@test ForwardDiff.derivative(g, p) 1 / (2 * sqrt(p))
54+
end
4155
end
4256

4357
# Scalar
4458
f, u0 = (u, p) -> u * u - p, 1.0
59+
for alg in [SimpleNewtonRaphson(), Broyden()]
60+
g = function (p)
61+
probN = NonlinearProblem{false}(f, oftype(p, u0), p)
62+
sol = solve(probN, alg)
63+
return sol.u
64+
end
4565

46-
# SimpleNewtonRaphson
47-
g = function (p)
48-
probN = NonlinearProblem{false}(f, oftype(p, u0), p)
49-
sol = solve(probN, SimpleNewtonRaphson())
50-
return sol.u
51-
end
52-
53-
@test ForwardDiff.derivative(g, 1.0) 0.5
54-
55-
for p in 1.1:0.1:100.0
56-
@test g(p) sqrt(p)
57-
@test ForwardDiff.derivative(g, p) 1 / (2 * sqrt(p))
66+
for p in 1.1:0.1:100.0
67+
@test g(p) sqrt(p)
68+
@test ForwardDiff.derivative(g, p) 1 / (2 * sqrt(p))
69+
end
5870
end
5971

6072
tspan = (1.0, 20.0)
@@ -77,24 +89,26 @@ for alg in [Bisection(), Falsi()]
7789
global g, p
7890
g = function (p)
7991
probN = IntervalNonlinearProblem{false}(f, tspan, p)
80-
sol = solve(probN, Bisection())
92+
sol = solve(probN, alg)
8193
return [sol.left]
8294
end
8395

8496
@test g(p) [sqrt(p[2] / p[1])]
8597
@test ForwardDiff.jacobian(g, p) ForwardDiff.jacobian(t, p)
8698
end
8799

88-
gnewton = function (p)
89-
probN = NonlinearProblem{false}(f, 0.5, p)
90-
sol = solve(probN, SimpleNewtonRaphson())
91-
return [sol.u]
100+
for alg in [SimpleNewtonRaphson(), Broyden()]
101+
global g, p
102+
g = function (p)
103+
probN = NonlinearProblem{false}(f, 0.5, p)
104+
sol = solve(probN, alg)
105+
return [sol.u]
106+
end
107+
@test g(p) [sqrt(p[2] / p[1])]
108+
@test ForwardDiff.jacobian(g, p) ForwardDiff.jacobian(t, p)
92109
end
93-
@test gnewton(p) [sqrt(p[2] / p[1])]
94-
@test ForwardDiff.jacobian(gnewton, p) ForwardDiff.jacobian(t, p)
95110

96111
# Error Checks
97-
98112
f, u0 = (u, p) -> u .* u .- 2.0, @SVector[1.0, 1.0]
99113
probN = NonlinearProblem(f, u0)
100114

@@ -103,6 +117,9 @@ probN = NonlinearProblem(f, u0)
103117
@test solve(probN, SimpleNewtonRaphson(; autodiff = false)).u[end] sqrt(2.0)
104118
@test solve(probN, SimpleNewtonRaphson(; autodiff = false)).u[end] sqrt(2.0)
105119

120+
@test solve(probN, Broyden()).u[end] sqrt(2.0)
121+
@test solve(probN, Broyden(); immutable = false).u[end] sqrt(2.0)
122+
106123
for u0 in [1.0, [1, 1.0]]
107124
local f, probN, sol
108125
f = (u, p) -> u .* u .- 2.0
@@ -112,6 +129,8 @@ for u0 in [1.0, [1, 1.0]]
112129
@test solve(probN, SimpleNewtonRaphson()).u sol
113130
@test solve(probN, SimpleNewtonRaphson()).u sol
114131
@test solve(probN, SimpleNewtonRaphson(; autodiff = false)).u sol
132+
133+
@test solve(probN, Broyden()).u sol
115134
end
116135

117136
# Bisection Tests

0 commit comments

Comments
 (0)