Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix precision mixing in RodasTableau and better test precision mix #2525

Merged
merged 2 commits into from
Nov 12, 2024
Merged
Show file tree
Hide file tree
Changes from all 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
2 changes: 1 addition & 1 deletion lib/OrdinaryDiffEqRosenbrock/src/rosenbrock_tableaus.jl
Original file line number Diff line number Diff line change
Expand Up @@ -225,7 +225,7 @@ end
struct RodasTableau{T, T2}
A::Matrix{T}
C::Matrix{T}
gamma::T
gamma::T2
c::Vector{T2}
d::Vector{T}
H::Matrix{T}
Expand Down
22 changes: 22 additions & 0 deletions test/interface/precision_mixing.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,22 @@
using OrdinaryDiffEq
function ODE(du, u, t, R, K)
du .= u
end
params = BigFloat[1. 0.91758707304098; 1.48439909482661 1.]
u0 = BigFloat[0.1, 0.1]
tspan = (1.0, 31.0)
R = BigFloat[0.443280390004304303, 1.172917082211452]
K = BigFloat[13.470600276901400604, 12.52980757005]
ODE_ = (du, u, params, t) -> ODE(du, u, t, R, K)
odeProblem = ODEProblem(ODE_, u0, tspan, params)
for alg in [AutoVern8(Rodas5(), nonstifftol = 11 / 10)
FBDF()
QNDF()
Tsit5()
Rodas5P()
TRBDF2()
KenCarp4()
RadauIIA5()
]
Solution = solve(odeProblem, alg, saveat = 1, abstol = 1.e-12, reltol = 1.e-6)
end
1 change: 1 addition & 0 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -85,6 +85,7 @@ end
if !is_APPVEYOR && (GROUP == "All" || GROUP == "InterfaceIV" || GROUP == "Interface")
@time @safetestset "Autodiff Error Tests" include("interface/autodiff_error_tests.jl")
@time @safetestset "Ambiguity Tests" include("interface/ambiguity_tests.jl")
@time @safetestset "Precision Mixing Tests" include("interface/precision_mixing.jl")
@time @safetestset "Sized Matrix Tests" include("interface/sized_matrix_tests.jl")
@time @safetestset "Second Order with First Order Solver Tests" include("interface/second_order_with_first_order_solvers.jl")
end
Expand Down
Loading