Skip to content

Commit d33f86b

Browse files
Merge pull request #1178 from AayushSabharwal/as/cse-bench
feat: add benchmark for effect of CSE on jacobian computation/compilation
2 parents 4bf0443 + 386dc8a commit d33f86b

File tree

5 files changed

+4250
-0
lines changed

5 files changed

+4250
-0
lines changed

benchmarks/Symbolics/BCR.jmd

Lines changed: 241 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,241 @@
1+
---
2+
title: BCR Symbolic Jacobian
3+
author: Aayush Sabharwal, Bowen Zhu, Chris Rackauckas
4+
---
5+
6+
The following benchmark is of 1122 ODEs with 24388 terms that describe a stiff
7+
chemical reaction network modeling the BCR signaling network from [Barua et
8+
al.](https://doi.org/10.4049/jimmunol.1102003). We use
9+
[`ReactionNetworkImporters`](https://github.com/isaacsas/ReactionNetworkImporters.jl)
10+
to load the BioNetGen model files as a
11+
[Catalyst](https://github.com/SciML/Catalyst.jl) model, and then use
12+
[ModelingToolkit](https://github.com/SciML/ModelingToolkit.jl) to convert the
13+
Catalyst network model to ODEs.
14+
15+
The resultant large model is used to benchmark the time taken to compute a symbolic
16+
jacobian, generate a function to calculate it and call the function.
17+
18+
19+
```julia
20+
using Catalyst, ReactionNetworkImporters,
21+
TimerOutputs, LinearAlgebra, ModelingToolkit, Chairmarks,
22+
LinearSolve, Symbolics, SymbolicUtils, SymbolicUtils.Code, SparseArrays, CairoMakie,
23+
PrettyTables
24+
25+
datadir = joinpath(dirname(pathof(ReactionNetworkImporters)),"../data/bcr")
26+
const to = TimerOutput()
27+
tf = 100000.0
28+
29+
# generate ModelingToolkit ODEs
30+
prnbng = loadrxnetwork(BNGNetwork(), joinpath(datadir, "bcr.net"))
31+
show(to)
32+
rn = complete(prnbng.rn; split = false)
33+
obs = [eq.lhs for eq in observed(rn)]
34+
osys = convert(ODESystem, rn)
35+
36+
rhs = [eq.rhs for eq in full_equations(osys)]
37+
vars = unknowns(osys)
38+
pars = parameters(osys)
39+
```
40+
41+
The `sparsejacobian` function in Symbolics.jl is optimized for hashconsing and caching, and as such
42+
performs very poorly without either of those features. We use the old implementation, optimized without
43+
hashconsing, to benchmark performance without hashconsing and without caching to avoid biasing the results.
44+
45+
```julia
46+
include("old_sparse_jacobian.jl")
47+
```
48+
49+
```julia
50+
SymbolicUtils.ENABLE_HASHCONSING[] = false
51+
@timeit to "Calculate jacobian - without hashconsing" jac_nohc = old_sparsejacobian(rhs, vars);
52+
SymbolicUtils.ENABLE_HASHCONSING[] = true
53+
Symbolics.toggle_derivative_caching!(false)
54+
Symbolics.clear_derivative_caches!()
55+
@timeit to "Calculate jacobian - hashconsing, without caching" jac_hc_nocache = old_sparsejacobian(rhs, vars);
56+
Symbolics.toggle_derivative_caching!(true)
57+
for fn in Symbolics.cached_derivative_functions()
58+
stats = SymbolicUtils.get_stats(fn)
59+
@assert stats.hits == stats.misses == 0
60+
end
61+
Symbolics.clear_derivative_caches!()
62+
@timeit to "Calculate jacobian - hashconsing and caching" jac_hc_cache = Symbolics.sparsejacobian(rhs, vars);
63+
64+
@assert isequal(jac_nohc, jac_hc_nocache)
65+
@assert isequal(jac_hc_nocache, jac_hc_cache)
66+
67+
jac = jac_hc_cache
68+
args = (vars, pars, ModelingToolkit.get_iv(osys))
69+
# out of place versions run into an error saying the expression is too large
70+
# due to the `SymbolicUtils.Code.create_array` call. `iip_config` prevents it
71+
# from trying to build the function.
72+
kwargs = (; iip_config = (false, true), expression = Val{true})
73+
@timeit to "Build jacobian - no CSE" _, jac_nocse_iip = build_function(jac, args...; cse = false, kwargs...);
74+
@timeit to "Build jacobian - CSE" _, jac_cse_iip = build_function(jac, args...; cse = true, kwargs...);
75+
76+
jac_nocse_iip = eval(jac_nocse_iip)
77+
jac_cse_iip = eval(jac_cse_iip)
78+
79+
defs = defaults(osys)
80+
u = Float64[Symbolics.fixpoint_sub(var, defs) for var in vars]
81+
buffer_cse = similar(jac, Float64)
82+
buffer_nocse = similar(jac, Float64)
83+
p = Float64[Symbolics.fixpoint_sub(par, defs) for par in pars]
84+
tt = 0.0
85+
86+
@timeit to "Compile jacobian - CSE" jac_cse_iip(buffer_cse, u, p, tt)
87+
@timeit to "Compute jacobian - CSE" jac_cse_iip(buffer_cse, u, p, tt)
88+
89+
@timeit to "Compile jacobian - no CSE" jac_nocse_iip(buffer_nocse, u, p, tt)
90+
@timeit to "Compute jacobian - no CSE" jac_nocse_iip(buffer_nocse, u, p, tt)
91+
92+
@assert isapprox(buffer_cse, buffer_nocse, rtol = 1e-10)
93+
94+
show(to)
95+
```
96+
97+
We'll also measure scaling.
98+
99+
100+
```julia
101+
function run_and_time_construct!(rhs, vars, pars, iv, N, i, jac_times, jac_allocs, build_times, functions)
102+
outputs = rhs[1:N]
103+
SymbolicUtils.ENABLE_HASHCONSING[] = false
104+
jac_result = @be old_sparsejacobian(outputs, vars)
105+
jac_times[1][i] = minimum(x -> x.time, jac_result.samples)
106+
jac_allocs[1][i] = minimum(x -> x.bytes, jac_result.samples)
107+
108+
SymbolicUtils.ENABLE_HASHCONSING[] = true
109+
jac_result = @be (Symbolics.clear_derivative_caches!(); Symbolics.sparsejacobian(outputs, vars))
110+
jac_times[2][i] = minimum(x -> x.time, jac_result.samples)
111+
jac_allocs[2][i] = minimum(x -> x.bytes, jac_result.samples)
112+
113+
jac = Symbolics.sparsejacobian(outputs, vars)
114+
args = (vars, pars, iv)
115+
kwargs = (; iip_config = (false, true), expression = Val{true})
116+
117+
build_result = @be build_function(jac, args...; cse = false, kwargs...);
118+
build_times[1][i] = minimum(x -> x.time, build_result.samples)
119+
jacfn_nocse = eval(build_function(jac, args...; cse = false, kwargs...)[2])
120+
121+
build_result = @be build_function(jac, args...; cse = true, kwargs...);
122+
build_times[2][i] = minimum(x -> x.time, build_result.samples)
123+
jacfn_cse = eval(build_function(jac, args...; cse = true, kwargs...)[2])
124+
125+
functions[1][i] = let buffer = similar(jac, Float64), fn = jacfn_nocse
126+
function nocse(u, p, t)
127+
fn(buffer, u, p, t)
128+
buffer
129+
end
130+
end
131+
functions[2][i] = let buffer = similar(jac, Float64), fn = jacfn_cse
132+
function cse(u, p, t)
133+
fn(buffer, u, p, t)
134+
buffer
135+
end
136+
end
137+
138+
return nothing
139+
end
140+
141+
function run_and_time_call!(i, u, p, tt, functions, first_call_times, second_call_times)
142+
jacfn_nocse = functions[1][i]
143+
jacfn_cse = functions[2][i]
144+
145+
call_result = @timed jacfn_nocse(u, p, tt)
146+
first_call_times[1][i] = call_result.time
147+
call_result = @timed jacfn_cse(u, p, tt)
148+
first_call_times[2][i] = call_result.time
149+
150+
call_result = @be jacfn_nocse(u, p, tt)
151+
second_call_times[1][i] = minimum(x -> x.time, call_result.samples)
152+
call_result = @be jacfn_cse(u, p, tt)
153+
second_call_times[2][i] = minimum(x -> x.time, call_result.samples)
154+
end
155+
```
156+
157+
# Run benchmark
158+
159+
```julia
160+
Chairmarks.DEFAULTS.seconds = 15.0
161+
N = [10, 20, 40, 80, 160, 320]
162+
jacobian_times = [zeros(Float64, length(N)), zeros(Float64, length(N))]
163+
jacobian_allocs = copy.(jacobian_times)
164+
functions = [Vector{Any}(undef, length(N)), Vector{Any}(undef, length(N))]
165+
# [without_cse_times, with_cse_times]
166+
build_times = copy.(jacobian_times)
167+
first_call_times = copy.(jacobian_times)
168+
second_call_times = copy.(jacobian_times)
169+
170+
iv = ModelingToolkit.get_iv(osys)
171+
run_and_time_construct!(rhs, vars, pars, iv, 10, 1, jacobian_times, jacobian_allocs, build_times, functions)
172+
run_and_time_call!(1, u, p, tt, functions, first_call_times, second_call_times)
173+
for (i, n) in enumerate(N)
174+
@info i n
175+
run_and_time_construct!(rhs, vars, pars, iv, n, i, jacobian_times, jacobian_allocs, build_times, functions)
176+
end
177+
for (i, n) in enumerate(N)
178+
@info i n
179+
run_and_time_call!(i, u, p, tt, functions, first_call_times, second_call_times)
180+
end
181+
```
182+
183+
# Plot figures
184+
185+
```julia
186+
tabledata = hcat(N, jacobian_times..., jacobian_allocs..., build_times..., first_call_times..., second_call_times...)
187+
header = ["N", "Jacobian time (no hashconsing)", "Jacobian time (hashconsing)", "Jacobian allocated memory (no hashconsing) (B)", "Jacobian allocated memory (hashconsing) (B)", "`build_function` time (no CSE)", "`build_function` time (CSE)", "First call time (no CSE)", "First call time (CSE)", "Second call time (no CSE)", "Second call time (CSE)"]
188+
pretty_table(tabledata; header, backend = Val(:html))
189+
```
190+
191+
```julia
192+
f = Figure(size = (750, 400))
193+
titles = [
194+
"Jacobian symbolic computation", "Jacobian symbolic computation", "Code generation",
195+
"Numerical function compilation", "Numerical function evaluation"]
196+
labels = ["Time (seconds)", "Allocated memory (bytes)",
197+
"Time (seconds)", "Time (seconds)", "Time (seconds)"]
198+
times = [jacobian_times, jacobian_allocs, build_times, first_call_times, second_call_times]
199+
axes = Axis[]
200+
for i in 1:2
201+
label = labels[i]
202+
data = times[i]
203+
ax = Axis(f[1, i], xscale = log10, yscale = log10, xlabel = "model size",
204+
xlabelsize = 10, ylabel = label, ylabelsize = 10, xticks = N,
205+
title = titles[i], titlesize = 12, xticklabelsize = 10, yticklabelsize = 10)
206+
push!(axes, ax)
207+
l1 = scatterlines!(ax, N, data[1], label = "without hashconsing")
208+
l2 = scatterlines!(ax, N, data[2], label = "with hashconsing")
209+
end
210+
Legend(f[1, 3], axes[1], "Methods", tellwidth = false, labelsize = 12, titlesize = 15)
211+
axes2 = Axis[]
212+
# make equal y-axis unit length
213+
mn3, mx3 = extrema(reduce(vcat, times[3]))
214+
xn3 = log10(mx3 / mn3)
215+
mn4, mx4 = extrema(reduce(vcat, times[4]))
216+
xn4 = log10(mx4 / mn4)
217+
mn5, mx5 = extrema(reduce(vcat, times[5]))
218+
xn5 = log10(mx5 / mn5)
219+
xn = max(xn3, xn4, xn5)
220+
xn += 0.2
221+
hxn = xn / 2
222+
hxn3 = (log10(mx3) + log10(mn3)) / 2
223+
hxn4 = (log10(mx4) + log10(mn4)) / 2
224+
hxn5 = (log10(mx5) + log10(mn5)) / 2
225+
ylims = [(exp10(hxn3 - hxn), exp10(hxn3 + hxn)), (exp10(hxn4 - hxn), exp10(hxn4 + hxn)),
226+
(exp10(hxn5 - hxn), exp10(hxn5 + hxn))]
227+
for i in 1:3
228+
ir = i + 2
229+
label = labels[ir]
230+
data = times[ir]
231+
ax = Axis(f[2, i], xscale = log10, yscale = log10, xlabel = "model size",
232+
xlabelsize = 10, ylabel = label, ylabelsize = 10, xticks = N,
233+
title = titles[ir], titlesize = 12, xticklabelsize = 10, yticklabelsize = 10)
234+
ylims!(ax, ylims[i]...)
235+
push!(axes2, ax)
236+
l1 = scatterlines!(ax, N, data[1], label = "without hashconsing")
237+
l2 = scatterlines!(ax, N, data[2], label = "with hashconsing")
238+
end
239+
save("bcr.pdf", f)
240+
f
241+
```

0 commit comments

Comments
 (0)