The ode_def macro

ParameterizedFunctions.@ode_defMacro
@ode_def name begin
    differential equation
end parameters :: ODEFunction

Definition of the Domain-Specific Language (DSL)

A helper macro is provided to make it easier to define a ParameterizedFunction, and it will symbolically compute a bunch of extra functions to make the differential equation solvers run faster. For example, to define the previous LotkaVolterra, you can use the following command:

f = @ode_def LotkaVolterra begin
    dx = a * x - b * x * y
    dy = -c * y + d * x * y
end a b c d

or you can define it anonymously:

f = @ode_def begin
    dx = a * x - b * x * y
    dy = -c * y + d * x * y
end a b c d

@ode_def uses ModelingToolkit.jl internally and returns an ODEFunction with the extra definitions (Jacobian, parameter Jacobian, etc.) defined through the MTK symbolic tools.

source
ParameterizedFunctions.@ode_def_bareMacro
@ode_def_bare name begin
    differential equation
end parameters :: ODEFunction

Like @ode_def but the opts options are set so that no symbolic functions are generated. See the @ode_def docstring for more details.

source
ParameterizedFunctions.@ode_def_allMacro
@ode_def_all name begin
    differential equation
end parameters :: ODEFunction

Like @ode_def but the opts options are set so that all possible symbolic functions are generated. See the @ode_def docstring for more details.

source

Internal API

ParameterizedFunctions.ode_def_optsFunction
ode_def_opts(name::Symbol, opts::Dict{Symbol, Bool}, curmod, ex::Expr, params...;
    depvar = :t)

The core internal. Users should only interact with this through the @ode_def_* macros.

Options are self-explanatory by name mapping to ODEFunction:

  • build_tgrad
  • build_jac
  • build_expjac
  • build_invjac
  • build_invW
  • buildinvWt
  • build_hes
  • build_invhes
  • build_dpfuncs

depvar sets the symbol for the dependent variable.

Example:

opts = Dict{Symbol, Bool}(:build_tgrad => true,
    :build_jac => true,
    :build_expjac => false,
    :build_invjac => true,
    :build_invW => true,
    :build_invW_t => true,
    :build_hes => false,
    :build_invhes => false,
    :build_dpfuncs => true)
source
ParameterizedFunctions.ParameterizedFunctionsModule

ParameterizedFunctions.jl

Join the chat at https://julialang.zulipchat.com #sciml-bridgedGlobal Docs

codecovBuild Status

ColPrac: Contributor's Guide on Collaborative Practices for Community PackagesSciML Code Style

ParameterizedFunctions.jl is a component of the SciML ecosystem which allows for easily defining parameterized ODE models in a simple syntax.

Tutorials and Documentation

For information on using the package, see the stable documentation. Use the in-development documentation for the version of the documentation, which contains the unreleased features.

Example

The following are valid ODE definitions.

using DifferentialEquations, ParameterizedFunctions

# Non-Stiff ODE

lotka_volterra = @ode_def begin
    d🐁 = α * 🐁 - β * 🐁 * 🐈
    d🐈 = -γ * 🐈 + δ * 🐁 * 🐈
end α β γ δ

p = [1.5, 1.0, 3.0, 1.0];
u0 = [1.0; 1.0];
prob = ODEProblem(lotka_volterra, u0, (0.0, 10.0), p)
sol = solve(prob, Tsit5(), reltol = 1e-6, abstol = 1e-6)

# Stiff ODE

rober = @ode_def begin
    dy₁ = -k₁ * y₁ + k₃ * y₂ * y₃
    dy₂ = k₁ * y₁ - k₂ * y₂^2 - k₃ * y₂ * y₃
    dy₃ = k₂ * y₂^2
end k₁ k₂ k₃

prob = ODEProblem(rober, [1.0, 0.0, 0.0], (0.0, 1e5), [0.04, 3e7, 1e4])
sol = solve(prob)
source