Skip to content

Commit f5d3ef7

Browse files
authored
Add Muller's method (#568)
* Add Muller's method * Modify Muller to be an IntervalNonlinearProblem * Fix bug with complex roots * Incorporate Muller into test suite Add `left` and `right` fields to solution * Add functionality to specify `middle` guess
1 parent d56a975 commit f5d3ef7

File tree

4 files changed

+168
-5
lines changed

4 files changed

+168
-5
lines changed

lib/BracketingNonlinearSolve/src/BracketingNonlinearSolve.jl

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -17,6 +17,7 @@ include("bisection.jl")
1717
include("brent.jl")
1818
include("falsi.jl")
1919
include("itp.jl")
20+
include("muller.jl")
2021
include("ridder.jl")
2122

2223
# Default Algorithm
@@ -44,6 +45,6 @@ end
4445

4546
@reexport using SciMLBase, NonlinearSolveBase
4647

47-
export Alefeld, Bisection, Brent, Falsi, ITP, Ridder
48+
export Alefeld, Bisection, Brent, Falsi, ITP, Muller, Ridder
4849

4950
end
Lines changed: 66 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,66 @@
1+
"""
2+
Muller(; middle = nothing)
3+
4+
Muller's method for determining a root of a univariate, scalar function. The
5+
algorithm, described in Sec. 9.5.2 of
6+
[Press et al. (2007)](https://numerical.recipes/book.html), requires three
7+
initial guesses `(left, middle, right)` for the root.
8+
9+
### Keyword Arguments
10+
11+
- `middle`: the initial guess for the middle point. If not provided, the
12+
midpoint of the interval `(left, right)` is used.
13+
"""
14+
struct Muller{T} <: AbstractBracketingAlgorithm
15+
middle::T
16+
end
17+
18+
Muller() = Muller(nothing)
19+
20+
function CommonSolve.solve(prob::IntervalNonlinearProblem, alg::Muller, args...;
21+
abstol = nothing, maxiters = 1000, kwargs...)
22+
@assert !SciMLBase.isinplace(prob) "`Muller` only supports out-of-place problems."
23+
xᵢ₋₂, xᵢ = prob.tspan
24+
xᵢ₋₁ = isnothing(alg.middle) ? (xᵢ₋₂ + xᵢ) / 2 : alg.middle
25+
xᵢ₋₂, xᵢ₋₁, xᵢ = promote(xᵢ₋₂, xᵢ₋₁, xᵢ)
26+
@assert xᵢ₋₂ xᵢ₋₁ xᵢ xᵢ₋₂
27+
f = Base.Fix2(prob.f, prob.p)
28+
fxᵢ₋₂, fxᵢ₋₁, fxᵢ = f(xᵢ₋₂), f(xᵢ₋₁), f(xᵢ)
29+
30+
xᵢ₊₁, fxᵢ₊₁ = xᵢ₋₂, fxᵢ₋₂
31+
32+
abstol = abs(NonlinearSolveBase.get_tolerance(
33+
xᵢ₋₂, abstol, promote_type(eltype(xᵢ₋₂), eltype(xᵢ))))
34+
35+
for _ 1:maxiters
36+
q = (xᵢ - xᵢ₋₁)/(xᵢ₋₁ - xᵢ₋₂)
37+
A = q*fxᵢ - q*(1 + q)*fxᵢ₋₁ + q^2*fxᵢ₋₂
38+
B = (2*q + 1)*fxᵢ - (1 + q)^2*fxᵢ₋₁ + q^2*fxᵢ₋₂
39+
C = (1 + q)*fxᵢ
40+
41+
denom₊ = B + (B^2 - 4*A*C)
42+
denom₋ = B - (B^2 - 4*A*C)
43+
44+
if abs(denom₊) abs(denom₋)
45+
xᵢ₊₁ = xᵢ - (xᵢ - xᵢ₋₁)*2*C/denom₊
46+
else
47+
xᵢ₊₁ = xᵢ - (xᵢ - xᵢ₋₁)*2*C/denom₋
48+
end
49+
50+
fxᵢ₊₁ = f(xᵢ₊₁)
51+
52+
# Termination Check
53+
if abstol abs(fxᵢ₊₁)
54+
return SciMLBase.build_solution(prob, alg, xᵢ₊₁, fxᵢ₊₁;
55+
retcode = ReturnCode.Success,
56+
left = xᵢ₊₁, right = xᵢ₊₁)
57+
end
58+
59+
xᵢ₋₂, xᵢ₋₁, xᵢ = xᵢ₋₁, xᵢ, xᵢ₊₁
60+
fxᵢ₋₂, fxᵢ₋₁, fxᵢ = fxᵢ₋₁, fxᵢ, fxᵢ₊₁
61+
end
62+
63+
return SciMLBase.build_solution(prob, alg, xᵢ₊₁, fxᵢ₊₁;
64+
retcode = ReturnCode.MaxIters,
65+
left = xᵢ₊₁, right = xᵢ₊₁)
66+
end
Lines changed: 96 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,96 @@
1+
@testitem "Muller" begin
2+
f(u, p) = u^2 - p
3+
g(u, p) = sin(u)
4+
h(u, p) = exp(-u)*sin(u)
5+
i(u, p) = u^3 - 1
6+
7+
@testset "Quadratic function" begin
8+
tspan = (10.0, 30.0)
9+
p = 612.0
10+
prob = IntervalNonlinearProblem{false}(f, tspan, p)
11+
sol = solve(prob, Muller())
12+
13+
@test sol.u 612
14+
15+
tspan = (-10.0, -30.0)
16+
prob = IntervalNonlinearProblem{false}(f, tspan, p)
17+
sol = solve(prob, Muller())
18+
19+
@test sol.u -√612
20+
end
21+
22+
@testset "Sine function" begin
23+
tspan = (1.0, 3.0)
24+
prob = IntervalNonlinearProblem{false}(g, tspan)
25+
sol = solve(prob, Muller())
26+
27+
@test sol.u π
28+
29+
tspan = (2.0, 6.0)
30+
prob = IntervalNonlinearProblem{false}(g, tspan)
31+
sol = solve(prob, Muller())
32+
33+
@test sol.u 2*π
34+
end
35+
36+
@testset "Exponential-sine function" begin
37+
tspan = (-2.0, -4.0)
38+
prob = IntervalNonlinearProblem{false}(h, tspan)
39+
sol = solve(prob, Muller())
40+
41+
@test sol.u -π
42+
43+
tspan = (-3.0, 1.0)
44+
prob = IntervalNonlinearProblem{false}(h, tspan)
45+
sol = solve(prob, Muller())
46+
47+
@test sol.u 0 atol = 1e-15
48+
49+
tspan = (-1.0, 1.0)
50+
prob = IntervalNonlinearProblem{false}(h, tspan)
51+
sol = solve(prob, Muller())
52+
53+
@test sol.u π
54+
end
55+
56+
@testset "Complex roots" begin
57+
tspan = (-1.0, 1.0*im)
58+
prob = IntervalNonlinearProblem{false}(i, tspan)
59+
sol = solve(prob, Muller())
60+
61+
@test sol.u (-1 + 3*im)/2
62+
63+
tspan = (-1.0, -1.0*im)
64+
prob = IntervalNonlinearProblem{false}(i, tspan)
65+
sol = solve(prob, Muller())
66+
67+
@test sol.u (-1 - 3*im)/2
68+
end
69+
70+
@testset "Middle" begin
71+
tspan = (10.0, 30.0)
72+
p = 612.0
73+
prob = IntervalNonlinearProblem{false}(f, tspan, p)
74+
sol = solve(prob, Muller(20.0))
75+
76+
@test sol.u 612
77+
78+
tspan = (1.0, 3.0)
79+
prob = IntervalNonlinearProblem{false}(g, tspan)
80+
sol = solve(prob, Muller(2.0))
81+
82+
@test sol.u π
83+
84+
tspan = (-2.0, -4.0)
85+
prob = IntervalNonlinearProblem{false}(h, tspan)
86+
sol = solve(prob, Muller(-3.0))
87+
88+
@test sol.u -π
89+
90+
tspan = (-1.0, 1.0*im)
91+
prob = IntervalNonlinearProblem{false}(i, tspan)
92+
sol = solve(prob, Muller(0.0))
93+
94+
@test sol.u (-1 + 3*im)/2
95+
end
96+
end

