Skip to content

Commit 69ba15c

Browse files
committed
Add Muller's method
1 parent b182ae1 commit 69ba15c

File tree

3 files changed

+110
-0
lines changed

3 files changed

+110
-0
lines changed

lib/SimpleNonlinearSolve/src/SimpleNonlinearSolve.jl

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -54,6 +54,7 @@ include("klement.jl")
5454
include("lbroyden.jl")
5555
include("raphson.jl")
5656
include("trust_region.jl")
57+
include("muller.jl")
5758

5859
# By Pass the highlevel checks for NonlinearProblem for Simple Algorithms
5960
function CommonSolve.solve(
@@ -166,6 +167,7 @@ export SimpleBroyden, SimpleKlement, SimpleLimitedMemoryBroyden
166167
export SimpleDFSane
167168
export SimpleGaussNewton, SimpleNewtonRaphson, SimpleTrustRegion
168169
export SimpleHalley
170+
export SimpleMuller
169171

170172
export solve
171173

Lines changed: 52 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,52 @@
1+
"""
2+
SimpleMuller()
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 `(xᵢ₋₂, xᵢ₋₁, xᵢ)` for the root.
8+
"""
9+
struct SimpleMuller <: AbstractSimpleNonlinearSolveAlgorithm end
10+
11+
function SciMLBase.solve(prob::NonlinearProblem, alg::SimpleMuller, args...;
12+
abstol = 1e-3, maxiters = 1000, kwargs...)
13+
@assert !isinplace(prob) "`SimpleMuller` only supports OOP problems."
14+
@assert length(prob.u0) == 3 "`SimpleMuller` requires three initial guesses."
15+
xᵢ₋₂, xᵢ₋₁, xᵢ = prob.u0
16+
xᵢ₋₂, xᵢ₋₁, xᵢ = promote(xᵢ₋₂, xᵢ₋₁, xᵢ)
17+
@assert xᵢ₋₂ xᵢ₋₁ xᵢ xᵢ₋₂
18+
f = Base.Fix2(prob.f, prob.p)
19+
fxᵢ₋₂, fxᵢ₋₁, fxᵢ = f(xᵢ₋₂), f(xᵢ₋₁), f(xᵢ)
20+
21+
xᵢ₊₁, fxᵢ₊₁ = xᵢ₋₂, fxᵢ₋₂
22+
23+
for _ 1:maxiters
24+
q = (xᵢ - xᵢ₋₁)/(xᵢ₋₁ - xᵢ₋₂)
25+
A = q*fxᵢ - q*(1 + q)*fxᵢ₋₁ + q^2*fxᵢ₋₂
26+
B = (2*q + 1)*fxᵢ - (1 + q)^2*fxᵢ₋₁ + q^2*fxᵢ₋₂
27+
C = (1 + q)*fxᵢ
28+
29+
denom₊ = B + (B^2 - 4*A*C)
30+
denom₋ = B - (B^2 - 4*A*C)
31+
32+
if abs(denom₊) abs(denom₋)
33+
xᵢ₊₁ = xᵢ - (xᵢ - xᵢ₋₁)*2*C/denom₊
34+
else
35+
xᵢ₊₁ = xᵢ - (xᵢ - xᵢ₋₁)*2*C/denom₋
36+
end
37+
38+
fxᵢ₊₁ = f(xᵢ₊₁)
39+
40+
# Termination Check
41+
if abstol abs(fxᵢ₊₁)
42+
return SciMLBase.build_solution(prob, alg, xᵢ₊₁, fxᵢ₊₁;
43+
retcode = ReturnCode.Success)
44+
end
45+
46+
xᵢ₋₂, xᵢ₋₁, xᵢ = xᵢ₋₁, xᵢ, xᵢ₊₁
47+
fxᵢ₋₂, fxᵢ₋₁, fxᵢ = fxᵢ₋₁, fxᵢ, fxᵢ₊₁
48+
end
49+
50+
return SciMLBase.build_solution(prob, alg, xᵢ₊₁, fxᵢ₊₁;
51+
retcode = ReturnCode.MaxIters)
52+
end
Lines changed: 56 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,56 @@
1+
@testitem "SimpleMuller" begin
2+
@testset "Quadratic function" begin
3+
f(u, p) = u^2 - p
4+
5+
u0 = (10.0, 20.0, 30.0)
6+
p = 612.0
7+
prob = NonlinearProblem{false}(f, u0, p)
8+
sol = solve(prob, SimpleMuller())
9+
10+
@test sol.u 612
11+
12+
u0 = (-10.0, -20.0, -30.0)
13+
prob = NonlinearProblem{false}(f, u0, p)
14+
sol = solve(prob, SimpleMuller())
15+
16+
@test sol.u -√612
17+
end
18+
19+
@testset "Sine function" begin
20+
f(u, p) = sin(u)
21+
22+
u0 = (1.0, 2.0, 3.0)
23+
prob = NonlinearProblem{false}(f, u0)
24+
sol = solve(prob, SimpleMuller())
25+
26+
@test sol.u π
27+
28+
u0 = (2.0, 4.0, 6.0)
29+
prob = NonlinearProblem{false}(f, u0)
30+
sol = solve(prob, SimpleMuller())
31+
32+
@test sol.u 2*π
33+
end
34+
35+
@testset "Exponential-sine function" begin
36+
f(u, p) = exp(-u)*sin(u)
37+
38+
u0 = (-2.0, -3.0, -4.0)
39+
prob = NonlinearProblem{false}(f, u0)
40+
sol = solve(prob, SimpleMuller())
41+
42+
@test sol.u -π
43+
44+
u0 = (-1.0, 0.0, 1/2)
45+
prob = NonlinearProblem{false}(f, u0)
46+
sol = solve(prob, SimpleMuller())
47+
48+
@test sol.u 0
49+
50+
u0 = (-1.0, 0.0, 1.0)
51+
prob = NonlinearProblem{false}(f, u0)
52+
sol = solve(prob, SimpleMuller())
53+
54+
@test sol.u π
55+
end
56+
end

0 commit comments

Comments
 (0)