Parallel Ensemble Simulations
Performing Monte Carlo simulations, solving with a predetermined set of initial conditions, and GPU-parallelizing a parameter search all fall under the ensemble simulation interface. This interface allows one to declare a template DEProblem to parallelize, tweak the template in trajectories
for many trajectories, solve each in parallel batches, reduce the solutions down to specific answers, and compute summary statistics on the results.
Performing an Ensemble Simulation
Building a Problem
To perform a simulation on an ensemble of trajectories, define a EnsembleProblem
. The constructor is:
EnsembleProblem(prob::DEProblem;
output_func = (sol, i) -> (sol, false),
prob_func = (prob, i, repeat) -> (prob),
reduction = (u, data, I) -> (append!(u, data), false),
u_init = [], safetycopy = prob_func !== DEFAULT_PROB_FUNC)
-
output_func
: The function determines what is saved from the solution to the output array. Defaults to saving the solution itself. The output is(out,rerun)
whereout
is the output andrerun
is a boolean which designates whether to rerun. -
prob_func
: The function by which the problem is to be modified.prob
is the problem,i
is the unique id1:trajectories
for the problem, andrepeat
is the iteration of the repeat. At first, it is1
, but ifrerun
was true this will be2
,3
, etc. counting the number of times problemi
has been repeated. -
reduction
: This function determines how to reduce the data in each batch. Defaults to appending thedata
intou
, initialised viau_data
, from the batches.I
is a range of indices giving the trajectories corresponding to the batches. The second part of the output determines whether the simulation has converged. Iftrue
, the simulation will exit early. By default, this is alwaysfalse
. -
u_init
: The initial form of the object that gets updated in-place inside thereduction
function. -
safetycopy
: Determines whether a safetydeepcopy
is called on theprob
before theprob_func
. By default, this is true for any user-givenprob_func
, as without this, modifying the arguments of something in theprob_func
, such as parameters or caches stored within the user function, are not necessarily thread-safe. If you know that your function is thread-safe, then setting this tofalse
can improve performance when used with threads. For nested problems, e.g., SDE problems with custom noise processes,deepcopy
might be insufficient. In such cases, use a customprob_func
.
One can specify a function prob_func
which changes the problem. For example:
function prob_func(prob, i, repeat)
@. prob.u0 = randn() * prob.u0
prob
end
modifies the initial condition for all of the problems by a standard normal random number (a different random number per simulation). Notice that since problem types are immutable, it uses .=
. Otherwise, one can just create a new problem type:
function prob_func(prob, i, repeat)
@. prob.u0 = u0_arr[i]
prob
end
If your function is a ParameterizedFunction
, you can do similar modifications to prob.f
to perform a parameter search. The output_func
is a reduction function. Its arguments are the generated solution and the unique index for the run. For example, if we wish to only save the 2nd coordinate at the end of each solution, we can do:
output_func(sol, i) = (sol[end, 2], false)
Thus, the ensemble simulation would return as its data an array which is the end value of the 2nd dependent variable for each of the runs.
Solving the Problem
sim = solve(prob, alg, ensemblealg, kwargs...)
The keyword arguments take in the arguments for the common solver interface and will pass them to the differential equation solver. The ensemblealg
is optional, and will default to EnsembleThreads()
. The special keyword arguments to note are:
-
trajectories
: The number of simulations to run. This argument is required. -
batch_size
: The size of the batches on which the reductions are applies. Defaults totrajectories
. -
pmap_batch_size
: The size of thepmap
batches. Default isbatch_size÷100 > 0 ? batch_size÷100 : 1
EnsembleAlgorithms
The choice of ensemble algorithm allows for control over how the multiple trajectories are handled. Currently, the ensemble algorithm types are:
-
EnsembleSerial()
- No parallelism -
EnsembleThreads()
- The default. This uses multithreading. It’s local (single computer, shared memory) parallelism only. Fastest when the trajectories are quick. -
EnsembleDistributed()
- Usespmap
internally. It will use as many processors as you have Julia processes. To add more processes, useaddprocs(n)
. See Julia’s documentation for more details. Recommended for the case when each trajectory calculation isn’t “too quick” (at least about a millisecond each?). -
EnsembleSplitThreads()
- This uses threading on each process, splitting the problem intonprocs()
even parts. This is for solving many quick trajectories on a multi-node machine. It’s recommended you have one process on each node. -
EnsembleGPUArray()
- Requires installing andusing DiffEqGPU
. This uses a GPU for computing the ensemble with hyperparallelism. It will automatically recompile your Julia functions to the GPU. A standard GPU sees a 5x performance increase over a 16 core Xeon CPU. However, there are limitations on what functions can auto-compile in this fashion, please see the DiffEqGPU README for more details
For example, EnsembleThreads()
is invoked by:
solve(ensembleprob, alg, EnsembleThreads(); trajectories = 1000)
Plot Recipe
There is a plot recipe for a AbstractEnsembleSimulation
which composes all of the plot recipes for the component solutions. The keyword arguments are passed along. A useful argument to use is linealpha
which will change the transparency of the plots. An additional argument is idxs
which allows you to choose which components of the solution to plot. For example, if the differential equation is a vector of 9 values, idxs=1:2:9
will plot only the solutions of the odd components. Another additional argument is zcolors
(an alias of marker_z
) which allows you to pass a zcolor
for each series. For details about zcolor
see the Series documentation for Plots.jl.
Analyzing an Ensemble Experiment
Analysis tools are included for generating summary statistics and summary plots for a EnsembleSimulation
.
To use this functionality, import the analysis module via:
using DifferentialEquations.EnsembleAnalysis
(or more directly SciMLBase.EnsembleAnalysis
).
Time steps vs time points
For the summary statistics, there are two types. You can either summarize by time steps or by time points. Summarizing by time steps assumes that the time steps are all the same time point, i.e. the integrator used a fixed dt
or the values were saved using saveat
. Summarizing by time points requires interpolating the solution.
Analysis at a time step or time point
get_timestep(sim, i) # Returns an iterator of each simulation at time step i
get_timepoint(sim, t) # Returns an iterator of each simulation at time point t
componentwise_vectors_timestep(sim, i) # Returns a vector of each simulation at time step i
componentwise_vectors_timepoint(sim, t) # Returns a vector of each simulation at time point t
Summary Statistics
Single Time Statistics
The available functions for time steps are:
timestep_mean(sim, i) # Computes the mean of each component at time step i
timestep_median(sim, i) # Computes the median of each component at time step i
timestep_quantile(sim, q, i) # Computes the quantile q of each component at time step i
timestep_meanvar(sim, i) # Computes the mean and variance of each component at time step i
timestep_meancov(sim, i, j) # Computes the mean at i and j, and the covariance, for each component
timestep_meancor(sim, i, j) # Computes the mean at i and j, and the correlation, for each component
timestep_weighted_meancov(sim, W, i, j) # Computes the mean at i and j, and the weighted covariance W, for each component
The available functions for time points are:
timepoint_mean(sim, t) # Computes the mean of each component at time t
timepoint_median(sim, t) # Computes the median of each component at time t
timepoint_quantile(sim, q, t) # Computes the quantile q of each component at time t
timepoint_meanvar(sim, t) # Computes the mean and variance of each component at time t
timepoint_meancov(sim, t1, t2) # Computes the mean at t1 and t2, the covariance, for each component
timepoint_meancor(sim, t1, t2) # Computes the mean at t1 and t2, the correlation, for each component
timepoint_weighted_meancov(sim, W, t1, t2) # Computes the mean at t1 and t2, the weighted covariance W, for each component
Full Timeseries Statistics
Additionally, the following functions are provided for analyzing the full timeseries. The mean
and meanvar
versions return a DiffEqArray
which can be directly plotted. The meancov
and meancor
return a matrix of tuples, where the tuples are the (mean_t1,mean_t2,cov or cor)
.
The available functions for the time steps are:
timeseries_steps_mean(sim) # Computes the mean at each time step
timeseries_steps_median(sim) # Computes the median at each time step
timeseries_steps_quantile(sim, q) # Computes the quantile q at each time step
timeseries_steps_meanvar(sim) # Computes the mean and variance at each time step
timeseries_steps_meancov(sim) # Computes the covariance matrix and means at each time step
timeseries_steps_meancor(sim) # Computes the correlation matrix and means at each time step
timeseries_steps_weighted_meancov(sim) # Computes the weighted covariance matrix and means at each time step
The available functions for the time points are:
timeseries_point_mean(sim, ts) # Computes the mean at each time point in ts
timeseries_point_median(sim, ts) # Computes the median at each time point in ts
timeseries_point_quantile(sim, q, ts) # Computes the quantile q at each time point in ts
timeseries_point_meanvar(sim, ts) # Computes the mean and variance at each time point in ts
timeseries_point_meancov(sim, ts) # Computes the covariance matrix and means at each time point in ts
timeseries_point_meancor(sim, ts) # Computes the correlation matrix and means at each time point in ts
timeseries_point_weighted_meancov(sim, ts) # Computes the weighted covariance matrix and means at each time point in ts
EnsembleSummary
The EnsembleSummary
type is included to help with analyzing the general summary statistics. Two constructors are provided:
EnsembleSummary(sim; quantiles = [0.05, 0.95])
EnsembleSummary(sim, ts; quantiles = [0.05, 0.95])
The first produces a (mean,var)
summary at each time step. As with the summary statistics, this assumes that the time steps are all the same. The second produces a (mean,var)
summary at each time point t
in ts
. This requires the ability to interpolate the solution. Quantile is used to determine the qlow
and qhigh
quantiles at each timepoint. It defaults to the 5% and 95% quantiles.
Plot Recipe
The EnsembleSummary
comes with a plot recipe for visualizing the summary statistics. The extra keyword arguments are:
-
idxs
: the solution components to plot. Defaults to plotting all components. -
error_style
: The style for plotting the error. Defaults toribbon
. Other choices are:bars
for error bars and:none
for no error bars. -
ci_type
: Defaults to:quantile
which has(qlow,qhigh)
quantiles whose limits were determined when constructing theEnsembleSummary
. Gaussian CI1.96*(standard error of the mean)
can be set usingci_type=:SEM
.
One useful argument is fillalpha
which controls the transparency of the ribbon around the mean.
Example 1: Solving an ODE With Different Initial Conditions
Random Initial Conditions
Let’s test the sensitivity of the linear ODE to its initial condition. To do this, we would like to solve the linear ODE 100 times and plot what the trajectories look like. Let’s start by opening up some extra processes so that way the computation will be parallelized. Here we will choose to use distributed parallelism, which means that the required functions must be made available to all processes. This can be achieved with @everywhere
macro:
using Distributed
using DifferentialEquations
using Plots
addprocs()
@everywhere using DifferentialEquations
Now let’s define the linear ODE, which is our base problem:
# Linear ODE which starts at 0.5 and solves from t=0.0 to t=1.0
prob = ODEProblem((u, p, t) -> 1.01u, 0.5, (0.0, 1.0))
For our ensemble simulation, we would like to change the initial condition around. This is done through the prob_func
. This function takes in the base problem and modifies it to create the new problem that the trajectory actually solves. Here, we will take the base problem, multiply the initial condition by a rand()
, and use that for calculating the trajectory:
@everywhere function prob_func(prob, i, repeat)
remake(prob, u0 = rand() * prob.u0)
end
Now we build and solve the EnsembleProblem
with this base problem and prob_func
:
ensemble_prob = EnsembleProblem(prob, prob_func = prob_func)
sim = solve(ensemble_prob, Tsit5(), EnsembleDistributed(), trajectories = 10)
We can use the plot recipe to plot what the 10 ODEs look like:
plot(sim, linealpha = 0.4)
We note that if we wanted to find out what the initial condition was for a given trajectory, we can retrieve it from the solution. sim[i]
returns the i
th solution object. sim[i].prob
is the problem that specific trajectory solved, and sim[i].prob.u0
would then be the initial condition used in the i
th trajectory.
If the problem has callbacks, the functions for the condition and affect! must be named functions (not anonymous functions).
|
Using multithreading
The previous ensemble simulation can also be parallelized using a multithreading approach, which will make use of the different cores within a single computer. Because the memory is shared across the different threads, it is not necessary to use the @everywhere
macro. Instead, the same problem can be implemented simply as:
using DifferentialEquations
prob = ODEProblem((u, p, t) -> 1.01u, 0.5, (0.0, 1.0))
function prob_func(prob, i, repeat)
remake(prob, u0 = rand() * prob.u0)
end
ensemble_prob = EnsembleProblem(prob, prob_func = prob_func)
sim = solve(ensemble_prob, Tsit5(), EnsembleThreads(), trajectories = 10)
using Plots;
plot(sim);
The number of threads to be used has to be defined outside of Julia, in the environmental variable JULIA_NUM_THREADS
(see Julia’s documentation for details).
Pre-Determined Initial Conditions
Often, you may already know what initial conditions you want to use. This can be specified by the i
argument of the prob_func
. This i
is the unique index of each trajectory. So, if we have trajectories=100
, then we have i
as some index in 1:100
, and it’s different for each trajectory.
So, if we wanted to use a grid of evenly spaced initial conditions from 0
to 1
, we could simply index the linspace
type:
initial_conditions = range(0, stop = 1, length = 100)
function prob_func(prob, i, repeat)
remake(prob, u0 = initial_conditions[i])
end
prob_func (generic function with 1 method)
It’s worth noting that if you run this code successfully, there will be no visible output.
Example 2: Solving an SDE with Different Parameters
Let’s solve the same SDE, but with varying parameters. Let’s create a Lotka-Volterra system with multiplicative noise. Our Lotka-Volterra system will have as its drift component:
function f(du, u, p, t)
du[1] = p[1] * u[1] - p[2] * u[1] * u[2]
du[2] = -3 * u[2] + u[1] * u[2]
end
f (generic function with 1 method)
For our noise function, we will use multiplicative noise:
function g(du, u, p, t)
du[1] = p[3] * u[1]
du[2] = p[4] * u[2]
end
g (generic function with 1 method)
Now we build the SDE with these functions:
using DifferentialEquations
p = [1.5, 1.0, 0.1, 0.1]
prob = SDEProblem(f, g, [1.0, 1.0], (0.0, 10.0), p)
SDEProblem with uType Vector{Float64} and tType Float64. In-place: true
timespan: (0.0, 10.0)
u0: 2-element Vector{Float64}:
1.0
1.0
This is the base problem for our study. What would like to do with this experiment is keep the same parameters in the deterministic component each time, but vary the parameters for the amount of noise using 0.3rand(2)
as our parameters. Once again, we do this with a prob_func
, and here we modify the parameters in prob.p
:
# `p` is a global variable, referencing it would be type unstable.
# Using a let block defines a small local scope in which we can
# capture that local `p` which isn't redefined anywhere in that local scope.
# This allows it to be type stable.
prob_func = let p = p
(prob, i, repeat) -> begin
x = 0.3rand(2)
remake(prob, p = [p[1], p[2], x[1], x[2]])
end
end
#1 (generic function with 1 method)
Now we solve the problem 10 times and plot all of the trajectories in phase space:
ensemble_prob = EnsembleProblem(prob, prob_func = prob_func)
sim = solve(ensemble_prob, SRIW1(), trajectories = 10)
using Plots;
plot(sim, linealpha = 0.6, color = :blue, idxs = (0, 1), title = "Phase Space Plot");
plot!(sim, linealpha = 0.6, color = :red, idxs = (0, 2), title = "Phase Space Plot")
We can then summarize this information with the mean/variance bounds using a EnsembleSummary
plot. We will take the mean/quantile at every 0.1
time units and directly plot the summary:
summ = EnsembleSummary(sim, 0:0.1:10)
plot(summ, fillalpha = 0.5)
Note that here we used the quantile bounds, which default to [0.05,0.95]
in the EnsembleSummary
constructor. We can change to standard error of the mean bounds using ci_type=:SEM
in the plot recipe.
Example 3: Using the Reduction to Halt When Estimator is Within Tolerance
In this problem, we will solve the equation just as many times as needed to get the standard error of the mean for the final time point below our tolerance 0.5
. Since we only care about the endpoint, we can tell the output_func
to discard the rest of the data.
function output_func(sol, i)
last(sol), false
end
output_func (generic function with 1 method)
Our prob_func
will simply randomize the initial condition:
using DifferentialEquations
# Linear ODE which starts at 0.5 and solves from t=0.0 to t=1.0
prob = ODEProblem((u, p, t) -> 1.01u, 0.5, (0.0, 1.0))
function prob_func(prob, i, repeat)
remake(prob, u0 = rand() * prob.u0)
end
prob_func (generic function with 1 method)
Our reduction function will append the data from the current batch to the previous batch, and declare convergence if the standard error of the mean is calculated as sufficiently small:
using Statistics
function reduction(u, batch, I)
u = append!(u, batch)
finished = (var(u) / sqrt(last(I))) / mean(u) < 0.5
u, finished
end
reduction (generic function with 1 method)
Then we can define and solve the problem:
prob2 = EnsembleProblem(prob, prob_func = prob_func, output_func = output_func,
reduction = reduction, u_init = Vector{Float64}())
sim = solve(prob2, Tsit5(), trajectories = 10000, batch_size = 20)
EnsembleSolution Solution of length 20 with uType:
Float64
Since batch_size=20
, this means that every 20 simulations, it will take this batch, append the results to the previous batch, calculate (var(u)/sqrt(last(I)))/mean(u)
, and if that’s small enough, exit the simulation. In this case, the simulation exits only after 20 simulations (i.e. after calculating the first batch). This can save a lot of time!
In addition to saving time by checking convergence, we can save memory by reducing between batches. For example, say we only care about the mean at the end once again. Instead of saving the solution at the end for each trajectory, we can instead save the running summation of the endpoints:
function reduction(u, batch, I)
u + sum(batch), false
end
prob2 = EnsembleProblem(prob, prob_func = prob_func, output_func = output_func,
reduction = reduction, u_init = 0.0)
sim2 = solve(prob2, Tsit5(), trajectories = 100, batch_size = 20)
EnsembleSolution Solution of length 1 with uType:
Float64
this will sum up the endpoints after every 20 solutions, and save the running sum. The final result will have sim2.u
as simply a number, and thus sim2.u/100
would be the mean.
Example 4: Using the Analysis Tools
In this example, we will show how to analyze a EnsembleSolution
. First, let’s generate a 10 solution Monte Carlo experiment. For our problem, we will use a 4x2
system of linear stochastic differential equations:
function f(du, u, p, t)
for i in 1:length(u)
du[i] = 1.01 * u[i]
end
end
function σ(du, u, p, t)
for i in 1:length(u)
du[i] = 0.87 * u[i]
end
end
using DifferentialEquations
prob = SDEProblem(f, σ, ones(4, 2) / 2, (0.0, 1.0)) #prob_sde_2Dlinear
SDEProblem with uType Matrix{Float64} and tType Float64. In-place: true
timespan: (0.0, 1.0)
u0: 4×2 Matrix{Float64}:
0.5 0.5
0.5 0.5
0.5 0.5
0.5 0.5
To solve this 10 times, we use the EnsembleProblem
constructor and solve with trajectories=10
. Since we wish to compare values at the timesteps, we need to make sure the steps all hit the same times. We thus set adaptive=false
and explicitly give a dt
.
prob2 = EnsembleProblem(prob)
sim = solve(prob2, SRIW1(), dt = 1 // 2^(3), trajectories = 10, adaptive = false)
EnsembleSolution Solution of length 10 with uType:
RODESolution{Float64, 3, Vector{Matrix{Float64}}, Nothing, Nothing, Vector{Float64}, NoiseProcess{Float64, 3, Float64, Matrix{Float64}, Matrix{Float64}, Vector{Matrix{Float64}}, typeof(DiffEqNoiseProcess.INPLACE_WHITE_NOISE_DIST), typeof(DiffEqNoiseProcess.INPLACE_WHITE_NOISE_BRIDGE), true, ResettableStacks.ResettableStack{Tuple{Float64, Matrix{Float64}, Matrix{Float64}}, true}, ResettableStacks.ResettableStack{Tuple{Float64, Matrix{Float64}, Matrix{Float64}}, true}, RSWM{Float64}, Nothing, RandomNumbers.Xorshifts.Xoroshiro128Plus}, SDEProblem{Matrix{Float64}, Tuple{Float64, Float64}, true, SciMLBase.NullParameters, Nothing, SDEFunction{true, SciMLBase.FullSpecialize, typeof(Main.f), typeof(Main.σ), UniformScaling{Bool}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing}, typeof(Main.σ), Base.Pairs{Symbol, Union{}, Tuple{}, NamedTuple{(), Tuple{}}}, Nothing}, SRIW1, StochasticDiffEq.LinearInterpolationData{Vector{Matrix{Float64}}, Vector{Float64}}, DiffEqBase.Stats, Nothing}
Note that if you don’t do the timeseries_steps
calculations, this code is compatible with adaptive timestepping. Using adaptivity is usually more efficient!
We can compute the mean and the variance at the 3rd timestep using:
using DifferentialEquations.EnsembleAnalysis
m, v = timestep_meanvar(sim, 3)
([0.7059438825839033 0.507097988064184; 0.8473751737803525 0.7024607895232441; 0.685934978893289 0.5487963290163301; 0.6326564881228885 0.8246258853792245], [0.08220192493826829 0.04094223829627257; 0.15357597542142784 0.07643987057982006; 0.05241943872642733 0.037769230968575186; 0.0903645667625475 0.07902457497498139])
or we can compute the mean and the variance at the t=0.5
using:
m, v = timepoint_meanvar(sim, 0.5)
([0.9515887327672015 0.6600224705927742; 1.0559274818316131 0.8219629156937067; 0.6319679185906485 0.7997168492442293; 0.7880845795244602 0.9886773300584073], [0.18470071197029475 0.22807965935937002; 0.556424490399348 0.2086645242977786; 0.09824942288195918 0.6724025460848169; 0.22558854357279806 0.2526245167702519])
We can get a series for the mean and the variance at each time step using:
m_series, v_series = timeseries_steps_meanvar(sim)
(DiffEqArray{Float64, 3, Vector{Matrix{Float64}}, Vector{Float64}, SymbolicIndexingInterface.SymbolCache{Nothing, Nothing, Nothing}, Nothing, Nothing}([[0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [0.5907422480870421 0.47670059358761757; 0.620331566216781 0.5617611392227679; 0.6344684552006176 0.5627747292169305; 0.5477639822273335 0.6922339178643483], [0.7059438825839033 0.507097988064184; 0.8473751737803525 0.7024607895232441; 0.685934978893289 0.5487963290163301; 0.6326564881228885 0.8246258853792245], [0.7341285196161411 0.5345032287984489; 1.0208464804277884 0.6707094238040601; 0.7242417684224293 0.552797007551382; 0.6742516122474017 0.8597673851846357], [0.9515887327672015 0.6600224705927742; 1.0559274818316131 0.8219629156937067; 0.6319679185906485 0.7997168492442293; 0.7880845795244602 0.9886773300584073], [1.0428237451034787 0.667282003916669; 1.1378286087828915 0.938346903955284; 0.6427318009602387 0.7313467864263105; 0.840679501086577 1.149008662177291], [1.1933448270652285 0.8619738643156023; 1.141007940476052 0.977291388985809; 0.8737169860152547 0.7710177755726961; 1.0955263068343075 1.1723661063580615], [1.4917319621735583 0.9803049865605309; 1.491562258522121 1.1865468788615872; 0.9935685698835602 0.9457012518264655; 1.2596297188646093 1.5185848360656915], [1.9830491482389645 1.5024128142386555; 1.6296115539269103 1.6983611117680775; 1.0588790764004232 1.048712234611464; 1.3746252663438125 1.632360373184002]], [0.0, 0.125, 0.25, 0.375, 0.5, 0.625, 0.75, 0.875, 1.0], SymbolicIndexingInterface.SymbolCache{Nothing, Nothing, Nothing}(nothing, nothing, nothing), nothing, nothing), DiffEqArray{Float64, 3, Vector{Matrix{Float64}}, Vector{Float64}, SymbolicIndexingInterface.SymbolCache{Nothing, Nothing, Nothing}, Nothing, Nothing}([[0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0], [0.01776004619856157 0.017150069440896233; 0.007464822899251548 0.010749563070496258; 0.0503692942338929 0.023990014895319062; 0.025565193824098633 0.050221687927936645], [0.08220192493826829 0.04094223829627257; 0.15357597542142784 0.07643987057982006; 0.05241943872642733 0.037769230968575186; 0.0903645667625475 0.07902457497498139], [0.08908890538434039 0.06738722755177229; 0.48632108007863617 0.09549009454811547; 0.09785735398242405 0.056367823651902085; 0.09846163931517155 0.13703898146914115], [0.18470071197029475 0.22807965935937002; 0.556424490399348 0.2086645242977786; 0.09824942288195918 0.6724025460848169; 0.22558854357279806 0.2526245167702519], [0.4014736122621857 0.29779995530296133; 0.6236246385491712 0.290802223589946; 0.1100272342836092 0.2877908832401255; 0.3479755076551748 0.6281619192815174], [0.6876532727219186 0.6121096811643733; 0.6209063336971192 0.45991376000562806; 0.27684505476158017 0.57802669702141; 0.7164492474489248 0.43014081986211994], [1.5815841149730185 0.9499215307307226; 2.338266758810592 0.8537427707209551; 0.38879029488916766 1.2119005556305449; 0.8389201622044382 0.7775070942591135], [3.10917079628439 6.148445778195559; 1.8175196656258836 2.8296611898538258; 0.46502349756642875 1.4111594581364308; 0.876508598760132 1.494240020464408]], [0.0, 0.125, 0.25, 0.375, 0.5, 0.625, 0.75, 0.875, 1.0], SymbolicIndexingInterface.SymbolCache{Nothing, Nothing, Nothing}(nothing, nothing, nothing), nothing, nothing))
or at chosen values of t
:
ts = 0:0.1:1
m_series = timeseries_point_mean(sim, ts)
t: 0.0:0.1:1.0
u: 11-element Vector{Matrix{Float64}}:
[0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5]
[0.5725937984696337 0.4813604748700941; 0.5962652529734249 0.5494089113782143; 0.6075747641604942 0.5502197833735446; 0.5382111857818669 0.6537871342914786]
[0.6598632287851589 0.49493903027355757; 0.756557730754924 0.6461809294030537; 0.6653483694162206 0.5543876890965702; 0.5986994857646666 0.7716690983732741]
[0.7172177373967985 0.5180600843578901; 0.9167636964393271 0.6897602432355704; 0.7012576947049451 0.5503966004303509; 0.6492945377726939 0.8386824853013891]
[0.7776205622463533 0.559607077157314; 1.0278626807085534 0.7009601221819894; 0.7057869984560732 0.6021809758899515; 0.6970182057028135 0.8855493741593901]
[0.9515887327672015 0.6600224705927742; 1.0559274818316131 0.8219629156937069; 0.6319679185906486 0.7997168492442291; 0.7880845795244602 0.9886773300584073]
[1.0245767426362233 0.66583009725189; 1.1214483833926359 0.9150701063029685; 0.6405790244863209 0.7450207989898943; 0.8301605167741535 1.116942395753514]
[1.1331363942805286 0.7840971201560288; 1.1397362077987878 0.9617135949735989; 0.7813229119932484 0.7551493799141419; 0.993587584535215 1.1630231286857533]
[1.3126996811085603 0.9093063132135736; 1.2812296676944797 1.0609935849361203; 0.921657619562577 0.840891166074204; 1.161167671646428 1.3108535982411138]
[1.5899953993866396 1.0847265520961558; 1.519172117603079 1.2889097254428852; 1.0066306711869328 0.9663034483834654; 1.28262882836045 1.5413399434893538]
[1.9830491482389647 1.5024128142386552; 1.6296115539269107 1.6983611117680772; 1.0588790764004234 1.0487122346114641; 1.3746252663438128 1.632360373184002]
Note that these mean and variance series can be directly plotted. We can compute covariance matrices similarly:
timeseries_steps_meancov(sim) # Use the time steps, assume fixed dt
timeseries_point_meancov(sim, 0:(1 // 2^(3)):1, 0:(1 // 2^(3)):1) # Use time points, interpolate
9×9 Matrix{Tuple{Matrix{Float64}, Matrix{Float64}, Matrix{Float64}}}:
([0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0]) … ([1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0])
([0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [0.590742 0.476701; 0.620332 0.561761; 0.634468 0.562775; 0.547764 0.692234], [0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0]) ([1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [0.590742 0.476701; 0.620332 0.561761; 0.634468 0.562775; 0.547764 0.692234], [0.0967682 -0.0223925; 0.0274006 0.0324591; 0.0235797 0.0275256; 0.0405064 -0.0538027])
([0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [0.705944 0.507098; 0.847375 0.702461; 0.685935 0.548796; 0.632656 0.824626], [0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0]) ([1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [0.705944 0.507098; 0.847375 0.702461; 0.685935 0.548796; 0.632656 0.824626], [0.162477 0.190289; 0.405919 0.344204; 0.0838079 0.0995001; 0.0918312 -0.0158649])
([0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [0.734129 0.534503; 1.02085 0.670709; 0.724242 0.552797; 0.674252 0.859767], [0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0]) ([1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [0.734129 0.534503; 1.02085 0.670709; 0.724242 0.552797; 0.674252 0.859767], [0.400285 0.489974; 0.773062 0.423663; 0.111857 0.241149; 0.12542 0.160866])
([0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [0.951589 0.660022; 1.05593 0.821963; 0.631968 0.799717; 0.788085 0.988677], [0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0]) ([1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [0.951589 0.660022; 1.05593 0.821963; 0.631968 0.799717; 0.788085 0.988677], [0.663387 0.98501; 0.916895 0.664875; 0.105157 0.916273; 0.256409 0.316028])
([0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [1.04282 0.667282; 1.13783 0.938347; 0.642732 0.731347; 0.84068 1.14901], [0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0]) … ([1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [1.04282 0.667282; 1.13783 0.938347; 0.642732 0.731347; 0.84068 1.14901], [0.967997 1.15583; 0.993903 0.69763; 0.181501 0.597774; 0.307125 0.290157])
([0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [1.19334 0.861974; 1.14101 0.977291; 0.873717 0.771018; 1.09553 1.17237], [0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0]) ([1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [1.19334 0.861974; 1.14101 0.977291; 0.873717 0.771018; 1.09553 1.17237], [1.26753 1.57324; 0.944408 0.923347; 0.30119 0.878651; 0.591724 0.423086])
([0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [1.49173 0.980305; 1.49156 1.18655; 0.993569 0.945701; 1.25963 1.51858], [0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0]) ([1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [1.49173 0.980305; 1.49156 1.18655; 0.993569 0.945701; 1.25963 1.51858], [2.1187 2.16082; 1.98277 1.46278; 0.398108 1.2802; 0.84249 0.824917])
([0.5 0.5; 0.5 0.5; 0.5 0.5; 0.5 0.5], [1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [0.0 0.0; 0.0 0.0; 0.0 0.0; 0.0 0.0]) ([1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [1.98305 1.50241; 1.62961 1.69836; 1.05888 1.04871; 1.37463 1.63236], [3.10917 6.14845; 1.81752 2.82966; 0.465023 1.41116; 0.876509 1.49424])
For general analysis, we can build a EnsembleSummary
type.
summ = EnsembleSummary(sim)
EnsembleSolution Solution of length 9 with uType:
Float64
will summarize at each time step, while
summ = EnsembleSummary(sim, 0.0:0.1:1.0)
EnsembleSolution Solution of length 11 with uType:
Float64
will summarize at the 0.1
time points using the interpolations. To visualize the results, we can plot it. Since there are 8 components to the differential equation, this can get messy, so let’s only plot the 3rd component:
using Plots;
plot(summ; idxs = 3);
We can change to errorbars instead of ribbons and plot two different indices:
plot(summ; idxs = (3, 5), error_style = :bars)
Or we can simply plot the mean of every component over time:
plot(summ; error_style = :none)