I/O: Saving and Loading Solution Data

The ability to save and load solutions is important for handling large datasets and analyzing the results over multiple Julia sessions. This page explains the existing functionality for doing so.

Tabular Data: IterableTables

An interface to IterableTables.jl is provided. This IterableTables link allows you to use a solution type as the data source to convert to other tabular data formats. For example, let's solve a 4x2 system of ODEs and get the DataFrame:

import OrdinaryDiffEq as ODE, DataFrames
f_2dlinear = (du, u, p, t) -> du .= 1.01u;
tspan = (0.0, 1.0)
prob = ODE.ODEProblem(f_2dlinear, rand(2, 2), tspan);
sol = ODE.solve(prob, ODE.Euler(); dt = 1 // 2^(4));
df = DataFrames.DataFrame(sol)
17×5 DataFrame
Rowtimestampvalue1value2value3value4
Float64Float64Float64Float64Float64
10.00.7760810.9568580.5722260.538778
20.06250.8250721.017260.6083480.572788
30.1250.8771541.081470.646750.608946
40.18750.9325251.149740.6875760.647385
50.250.991391.222320.7309790.688251
60.31251.053971.299480.7771220.731697
70.3751.12051.381510.8261780.777886
80.43751.191241.468720.878330.82699
90.51.266431.561430.9337750.879194
100.56251.346381.659990.9927190.934693
110.6251.431371.764781.055380.993695
120.68751.521721.876181.122011.05642
130.751.617781.994621.192831.12311
140.81251.71992.120531.268131.194
150.8751.828472.254381.348181.26938
160.93751.943892.396691.433281.34951
171.02.06662.547981.523761.43469

If we set syms in the DiffEqFunction, then those names will be used:

f = ODE.ODEFunction(f_2dlinear, syms = [:a, :b, :c, :d])
prob = ODE.ODEProblem(f, rand(2, 2), (0.0, 1.0));
sol = ODE.solve(prob, ODE.Euler(); dt = 1 // 2^(4));
df = DataFrames.DataFrame(sol)
17×5 DataFrame
Rowtimestampabcd
Float64Float64Float64Float64Float64
10.00.5892980.4423730.6862050.0995601
20.06250.6264970.4702970.7295210.105845
30.1250.6660450.4999850.7755720.112526
40.18750.7080890.5315460.824530.11963
50.250.7527870.56510.8765790.127181
60.31250.8003070.6007720.9319130.135209
70.3750.8508260.6386960.990740.143745
80.43750.9045350.6790141.053280.152818
90.50.9616330.7218761.119770.162465
100.56251.022340.7674451.190450.172721
110.6251.086870.815891.26560.183624
120.68751.155480.8673931.345490.195215
130.751.228420.9221471.430430.207538
140.81251.305960.9803581.520720.220639
150.8751.38841.042241.616720.234567
160.93751.476051.108031.718770.249374
171.01.569221.177981.827270.265115

Many modeling frameworks will automatically set syms for this feature. Additionally, this data can be saved to a CSV:

import CSV
CSV.write("out.csv", df)
"out.csv"

JLD2 and BSON.jl

JLD2.jl and BSON.jl will work with the full solution type if you bring the required functions back into scope before loading. For example, if we save the solution:

sol = ODE.solve(prob, ODE.Euler(); dt = 1 // 2^(4))
import JLD2
JLD2.@save "out.jld2" sol
┌ Warning: Attempting to store Main.var"#1#2".
JLD2 only stores functions by name.
 This may not be useful for anonymous functions.
@ JLD2 ~/.cache/julia-buildkite-plugin/depots/0185fce3-4489-413a-a934-123dd653ef61/packages/JLD2/SgtOb/src/data/writing_datatypes.jl:447
┌ Warning: Attempting to store Main.var"#1#2".
JLD2 only stores functions by name.
 This may not be useful for anonymous functions.
@ JLD2 ~/.cache/julia-buildkite-plugin/depots/0185fce3-4489-413a-a934-123dd653ef61/packages/JLD2/SgtOb/src/data/writing_datatypes.jl:447

then we can get the full solution type back, interpolations and all, if we load the dependent functions first:

# New session
import JLD2
import OrdinaryDiffEq as ODE
JLD2.@load "out.jld2" sol
1-element Vector{Symbol}:
 :sol

The example with BSON.jl is:

sol = ODE.solve(prob, ODE.Euler(); dt = 1 // 2^(4))
import BSON
BSON.bson("test.bson", Dict(:sol => sol))
# New session
import OrdinaryDiffEq as ODE
import BSON
# BSON.load("test.bson") # currently broken: https://github.com/JuliaIO/BSON.jl/issues/109

If you load it without the DE function then for some algorithms the interpolation may not work, and for all algorithms you'll need at least a solver package or SciMLBase.jl in scope in order for the solution interface (plot recipes, array indexing, etc.) to work. If none of these are put into scope, the solution type will still load and hold all of the values (so sol.u and sol.t will work), but none of the interface will be available.

If you want a copy of the solution that contains no function information you can use the function SciMLBase.strip_solution(sol). This will return a copy of the solution that doesn't have any functions, which you can serialize and deserialize without having any of the problems that typically come with serializing functions.

JLD

Don't use JLD. It's dead. Julia types can be saved via JLD.jl. However, they cannot save types which have functions, which means that the solution type is currently not compatible with JLD.

import JLD
JLD.save("out.jld", "sol", sol)