using Turing
using DifferentialEquations
# Load StatsPlots for visualizations and diagnostics.
using StatsPlots
using LinearAlgebra
# Set a seed for reproducibility.
using Random
Random.seed!(14);
Bayesian Estimation of Differential Equations
Most of the scientific community deals with the basic problem of trying to mathematically model the reality around them and this often involves dynamical systems. The general trend to model these complex dynamical systems is through the use of differential equations. Differential equation models often have non-measurable parameters. The popular “forward-problem” of simulation consists of solving the differential equations for a given set of parameters, the “inverse problem” to simulation, known as parameter estimation, is the process of utilizing data to determine these model parameters. Bayesian inference provides a robust approach to parameter estimation with quantified uncertainty.
The Lotka-Volterra Model
The Lotka–Volterra equations, also known as the predator–prey equations, are a pair of first-order nonlinear differential equations. These differential equations are frequently used to describe the dynamics of biological systems in which two species interact, one as a predator and the other as prey. The populations change through time according to the pair of equations
\[ \begin{aligned} \frac{\mathrm{d}x}{\mathrm{d}t} &= (\alpha - \beta y(t))x(t), \\ \frac{\mathrm{d}y}{\mathrm{d}t} &= (\delta x(t) - \gamma)y(t) \end{aligned} \]
where \(x(t)\) and \(y(t)\) denote the populations of prey and predator at time \(t\), respectively, and \(\alpha, \beta, \gamma, \delta\) are positive parameters.
We implement the Lotka-Volterra model and simulate it with parameters \(\alpha = 1.5\), \(\beta = 1\), \(\gamma = 3\), and \(\delta = 1\) and initial conditions \(x(0) = y(0) = 1\).
# Define Lotka-Volterra model.
function lotka_volterra(du, u, p, t)
# Model parameters.
= p
α, β, γ, δ # Current state.
= u
x, y
# Evaluate differential equations.
1] = (α - β * y) * x # prey
du[2] = (δ * x - γ) * y # predator
du[
return nothing
end
# Define initial-value problem.
= [1.0, 1.0]
u0 = [1.5, 1.0, 3.0, 1.0]
p = (0.0, 10.0)
tspan = ODEProblem(lotka_volterra, u0, tspan, p)
prob
# Plot simulation.
plot(solve(prob, Tsit5()))
We generate noisy observations to use for the parameter estimation tasks in this tutorial. With the saveat
argument we specify that the solution is stored only at 0.1
time units. To make the example more realistic we add random normally distributed noise to the simulation.
= solve(prob, Tsit5(); saveat=0.1)
sol = Array(sol) + 0.8 * randn(size(Array(sol)))
odedata
# Plot simulation and noisy observations.
plot(sol; alpha=0.3)
scatter!(sol.t, odedata'; color=[1 2], label="")
Alternatively, we can use real-world data from Hudson’s Bay Company records (an Stan implementation with slightly different priors can be found here: https://mc-stan.org/users/documentation/case-studies/lotka-volterra-predator-prey.html).
Direct Handling of Bayesian Estimation with Turing
Previously, functions in Turing and DifferentialEquations were not inter-composable, so Bayesian inference of differential equations needed to be handled by another package called DiffEqBayes.jl (note that DiffEqBayes works also with CmdStan.jl, Turing.jl, DynamicHMC.jl and ApproxBayes.jl - see the DiffEqBayes docs for more info).
Nowadays, however, Turing and DifferentialEquations are completely composable and we can just simulate differential equations inside a Turing @model
. Therefore, we write the Lotka-Volterra parameter estimation problem using the Turing @model
macro as below:
@model function fitlv(data, prob)
# Prior distributions.
~ InverseGamma(2, 3)
σ ~ truncated(Normal(1.5, 0.5); lower=0.5, upper=2.5)
α ~ truncated(Normal(1.2, 0.5); lower=0, upper=2)
β ~ truncated(Normal(3.0, 0.5); lower=1, upper=4)
γ ~ truncated(Normal(1.0, 0.5); lower=0, upper=2)
δ
# Simulate Lotka-Volterra model.
= [α, β, γ, δ]
p = solve(prob, Tsit5(); p=p, saveat=0.1)
predicted
# Observations.
for i in 1:length(predicted)
:, i] ~ MvNormal(predicted[i], σ^2 * I)
data[end
return nothing
end
= fitlv(odedata, prob)
model
# Sample 3 independent chains with forward-mode automatic differentiation (the default).
= sample(model, NUTS(), MCMCSerial(), 1000, 3; progress=false) chain
┌ Info: Found initial step size
└ ϵ = 0.2
┌ Info: Found initial step size
└ ϵ = 0.00625
┌ Info: Found initial step size
└ ϵ = 0.025
Chains MCMC chain (1000×17×3 Array{Float64, 3}):
Iterations = 501:1:1500
Number of chains = 3
Samples per chain = 1000
Wall duration = 50.3 seconds
Compute duration = 47.5 seconds
parameters = σ, α, β, γ, δ
internals = lp, n_steps, is_accept, acceptance_rate, log_density, hamiltonian_energy, hamiltonian_energy_error, max_hamiltonian_energy_error, tree_depth, numerical_error, step_size, nom_step_size
Summary Statistics
parameters mean std mcse ess_bulk ess_tail rhat e ⋯
Symbol Float64 Float64 Float64 Float64 Float64 Float64 ⋯
σ 1.1932 0.5938 0.2410 9.5537 104.5556 1.6631 ⋯
α 1.3888 0.1382 0.0496 9.7176 85.6723 1.6497 ⋯
β 0.9570 0.0970 0.0274 15.7758 73.3357 1.2876 ⋯
γ 2.4489 0.9443 0.3817 9.5715 100.7936 1.6664 ⋯
δ 0.8819 0.2312 0.0913 9.5719 72.8664 1.6649 ⋯
1 column omitted
Quantiles
parameters 2.5% 25.0% 50.0% 75.0% 97.5%
Symbol Float64 Float64 Float64 Float64 Float64
σ 0.7101 0.7622 0.8020 1.9558 2.1757
α 1.0517 1.2899 1.4425 1.4886 1.5567
β 0.7299 0.9149 0.9826 1.0210 1.0911
γ 1.0090 1.1764 3.0053 3.1488 3.3965
δ 0.4824 0.6027 1.0017 1.0532 1.1454
The estimated parameters are close to the parameter values the observations were generated with. We can also check visually that the chains have converged.
plot(chain)
Data retrodiction
In Bayesian analysis it is often useful to retrodict the data, i.e. generate simulated data using samples from the posterior distribution, and compare to the original data (see for instance section 3.3.2 - model checking of McElreath’s book “Statistical Rethinking”). Here, we solve the ODE for 300 randomly picked posterior samples in the chain
. We plot the ensemble of solutions to check if the solution resembles the data. The 300 retrodicted time courses from the posterior are plotted in gray, the noisy observations are shown as blue and red dots, and the green and purple lines are the ODE solution that was used to generate the data.
plot(; legend=false)
= sample(chain[[:α, :β, :γ, :δ]], 300; replace=false)
posterior_samples for p in eachrow(Array(posterior_samples))
= solve(prob, Tsit5(); p=p, saveat=0.1)
sol_p plot!(sol_p; alpha=0.1, color="#BBBBBB")
end
# Plot simulation and noisy observations.
plot!(sol; color=[1 2], linewidth=1)
scatter!(sol.t, odedata'; color=[1 2])
We can see that, even though we added quite a bit of noise to the data the posterior distribution reproduces quite accurately the “true” ODE solution.
Lotka-Volterra model without data of prey
One can also perform parameter inference for a Lotka-Volterra model with incomplete data. For instance, let us suppose we have only observations of the predators but not of the prey. I.e., we fit the model only to the \(y\) variable of the system without providing any data for \(x\):
@model function fitlv2(data::AbstractVector, prob)
# Prior distributions.
~ InverseGamma(2, 3)
σ ~ truncated(Normal(1.5, 0.5); lower=0.5, upper=2.5)
α ~ truncated(Normal(1.2, 0.5); lower=0, upper=2)
β ~ truncated(Normal(3.0, 0.5); lower=1, upper=4)
γ ~ truncated(Normal(1.0, 0.5); lower=0, upper=2)
δ
# Simulate Lotka-Volterra model but save only the second state of the system (predators).
= [α, β, γ, δ]
p = solve(prob, Tsit5(); p=p, saveat=0.1, save_idxs=2)
predicted
# Observations of the predators.
~ MvNormal(predicted.u, σ^2 * I)
data
return nothing
end
= fitlv2(odedata[2, :], prob)
model2
# Sample 3 independent chains.
= sample(model2, NUTS(0.45), MCMCSerial(), 5000, 3; progress=false) chain2
┌ Info: Found initial step size
└ ϵ = 0.025
┌ Info: Found initial step size
└ ϵ = 0.025
┌ Info: Found initial step size
└ ϵ = 0.025
Chains MCMC chain (5000×17×3 Array{Float64, 3}):
Iterations = 1001:1:6000
Number of chains = 3
Samples per chain = 5000
Wall duration = 33.74 seconds
Compute duration = 33.29 seconds
parameters = σ, α, β, γ, δ
internals = lp, n_steps, is_accept, acceptance_rate, log_density, hamiltonian_energy, hamiltonian_energy_error, max_hamiltonian_energy_error, tree_depth, numerical_error, step_size, nom_step_size
Summary Statistics
parameters mean std mcse ess_bulk ess_tail rhat e ⋯
Symbol Float64 Float64 Float64 Float64 Float64 Float64 ⋯
σ 0.7918 0.0677 0.0087 56.4333 41.4030 1.1432 ⋯
α 1.6164 0.1967 0.0244 60.0994 50.2946 1.1635 ⋯
β 1.1273 0.1514 0.0184 62.8473 69.8384 1.1553 ⋯
γ 2.9995 0.3043 0.0399 62.5634 46.0728 1.1453 ⋯
δ 0.8701 0.2410 0.0323 60.5329 57.0321 1.1644 ⋯
1 column omitted
Quantiles
parameters 2.5% 25.0% 50.0% 75.0% 97.5%
Symbol Float64 Float64 Float64 Float64 Float64
σ 0.6707 0.7406 0.7889 0.8377 0.9311
α 1.2737 1.4785 1.6022 1.7491 2.0408
β 0.8725 1.0194 1.1083 1.2290 1.4556
γ 2.4474 2.7842 2.9893 3.1746 3.6879
δ 0.4554 0.6843 0.8521 1.0238 1.3894
Again we inspect the trajectories of 300 randomly selected posterior samples.
plot(; legend=false)
= sample(chain2[[:α, :β, :γ, :δ]], 300; replace=false)
posterior_samples for p in eachrow(Array(posterior_samples))
= solve(prob, Tsit5(); p=p, saveat=0.1)
sol_p plot!(sol_p; alpha=0.1, color="#BBBBBB")
end
# Plot simulation and noisy observations.
plot!(sol; color=[1 2], linewidth=1)
scatter!(sol.t, odedata'; color=[1 2])
Note that here the observations of the prey (blue dots) were not used in the parameter estimation! Yet, the model can predict the values of \(x\) relatively accurately, albeit with a wider distribution of solutions, reflecting the greater uncertainty in the prediction of the \(x\) values.
Inference of Delay Differential Equations
Here we show an example of inference with another type of differential equation: a Delay Differential Equation (DDE). DDEs are differential equations where derivatives are function of values at an earlier point in time. This is useful to model a delayed effect, like incubation time of a virus for instance.
Here is a delayed version of the Lokta-Voltera system:
\[ \begin{aligned} \frac{\mathrm{d}x}{\mathrm{d}t} &= \alpha x(t-\tau) - \beta y(t) x(t),\\ \frac{\mathrm{d}y}{\mathrm{d}t} &= - \gamma y(t) + \delta x(t) y(t), \end{aligned} \]
where \(\tau\) is a (positive) delay and \(x(t-\tau)\) is the variable \(x\) at an earlier time point \(t - \tau\).
The initial-value problem of the delayed system can be implemented as a DDEProblem
. As described in the DDE example, here the function h
is the history function that can be used to obtain a state at an earlier time point. Again we use parameters \(\alpha = 1.5\), \(\beta = 1\), \(\gamma = 3\), and \(\delta = 1\) and initial conditions \(x(0) = y(0) = 1\). Moreover, we assume \(x(t) = 1\) for \(t < 0\).
function delay_lotka_volterra(du, u, h, p, t)
# Model parameters.
= p
α, β, γ, δ
# Current state.
= u
x, y # Evaluate differential equations
1] = α * h(p, t - 1; idxs=1) - β * x * y
du[2] = -γ * y + δ * x * y
du[
return nothing
end
# Define initial-value problem.
= (1.5, 1.0, 3.0, 1.0)
p = [1.0; 1.0]
u0 = (0.0, 10.0)
tspan h(p, t; idxs::Int) = 1.0
= DDEProblem(delay_lotka_volterra, u0, h, tspan, p); prob_dde
We generate observations by adding normally distributed noise to the results of our simulations.
= solve(prob_dde; saveat=0.1)
sol_dde = Array(sol_dde) + 0.5 * randn(size(sol_dde))
ddedata
# Plot simulation and noisy observations.
plot(sol_dde)
scatter!(sol_dde.t, ddedata'; color=[1 2], label="")
Now we define the Turing model for the Lotka-Volterra model with delay and sample 3 independent chains.
@model function fitlv_dde(data, prob)
# Prior distributions.
~ InverseGamma(2, 3)
σ ~ truncated(Normal(1.5, 0.5); lower=0.5, upper=2.5)
α ~ truncated(Normal(1.2, 0.5); lower=0, upper=2)
β ~ truncated(Normal(3.0, 0.5); lower=1, upper=4)
γ ~ truncated(Normal(1.0, 0.5); lower=0, upper=2)
δ
# Simulate Lotka-Volterra model.
= [α, β, γ, δ]
p = solve(prob, MethodOfSteps(Tsit5()); p=p, saveat=0.1)
predicted
# Observations.
for i in 1:length(predicted)
:, i] ~ MvNormal(predicted[i], σ^2 * I)
data[end
end
= fitlv_dde(ddedata, prob_dde)
model_dde
# Sample 3 independent chains.
= sample(model_dde, NUTS(), MCMCSerial(), 300, 3; progress=false) chain_dde
┌ Info: Found initial step size
└ ϵ = 0.0125
┌ Info: Found initial step size
└ ϵ = 0.025
┌ Info: Found initial step size
└ ϵ = 0.2
Chains MCMC chain (300×17×3 Array{Float64, 3}):
Iterations = 151:1:450
Number of chains = 3
Samples per chain = 300
Wall duration = 15.08 seconds
Compute duration = 14.76 seconds
parameters = σ, α, β, γ, δ
internals = lp, n_steps, is_accept, acceptance_rate, log_density, hamiltonian_energy, hamiltonian_energy_error, max_hamiltonian_energy_error, tree_depth, numerical_error, step_size, nom_step_size
Summary Statistics
parameters mean std mcse ess_bulk ess_tail rhat e ⋯
Symbol Float64 Float64 Float64 Float64 Float64 Float64 ⋯
σ 0.5282 0.0250 0.0009 789.3204 564.1834 1.0037 ⋯
α 1.4958 0.0596 0.0037 261.2371 408.4717 1.0030 ⋯
β 1.0157 0.0518 0.0027 381.3145 476.9419 1.0040 ⋯
γ 3.0202 0.1447 0.0091 253.2131 368.8544 1.0061 ⋯
δ 1.0009 0.0500 0.0032 245.8068 314.7055 1.0046 ⋯
1 column omitted
Quantiles
parameters 2.5% 25.0% 50.0% 75.0% 97.5%
Symbol Float64 Float64 Float64 Float64 Float64
σ 0.4815 0.5110 0.5277 0.5440 0.5808
α 1.3860 1.4570 1.4915 1.5340 1.6147
β 0.9224 0.9823 1.0097 1.0468 1.1223
γ 2.7454 2.9211 3.0241 3.1170 3.3011
δ 0.9058 0.9651 1.0022 1.0342 1.1002
plot(chain_dde)
Finally, plot trajectories of 300 randomly selected samples from the posterior. Again, the dots indicate our observations, the colored lines are the “true” simulations without noise, and the gray lines are trajectories from the posterior samples.
plot(; legend=false)
= sample(chain_dde[[:α, :β, :γ, :δ]], 300; replace=false)
posterior_samples for p in eachrow(Array(posterior_samples))
= solve(prob_dde, MethodOfSteps(Tsit5()); p=p, saveat=0.1)
sol_p plot!(sol_p; alpha=0.1, color="#BBBBBB")
end
# Plot simulation and noisy observations.
plot!(sol_dde; color=[1 2], linewidth=1)
scatter!(sol_dde.t, ddedata'; color=[1 2])
The fit is pretty good even though the data was quite noisy to start.
Scaling to Large Models: Adjoint Sensitivities
DifferentialEquations.jl’s efficiency for large stiff models has been shown in multiple benchmarks. To learn more about how to optimize solving performance for stiff problems you can take a look at the docs.
Sensitivity analysis, or automatic differentiation (AD) of the solver, is provided by the DiffEq suite. The model sensitivities are the derivatives of the solution with respect to the parameters. Specifically, the local sensitivity of the solution to a parameter is defined by how much the solution would change by changes in the parameter. Sensitivity analysis provides a cheap way to calculate the gradient of the solution which can be used in parameter estimation and other optimization tasks.
The AD ecosystem in Julia allows you to switch between forward mode, reverse mode, source to source and other choices of AD and have it work with any Julia code. For a user to make use of this within SciML, high level interactions in solve
automatically plug into those AD systems to allow for choosing advanced sensitivity analysis (derivative calculation) methods.
More theoretical details on these methods can be found at: https://docs.sciml.ai/latest/extras/sensitivity_math/.
While these sensitivity analysis methods may seem complicated, using them is dead simple. Here is a version of the Lotka-Volterra model using adjoint sensitivities.
All we have to do is switch the AD backend to one of the adjoint-compatible backends (ReverseDiff or Zygote)! Notice that on this model adjoints are slower. This is because adjoints have a higher overhead on small parameter models and therefore we suggest using these methods only for models with around 100 parameters or more. For more details, see https://arxiv.org/abs/1812.01892.
using Zygote, SciMLSensitivity
# Sample a single chain with 1000 samples using Zygote.
sample(model, NUTS(;adtype=AutoZygote()), 1000; progress=false)
┌ Info: Found initial step size
└ ϵ = 0.025
Chains MCMC chain (1000×17×1 Array{Float64, 3}):
Iterations = 501:1:1500
Number of chains = 1
Samples per chain = 1000
Wall duration = 872.07 seconds
Compute duration = 872.07 seconds
parameters = σ, α, β, γ, δ
internals = lp, n_steps, is_accept, acceptance_rate, log_density, hamiltonian_energy, hamiltonian_energy_error, max_hamiltonian_energy_error, tree_depth, numerical_error, step_size, nom_step_size
Summary Statistics
parameters mean std mcse ess_bulk ess_tail rhat e ⋯
Symbol Float64 Float64 Float64 Float64 Float64 Float64 ⋯
σ 0.7771 0.0412 0.0018 500.6684 486.0897 1.0011 ⋯
α 1.4697 0.0462 0.0029 251.3591 288.1128 1.0108 ⋯
β 1.0010 0.0437 0.0022 382.1511 402.0423 1.0070 ⋯
γ 3.1182 0.1435 0.0090 258.4142 310.0951 1.0104 ⋯
δ 1.0429 0.0525 0.0033 253.0956 285.3338 1.0089 ⋯
1 column omitted
Quantiles
parameters 2.5% 25.0% 50.0% 75.0% 97.5%
Symbol Float64 Float64 Float64 Float64 Float64
σ 0.7026 0.7489 0.7740 0.8031 0.8623
α 1.3854 1.4377 1.4693 1.4984 1.5699
β 0.9236 0.9693 1.0010 1.0300 1.0902
γ 2.8291 3.0239 3.1215 3.2071 3.3981
δ 0.9394 1.0081 1.0424 1.0779 1.1510
If desired, we can control the sensitivity analysis method that is used by providing the sensealg
keyword argument to solve
. Here we will not choose a sensealg
and let it use the default choice:
@model function fitlv_sensealg(data, prob)
# Prior distributions.
~ InverseGamma(2, 3)
σ ~ truncated(Normal(1.5, 0.5); lower=0.5, upper=2.5)
α ~ truncated(Normal(1.2, 0.5); lower=0, upper=2)
β ~ truncated(Normal(3.0, 0.5); lower=1, upper=4)
γ ~ truncated(Normal(1.0, 0.5); lower=0, upper=2)
δ
# Simulate Lotka-Volterra model and use a specific algorithm for computing sensitivities.
= [α, β, γ, δ]
p = solve(prob; p=p, saveat=0.1)
predicted
# Observations.
for i in 1:length(predicted)
:, i] ~ MvNormal(predicted[i], σ^2 * I)
data[end
return nothing
end;
= fitlv_sensealg(odedata, prob)
model_sensealg
# Sample a single chain with 1000 samples using Zygote.
sample(model_sensealg, NUTS(;adtype=AutoZygote()), 1000; progress=false)
┌ Info: Found initial step size
└ ϵ = 0.00625
Chains MCMC chain (1000×17×1 Array{Float64, 3}):
Iterations = 501:1:1500
Number of chains = 1
Samples per chain = 1000
Wall duration = 442.08 seconds
Compute duration = 442.08 seconds
parameters = σ, α, β, γ, δ
internals = lp, n_steps, is_accept, acceptance_rate, log_density, hamiltonian_energy, hamiltonian_energy_error, max_hamiltonian_energy_error, tree_depth, numerical_error, step_size, nom_step_size
Summary Statistics
parameters mean std mcse ess_bulk ess_tail rhat e ⋯
Symbol Float64 Float64 Float64 Float64 Float64 Float64 ⋯
σ 2.1328 0.1063 0.0041 679.6424 564.1412 1.0003 ⋯
α 2.1556 0.1232 0.0085 204.1724 389.8030 1.0030 ⋯
β 1.8679 0.1043 0.0044 382.8599 361.4067 1.0001 ⋯
γ 3.5353 0.2447 0.0161 202.0714 121.5562 1.0038 ⋯
δ 1.7290 0.1408 0.0090 235.0671 160.0460 1.0072 ⋯
1 column omitted
Quantiles
parameters 2.5% 25.0% 50.0% 75.0% 97.5%
Symbol Float64 Float64 Float64 Float64 Float64
σ 1.9375 2.0579 2.1318 2.2024 2.3524
α 1.9521 2.0627 2.1464 2.2362 2.4000
β 1.6094 1.8127 1.8910 1.9498 1.9948
γ 3.0583 3.3534 3.5432 3.7145 3.9790
δ 1.4771 1.6185 1.7347 1.8333 1.9866
For more examples of adjoint usage on large parameter models, consult the DiffEqFlux documentation.