Skip to content
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
8 changes: 1 addition & 7 deletions lib/ODEProblemLibrary/Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -4,23 +4,17 @@ version = "0.1.12"

[deps]
DiffEqBase = "2b5f629d-d688-5b77-993f-72d75c75574e"
Latexify = "23fbe1c1-3f47-55db-b15f-69d7ec21a316"
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
Markdown = "d6f4376e-aef5-505a-96c1-9c027394607a"
ModelingToolkit = "961ee093-0014-501f-94e3-6117800e7a78"
Random = "9a3f8284-a2c9-5f02-9a11-845980a1fd5c"
RuntimeGeneratedFunctions = "7e49a35a-f44a-4d26-94aa-eba1b4ca6b47"

[compat]
Aqua = "0.5"
DiffEqBase = "6"
Latexify = "0.16"
ModelingToolkit = "10"
RuntimeGeneratedFunctions = "0.5"
julia = "1.10"

[extras]
Aqua = "4c88cf16-eb10-579e-8560-4a9242c79595"

[targets]
test = ["Aqua"]
test = ["Aqua"]
19 changes: 3 additions & 16 deletions lib/ODEProblemLibrary/src/ODEProblemLibrary.jl
Original file line number Diff line number Diff line change
@@ -1,17 +1,10 @@
module ODEProblemLibrary

using DiffEqBase
using Latexify
using ModelingToolkit
using ModelingToolkit: t_nounits as t, D_nounits as D

using LinearAlgebra
using Markdown
using Random

import RuntimeGeneratedFunctions
RuntimeGeneratedFunctions.init(@__MODULE__)

Random.seed!(100)

#ODE Example Problems
Expand All @@ -30,19 +23,13 @@ export prob_ode_linear, prob_ode_bigfloatlinear, prob_ode_2Dlinear,
prob_ode_chen, prob_ode_rossler, prob_ode_rabinovich_fabrikant, prob_ode_sprott,
prob_ode_hindmarsh_rose

# For MTKv9 compatibility
@static if !isdefined(ModelingToolkit, :mtkcompile)
function mtkcompile(args...; kwargs...)
structural_simplify(args...; kwargs...)
end
end

# Include all problem definitions
include("ode_linear_prob.jl")
include("nonlinchem.jl")
include("ode_simple_nonlinear_prob.jl")
include("brusselator_prob.jl")
include("pollution_prob.jl")
include("filament_prob.jl")
include("nonlinchem.jl")
include("strange_attractors.jl")

end # module
end # module
133 changes: 68 additions & 65 deletions lib/ODEProblemLibrary/src/ode_simple_nonlinear_prob.jl
Original file line number Diff line number Diff line change
Expand Up @@ -48,13 +48,15 @@ with initial condition ``v=w=1``
prob_ode_fitzhughnagumo = ODEProblem(fitz, [1.0; 1.0], (0.0, 1.0),
(0.7, 0.8, 1 / 12.5, 0.5))

#Van der Pol Equations
@parameters μ
@variables x(t) y(t)
## Van der Pol Equations

eqs = [D(y) ~ μ * ((1 - x^2) * y - x),
D(x) ~ y]
van = System(eqs, t; name = :van_der_pol) |> mtkcompile
function vanderpol(du, u, p, t)
x = u[1]
y = u[2]
μ = p[1]
du[1] = y
du[2] = μ * ((1 - x^2) * y - x)
end

"""
Van der Pol Equations
Expand All @@ -66,12 +68,11 @@ Van der Pol Equations
\\frac{dy}{dt} = μ((1-x^2)y -x)
```

with ``μ=1.0`` and ``u_0=[x => \\sqrt{3}, y => 0]``
with ``μ=1.0`` and ``u_0=[\\sqrt{3}, 0]`` (where ``u[1] = x``, ``u[2] = y``)

Non-stiff parameters.
"""
prob_ode_vanderpol = ODEProblem(van, [y => 0, x => sqrt(3), μ => 1.0], (0.0, 1.0),
jac = true, eval_module = @__MODULE__)
prob_ode_vanderpol = ODEProblem(vanderpol, [sqrt(3), 0.0], (0.0, 1.0), [1.0])

"""
Van der Pol Equations
Expand All @@ -83,21 +84,25 @@ Van der Pol Equations
\\frac{dy}{dt} = μ((1-x^2)y -x)
```

with ``μ=10^6`` and ``u_0=[x => \\sqrt{3}, y => 0]``
with ``μ=10^6`` and ``u_0=[\\sqrt{3}, 0]`` (where ``u[1] = x``, ``u[2] = y``)

Stiff parameters.
"""
prob_ode_vanderpol_stiff = ODEProblem(van, [y => 0, x => sqrt(3), μ => 1e6], (0.0, 1.0),
jac = true, eval_module = @__MODULE__)

