Skip to content

Commit fdafa9c

Browse files
Merge pull request #1048 from prbzrg/new-AlgorithmInterpretation
update for new `AlgorithmInterpretation`
2 parents 4baff6b + aed8d48 commit fdafa9c

File tree

4 files changed

+12
-8
lines changed

4 files changed

+12
-8
lines changed

Project.toml

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -34,9 +34,9 @@ NonlinearSolve = "3"
3434
OrdinaryDiffEq = "6.53"
3535
RecursiveArrayTools = "3"
3636
Reexport = "1.0"
37-
SciMLBase = "2"
37+
SciMLBase = "2.51"
3838
SteadyStateDiffEq = "2"
39-
StochasticDiffEq = "6.61"
39+
StochasticDiffEq = "6.69"
4040
Sundials = "4.19"
4141
julia = "1.9"
4242

src/DifferentialEquations.jl

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -21,6 +21,7 @@ using LinearAlgebra
2121
import DiffEqBase: solve
2222
@reexport using LinearSolve
2323
using NonlinearSolve
24+
import SciMLBase
2425

2526
include("default_solve.jl")
2627
include("default_arg_parsing.jl")

src/sde_default_alg.jl

Lines changed: 4 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -13,16 +13,17 @@ function default_algorithm(
1313
end
1414

1515
is_stiff = :stiff alg_hints
16-
is_stratonovich = :Stratonovich alg_hints
16+
is_stratonovich = :stratonovich alg_hints
1717
if is_stiff || prob.f.mass_matrix !== I
1818
alg = ImplicitRKMil(autodiff = false)
1919
end
2020

2121
if is_stratonovich
2222
if is_stiff || prob.f.mass_matrix !== I
23-
alg = ImplicitRKMil(autodiff = false, interpretation = :Stratonovich)
23+
alg = ImplicitRKMil(autodiff = false,
24+
interpretation = SciMLBase.AlgorithmInterpretation.Stratonovich)
2425
else
25-
alg = RKMil(interpretation = :Stratonovich)
26+
alg = RKMil(interpretation = SciMLBase.AlgorithmInterpretation.Stratonovich)
2627
end
2728
end
2829

test/default_sde_alg_test.jl

Lines changed: 5 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,5 @@
11
using DifferentialEquations, Test
2+
import SciMLBase
23

34
f_additive(u, p, t) = @. p[2] / sqrt(1 + t) - u / (2 * (1 + t))
45
σ_additive(u, p, t) = @. p[1] * p[2] / sqrt(1 + t)
@@ -17,8 +18,9 @@ sol = solve(prob, dt = 1 / 2^(3))
1718
sol = solve(prob, dt = 1 / 2^(3), alg_hints = [:additive])
1819
@test sol.alg isa SOSRA
1920

20-
sol = solve(prob, dt = 1 / 2^(3), alg_hints = [:Stratonovich])
21-
@test StochasticDiffEq.alg_interpretation(sol.alg) == :Stratonovich
21+
sol = solve(prob, dt = 1 / 2^(3), alg_hints = [:stratonovich])
22+
@test SciMLBase.alg_interpretation(sol.alg) ==
23+
SciMLBase.AlgorithmInterpretation.Stratonovich
2224
@test sol.alg isa RKMil
2325

2426
f = (du, u, p, t) -> du .= 1.01u
@@ -43,5 +45,5 @@ sol = solve(prob, dt = 1 / 2^(3), alg_hints = [:stiff])
4345
sol = solve(prob, dt = 1 / 2^(3), alg_hints = [:additive])
4446
@test sol.alg isa SOSRA
4547

46-
sol = solve(prob, dt = 1 / 2^(3), alg_hints = [:Stratonovich])
48+
sol = solve(prob, dt = 1 / 2^(3), alg_hints = [:stratonovich])
4749
@test sol.alg isa LambaEulerHeun

0 commit comments

Comments
 (0)