Skip to content
Merged
Show file tree
Hide file tree
Changes from 3 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
3 changes: 2 additions & 1 deletion lib/BracketingNonlinearSolve/src/BracketingNonlinearSolve.jl
Original file line number Diff line number Diff line change
Expand Up @@ -17,6 +17,7 @@ include("bisection.jl")
include("brent.jl")
include("falsi.jl")
include("itp.jl")
include("muller.jl")
include("ridder.jl")

# Default Algorithm
Expand Down Expand Up @@ -44,6 +45,6 @@ end

@reexport using SciMLBase, NonlinearSolveBase

export Alefeld, Bisection, Brent, Falsi, ITP, Ridder
export Alefeld, Bisection, Brent, Falsi, Muller, ITP, Ridder

end
55 changes: 55 additions & 0 deletions lib/BracketingNonlinearSolve/src/muller.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,55 @@
"""
Muller()

Muller's method for determining a root of a univariate, scalar function. The
algorithm, described in Sec. 9.5.2 of
[Press et al. (2007)](https://numerical.recipes/book.html), requires three
initial guesses `(xᵢ₋₂, xᵢ₋₁, xᵢ)` for the root.
"""
struct Muller <: AbstractBracketingAlgorithm end

function CommonSolve.solve(prob::IntervalNonlinearProblem, alg::Muller, args...;
abstol = nothing, maxiters = 1000, kwargs...)
@assert !SciMLBase.isinplace(prob) "`Muller` only supports out-of-place problems."
xᵢ₋₂, xᵢ = prob.tspan
xᵢ₋₁ = (xᵢ₋₂ + xᵢ) / 2 # Use midpoint for middle guess
xᵢ₋₂, xᵢ₋₁, xᵢ = promote(xᵢ₋₂, xᵢ₋₁, xᵢ)
@assert xᵢ₋₂ ≠ xᵢ₋₁ ≠ xᵢ ≠ xᵢ₋₂
f = Base.Fix2(prob.f, prob.p)
fxᵢ₋₂, fxᵢ₋₁, fxᵢ = f(xᵢ₋₂), f(xᵢ₋₁), f(xᵢ)

xᵢ₊₁, fxᵢ₊₁ = xᵢ₋₂, fxᵢ₋₂

abstol = abs(NonlinearSolveBase.get_tolerance(
xᵢ₋₂, abstol, promote_type(eltype(xᵢ₋₂), eltype(xᵢ))))

for _ ∈ 1:maxiters
q = (xᵢ - xᵢ₋₁)/(xᵢ₋₁ - xᵢ₋₂)
A = q*fxᵢ - q*(1 + q)*fxᵢ₋₁ + q^2*fxᵢ₋₂
B = (2*q + 1)*fxᵢ - (1 + q)^2*fxᵢ₋₁ + q^2*fxᵢ₋₂
C = (1 + q)*fxᵢ

denom₊ = B + √(B^2 - 4*A*C)
denom₋ = B - √(B^2 - 4*A*C)

if abs(denom₊) ≥ abs(denom₋)
xᵢ₊₁ = xᵢ - (xᵢ - xᵢ₋₁)*2*C/denom₊
else
xᵢ₊₁ = xᵢ - (xᵢ - xᵢ₋₁)*2*C/denom₋
end

fxᵢ₊₁ = f(xᵢ₊₁)

# Termination Check
if abstol ≥ abs(fxᵢ₊₁)
return SciMLBase.build_solution(prob, alg, xᵢ₊₁, fxᵢ₊₁;
retcode = ReturnCode.Success)
end

xᵢ₋₂, xᵢ₋₁, xᵢ = xᵢ₋₁, xᵢ, xᵢ₊₁
fxᵢ₋₂, fxᵢ₋₁, fxᵢ = fxᵢ₋₁, fxᵢ, fxᵢ₊₁
end

return SciMLBase.build_solution(prob, alg, xᵢ₊₁, fxᵢ₊₁;
retcode = ReturnCode.MaxIters)
end
72 changes: 72 additions & 0 deletions lib/BracketingNonlinearSolve/test/muller_tests.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,72 @@
@testitem "Muller" begin
@testset "Quadratic function" begin
f(u, p) = u^2 - p

tspan = (10.0, 30.0)
p = 612.0
prob = IntervalNonlinearProblem{false}(f, tspan, p)
sol = solve(prob, Muller())

@test sol.u ≈ √612

tspan = (-10.0, -30.0)
prob = IntervalNonlinearProblem{false}(f, tspan, p)
sol = solve(prob, Muller())

@test sol.u ≈ -√612
end

@testset "Sine function" begin
f(u, p) = sin(u)

tspan = (1.0, 3.0)
prob = IntervalNonlinearProblem{false}(f, tspan)
sol = solve(prob, Muller())

@test sol.u ≈ π

tspan = (2.0, 6.0)
prob = IntervalNonlinearProblem{false}(f, tspan)
sol = solve(prob, Muller())

@test sol.u ≈ 2*π
end

@testset "Exponential-sine function" begin
f(u, p) = exp(-u)*sin(u)

tspan = (-2.0, -4.0)
prob = IntervalNonlinearProblem{false}(f, tspan)
sol = solve(prob, Muller())

@test sol.u ≈ -π

tspan = (-3.0, 1.0)
prob = IntervalNonlinearProblem{false}(f, tspan)
sol = solve(prob, Muller())

@test sol.u ≈ 0 atol = 1e-15

tspan = (-1.0, 1.0)
prob = IntervalNonlinearProblem{false}(f, tspan)
sol = solve(prob, Muller())

@test sol.u ≈ π
end

@testset "Complex roots" begin
f(u, p) = u^3 - 1

tspan = (-1.0, 1.0*im)
prob = IntervalNonlinearProblem{false}(f, tspan)
sol = solve(prob, Muller())

@test sol.u ≈ (-1 + √3*im)/2

tspan = (-1.0, -1.0*im)
prob = IntervalNonlinearProblem{false}(f, tspan)
sol = solve(prob, Muller())

@test sol.u ≈ (-1 - √3*im)/2
end
end
Loading