Skip to content

Commit 44b110d

Browse files
Merge pull request #2525 from SciML/precision
Fix precision mixing in RodasTableau and better test precision mix
2 parents bdef41c + e379c46 commit 44b110d

File tree

3 files changed

+24
-1
lines changed

3 files changed

+24
-1
lines changed

lib/OrdinaryDiffEqRosenbrock/src/rosenbrock_tableaus.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -225,7 +225,7 @@ end
225225
struct RodasTableau{T, T2}
226226
A::Matrix{T}
227227
C::Matrix{T}
228-
gamma::T
228+
gamma::T2
229229
c::Vector{T2}
230230
d::Vector{T}
231231
H::Matrix{T}

test/interface/precision_mixing.jl

Lines changed: 22 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,22 @@
1+
using OrdinaryDiffEq
2+
function ODE(du, u, t, R, K)
3+
du .= u
4+
end
5+
params = BigFloat[1. 0.91758707304098; 1.48439909482661 1.]
6+
u0 = BigFloat[0.1, 0.1]
7+
tspan = (1.0, 31.0)
8+
R = BigFloat[0.443280390004304303, 1.172917082211452]
9+
K = BigFloat[13.470600276901400604, 12.52980757005]
10+
ODE_ = (du, u, params, t) -> ODE(du, u, t, R, K)
11+
odeProblem = ODEProblem(ODE_, u0, tspan, params)
12+
for alg in [AutoVern8(Rodas5(), nonstifftol = 11 / 10)
13+
FBDF()
14+
QNDF()
15+
Tsit5()
16+
Rodas5P()
17+
TRBDF2()
18+
KenCarp4()
19+
RadauIIA5()
20+
]
21+
Solution = solve(odeProblem, alg, saveat = 1, abstol = 1.e-12, reltol = 1.e-6)
22+
end

test/runtests.jl

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -85,6 +85,7 @@ end
8585
if !is_APPVEYOR && (GROUP == "All" || GROUP == "InterfaceIV" || GROUP == "Interface")
8686
@time @safetestset "Autodiff Error Tests" include("interface/autodiff_error_tests.jl")
8787
@time @safetestset "Ambiguity Tests" include("interface/ambiguity_tests.jl")
88+
@time @safetestset "Precision Mixing Tests" include("interface/precision_mixing.jl")
8889
@time @safetestset "Sized Matrix Tests" include("interface/sized_matrix_tests.jl")
8990
@time @safetestset "Second Order with First Order Solver Tests" include("interface/second_order_with_first_order_solvers.jl")
9091
end

0 commit comments

Comments
 (0)