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

Commit f8dc450

Browse files
committed
Implemented and tested the Klement solver
1 parent 325547e commit f8dc450

File tree

3 files changed

+86
-4
lines changed

3 files changed

+86
-4
lines changed

src/SimpleNonlinearSolve.jl

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -20,6 +20,7 @@ include("falsi.jl")
2020
include("raphson.jl")
2121
include("ad.jl")
2222
include("broyden.jl")
23+
include("klement.jl")
2324

2425
import SnoopPrecompile
2526

@@ -46,6 +47,6 @@ SnoopPrecompile.@precompile_all_calls begin for T in (Float32, Float64)
4647
end end
4748

4849
# DiffEq styled algorithms
49-
export Bisection, Broyden, Falsi, SimpleNewtonRaphson
50+
export Bisection, Broyden, Falsi, Klement, SimpleNewtonRaphson
5051

5152
end # module

src/klement.jl

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

test/basictests.jl

Lines changed: 19 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -35,13 +35,24 @@ sol = benchmark_scalar(sf, csu0)
3535
@test sol.u * sol.u - 2 < 1e-9
3636
@test (@ballocated benchmark_scalar(sf, csu0)) == 0
3737

38+
# Klement
39+
function benchmark_scalar(f, u0)
40+
probN = NonlinearProblem{false}(f, u0)
41+
sol = (solve(probN, Klement()))
42+
end
43+
44+
sol = benchmark_scalar(sf, csu0)
45+
@test sol.retcode === ReturnCode.Success
46+
@test sol.u * sol.u - 2 < 1e-9
47+
@test (@ballocated benchmark_scalar(sf, csu0)) == 0
48+
3849
# AD Tests
3950
using ForwardDiff
4051

4152
# Immutable
4253
f, u0 = (u, p) -> u .* u .- p, @SVector[1.0, 1.0]
4354

44-
for alg in [SimpleNewtonRaphson(), Broyden()]
55+
for alg in [SimpleNewtonRaphson(), Broyden(), Klement()]
4556
g = function (p)
4657
probN = NonlinearProblem{false}(f, csu0, p)
4758
sol = solve(probN, alg, tol = 1e-9)
@@ -56,7 +67,7 @@ end
5667

5768
# Scalar
5869
f, u0 = (u, p) -> u * u - p, 1.0
59-
for alg in [SimpleNewtonRaphson(), Broyden()]
70+
for alg in [SimpleNewtonRaphson(), Broyden(), Klement()]
6071
g = function (p)
6172
probN = NonlinearProblem{false}(f, oftype(p, u0), p)
6273
sol = solve(probN, alg)
@@ -97,7 +108,7 @@ for alg in [Bisection(), Falsi()]
97108
@test ForwardDiff.jacobian(g, p) ForwardDiff.jacobian(t, p)
98109
end
99110

100-
for alg in [SimpleNewtonRaphson(), Broyden()]
111+
for alg in [SimpleNewtonRaphson(), Broyden(), Klement()]
101112
global g, p
102113
g = function (p)
103114
probN = NonlinearProblem{false}(f, 0.5, p)
@@ -120,6 +131,9 @@ probN = NonlinearProblem(f, u0)
120131
@test solve(probN, Broyden()).u[end] sqrt(2.0)
121132
@test solve(probN, Broyden(); immutable = false).u[end] sqrt(2.0)
122133

134+
@test solve(probN, Klement()).u[end] sqrt(2.0)
135+
@test solve(probN, Klement(); immutable = false).u[end] sqrt(2.0)
136+
123137
for u0 in [1.0, [1, 1.0]]
124138
local f, probN, sol
125139
f = (u, p) -> u .* u .- 2.0
@@ -131,6 +145,8 @@ for u0 in [1.0, [1, 1.0]]
131145
@test solve(probN, SimpleNewtonRaphson(; autodiff = false)).u sol
132146

133147
@test solve(probN, Broyden()).u sol
148+
149+
@test solve(probN, Klement()).u sol
134150
end
135151

136152
# Bisection Tests

0 commit comments

Comments
 (0)