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.4053290.7862470.2970470.0910256
20.06250.4309160.8358790.3157990.0967715
30.1250.4581170.8886430.3357330.10288
40.18750.4870360.9447390.3569270.109375
50.250.517781.004380.3794580.116279
60.31250.5504651.067780.4034110.123619
70.3750.5852131.135180.4288760.131422
80.43750.6221551.206840.4559490.139718
90.50.6614281.283020.4847310.148538
100.56250.7031811.364010.5153290.157915
110.6250.7475691.450110.5478590.167883
120.68750.794761.541650.5824430.178481
130.750.8449291.638970.619210.189747
140.81250.8982651.742430.6582970.201725
150.8750.9549681.852420.6998520.214459
160.93751.015251.969350.7440310.227997
171.01.079342.093670.7909980.242389

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.2485020.5690230.5782610.413703
20.06250.2641880.6049420.6147630.439818
30.1250.2808650.6431290.653570.467581
40.18750.2985950.6837270.6948270.497097
50.250.3174440.7268870.7386880.528477
60.31250.3374820.7727720.7853170.561837
70.3750.3587860.8215530.8348910.597303
80.43750.3814340.8734130.8875930.635008
90.50.4055120.9285480.9436220.675092
100.56250.431110.9871621.003190.717708
110.6250.4583241.049481.066510.763013
120.68750.4872561.115731.133840.811178
130.750.5180141.186161.205410.862384
140.81250.5507131.261031.28150.916822
150.8750.5854771.340631.36240.974696
160.93750.6224351.425261.44841.03622
171.00.6617271.515231.539831.10164

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"#2#3".
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/hbsZG/src/data/writing_datatypes.jl:447
┌ Warning: Attempting to store Main.var"#2#3".
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/hbsZG/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)