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.5623230.4094580.1717850.558916
20.06250.597820.4353050.1826290.594197
30.1250.6355570.4627830.1941570.631706
40.18750.6756770.4919970.2064130.671582
50.250.7183290.5230540.2194430.713976
60.31250.7636730.5560720.2332960.759046
70.3750.811880.5911740.2480220.80696
80.43750.863130.6284910.2636790.8579
90.50.9176150.6681650.2803240.912055
100.56250.975540.7103430.2980190.969628
110.6251.037120.7551830.3168311.03084
120.68751.102590.8028540.3368311.09591
130.751.172190.8535340.3580941.16509
140.81251.246180.9074140.3806991.23863
150.8751.324850.9646940.404731.31682
160.93751.408481.025590.4302791.39995
171.01.497391.090330.457441.48832

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.4171220.5519110.4800570.418108
20.06250.4434520.586750.5103610.444502
30.1250.4714450.6237890.5425780.472561
40.18750.5012050.6631660.5768280.502391
50.250.5328440.7050280.613240.534104
60.31250.566480.7495330.6519510.56782
70.3750.6022390.7968470.6931050.603663
80.43750.6402550.8471480.7368570.64177
90.50.6806710.9006250.7833720.682281
100.56250.7236380.9574760.8328220.72535
110.6250.7693181.017920.8853940.771138
120.68750.8178811.082170.9412840.819816
130.750.869511.150491.00070.871567
140.81250.9243981.223111.063870.926585
150.8750.9827511.300321.131030.985076
160.93751.044791.38241.202431.04726
171.01.110741.469671.278331.11337

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)