lib/BracketingNonlinearSolve/test/rootfind_tests.jl

Lines changed: 4 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -7,7 +7,7 @@ end
77
@testitem "Interval Nonlinear Problems" setup=[RootfindingTestSnippet] tags=[:core] begin
88
using ForwardDiff
99

10-
@testset for alg in (Bisection(), Falsi(), Ridder(), Brent(), ITP(), Alefeld(), nothing)
10+
@testset for alg in (Alefeld(), Bisection(), Brent(), Falsi(), ITP(), Muller(), Ridder(), nothing)
1111
tspan = (1.0, 20.0)
1212

1313
function g(p)
@@ -52,7 +52,7 @@ end
5252
prob = IntervalNonlinearProblem(quadratic_f, (1.0, 20.0), 2.0)
5353
ϵ = eps(Float64) # least possible tol for all methods
5454

55-
@testset for alg in (Bisection(), Falsi(), ITP(), nothing)
55+
@testset for alg in (Bisection(), Falsi(), ITP(), Muller(), nothing)
5656
@testset for abstol in [0.1, 0.01, 0.001, 0.0001, 1e-5, 1e-6]
5757
sol = solve(prob, alg; abstol)
5858
result_tol = abs(sol.u - sqrt(2))
@@ -62,7 +62,7 @@ end
6262
end
6363
end
6464

65-
@testset for alg in (Ridder(), Brent())
65+
@testset for alg in (Brent(), Ridder())
6666
# Ridder and Brent converge rapidly so as we lower tolerance below 0.01, it
6767
# converges with max precision to the solution
6868
@testset for abstol in [0.1]
@@ -76,7 +76,7 @@ end
7676
end
7777

7878
@testitem "Flipped Signs and Reversed Tspan" setup=[RootfindingTestSnippet] tags=[:core] begin
79-
@testset for alg in (Alefeld(), Bisection(), Falsi(), Brent(), ITP(), Ridder(), nothing)
79+
@testset for alg in (Alefeld(), Bisection(), Brent(), Falsi(), ITP(), Muller(), Ridder(), nothing)
8080
f1(u, p) = u * u - p
8181
f2(u, p) = p - u * u
8282

0 commit comments

Comments
 (0)