# ROBER
@parameters k₁ k₂ k₃
@variables y₁(t) y₂(t) y₃(t)

eqs = [D(y₁) ~ -k₁ * y₁ + k₃ * y₂ * y₃,
D(y₂) ~ k₁ * y₁ - k₂ * y₂^2 - k₃ * y₂ * y₃,
D(y₃) ~ k₂ * y₂^2]
rober = System(eqs, t; name = :rober) |> mtkcompile
prob_ode_vanderpol_stiff = ODEProblem(vanderpol, [sqrt(3), 0.0], (0.0, 1.0), [1e6])

## ROBER

function rober(du, u, p, t)
y₁ = u[1]
y₂ = u[2]
y₃ = u[3]
k₁ = p[1]
k₂ = p[2]
k₃ = p[3]
du[1] = -k₁ * y₁ + k₃ * y₂ * y₃
du[2] = k₁ * y₁ - k₂ * y₂^2 - k₃ * y₂ * y₃
du[3] = k₂ * y₂^2
end

"""
The Robertson biochemical reactions: (Stiff)
Expand All @@ -118,9 +123,7 @@ Hairer Norsett Wanner Solving Ordinary Differential Equations I - Nonstiff Probl

Usually solved on ``[0,1e11]``
"""
prob_ode_rober = ODEProblem(
rober, [[y₁, y₂, y₃] .=> [1.0; 0.0; 0.0]; [k₁, k₂, k₃] .=> (0.04, 3e7, 1e4)],
(0.0, 1e11), jac = true, eval_module = @__MODULE__)
prob_ode_rober = ODEProblem(rober, [1.0, 0.0, 0.0], (0.0, 1e11), [0.04, 3e7, 1e4])

# Three Body
const threebody_μ = big(0.012277471);
Expand Down Expand Up @@ -174,15 +177,19 @@ prob_ode_threebody = ODEProblem(threebody,
[0.994, 0.0, 0.0, big(-2.00158510637908252240537862224)],
(big(0.0), big(17.0652165601579625588917206249)))

# Rigid Body Equations

@parameters I₁ I₂ I₃
@variables y₁(t) y₂(t) y₃(t)

eqs = [D(y₁) ~ I₁ * y₂ * y₃,
D(y₂) ~ I₂ * y₁ * y₃,
D(y₃) ~ I₃ * y₁ * y₂]
rigid = System(eqs, t; name = :rigid_body) |> mtkcompile
## Rigid Body Equations

function rigidbody(du, u, p, t)
y₁ = u[1]
y₂ = u[2]
y₃ = u[3]
I₁ = p[1]
I₂ = p[2]
I₃ = p[3]
du[1] = I₁ * y₂ * y₃
du[2] = I₂ * y₁ * y₃
du[3] = I₃ * y₁ * y₂
end

"""
Rigid Body Equations (Non-stiff)
Expand All @@ -207,10 +214,7 @@ or Hairer Norsett Wanner Solving Ordinary Differential Equations I - Nonstiff Pr

Usually solved from 0 to 20.
"""
prob_ode_rigidbody = ODEProblem(
rigid, [[y₁, y₂, y₃] .=> [1.0, 0.0, 0.9]; [I₁, I₂, I₃] .=> (-2.0, 1.25, -0.5)],
(0.0, 20.0),
jac = true, eval_module = @__MODULE__)
prob_ode_rigidbody = ODEProblem(rigidbody, [1.0, 0.0, 0.9], (0.0, 20.0), [-2.0, 1.25, -0.5])

# Pleiades Problem

Expand Down Expand Up @@ -356,28 +360,27 @@ mm_f = ODEFunction(mm_linear; analytic = (u0, p, t) -> exp(inv(MM_linear) * mm_A
mass_matrix = MM_linear)
prob_ode_mm_linear = ODEProblem(mm_f, rand(4), (0.0, 1.0))

@parameters p1 p2 p3 p4 p5 p6 p7 p8 p9 p10 p11 p12
@variables y1(t) y2(t) y3(t) y4(t) y5(t) y6(t) y7(t) y8(t)

eqs = [D(y1) ~ -p1 * y1 + p2 * y2 + p3 * y3 + p4,
D(y2) ~ p1 * y1 - p5 * y2,
D(y3) ~ -p6 * y3 + p2 * y4 + p7 * y5,
D(y4) ~ p3 * y2 + p1 * y3 - p8 * y4,
D(y5) ~ -p9 * y5 + p2 * y6 + p2 * y7,
D(y6) ~ -p10 * y6 * y8 + p11 * y4 + p1 * y5 -
p2 * y6 + p11 * y7,
D(y7) ~ p10 * y6 * y8 - p12 * y7,
D(y8) ~ -p10 * y6 * y8 + p12 * y7]
hires = System(eqs, t; name = :hires) |> mtkcompile
## Hires Problem

function hires(du, u, p, t)
y1, y2, y3, y4, y5, y6, y7, y8 = u
p1, p2, p3, p4, p5, p6, p7, p8, p9, p10, p11, p12 = p

du[1] = -p1 * y1 + p2 * y2 + p3 * y3 + p4
du[2] = p1 * y1 - p5 * y2
du[3] = -p6 * y3 + p2 * y4 + p7 * y5
du[4] = p3 * y2 + p1 * y3 - p8 * y4
du[5] = -p9 * y5 + p2 * y6 + p2 * y7
du[6] = -p10 * y6 * y8 + p11 * y4 + p1 * y5 - p2 * y6 + p11 * y7
du[7] = p10 * y6 * y8 - p12 * y7
du[8] = -p10 * y6 * y8 + p12 * y7
end

u0 = zeros(8)
u0[1] = 1
u0[8] = 0.0057

u0 = [y1, y2, y3, y4, y5, y6, y7, y8] .=> u0
p = [p1, p2, p3, p4, p5, p6, p7, p8, p9, p10, p11,
p12] .=> (1.71, 0.43, 8.32, 0.0007, 8.75,
10.03, 0.035, 1.12, 1.745, 280.0, 0.69, 1.81)
p = (1.71, 0.43, 8.32, 0.0007, 8.75, 10.03, 0.035, 1.12, 1.745, 280.0, 0.69, 1.81)

"""
Hires Problem (Stiff)
Expand All @@ -401,16 +404,18 @@ where ``f`` is defined by
Reference: [demohires.pdf](http://www.radford.edu/~thompson/vodef90web/problems/demosnodislin/Demos_Pitagora/DemoHires/demohires.pdf)
Notebook: [Hires.ipynb](http://nbviewer.jupyter.org/github/JuliaDiffEq/DiffEqBenchmarks.jl/blob/master/StiffODE/Hires.ipynb)
"""
prob_ode_hires = ODEProblem(
hires, [u0; p], (0.0, 321.8122), jac = true, eval_module = @__MODULE__)
prob_ode_hires = ODEProblem(hires, u0, (0.0, 321.8122), p)

@parameters p1 p2 p3
@variables y1(t) y2(t) y3(t)
## Orego Problem

eqs = [D(y1) ~ p1 * (y2 + y1 * (1 - p2 * y1 - y2)),
D(y2) ~ (y3 - (1 + y1) * y2) / p1,
D(y3) ~ p3 * (y1 - y3)]
orego = System(eqs, t; name = :orego) |> mtkcompile
function orego(du, u, p, t)
y1, y2, y3 = u
p1, p2, p3 = p

du[1] = p1 * (y2 + y1 * (1 - p2 * y1 - y2))
du[2] = (y3 - (1 + y1) * y2) / p1
du[3] = p3 * (y1 - y3)
end

"""
Orego Problem (Stiff)
Expand All @@ -430,6 +435,4 @@ where ``s=77.27``, ``w=0.161`` and ``q=8.375⋅10^{-6}``.
Reference: [demoorego.pdf](http://www.radford.edu/~thompson/vodef90web/problems/demosnodislin/Demos_Pitagora/DemoOrego/demoorego.pdf)
Notebook: [Orego.ipynb](http://nbviewer.jupyter.org/github/JuliaDiffEq/DiffEqBenchmarks.jl/blob/master/StiffODE/Orego.ipynb)
"""
prob_ode_orego = ODEProblem(
orego, [[y1, y2, y3] .=> [1.0, 2.0, 3.0]; [p1, p2, p3] .=> [77.27, 8.375e-6, 0.161]],
(0.0, 30.0), jac = true, eval_module = @__MODULE__)
prob_ode_orego = ODEProblem(orego, [1.0, 2.0, 3.0], (0.0, 30.0), [77.27, 8.375e-6, 0.161])
Loading
Loading