Optimizing DiffEq Code

Chris Rackauckas

In this notebook we will walk through some of the main tools for optimizing your code in order to efficiently solve DifferentialEquations.jl. User-side optimizations are important because, for sufficiently difficult problems, most of the time will be spent inside of your f function, the function you are trying to solve. "Efficient" integrators are those that reduce the required number of f calls to hit the error tolerance. The main ideas for optimizing your DiffEq code, or any Julia function, are the following:

  • Make it non-allocating

  • Use StaticArrays for small arrays

  • Use broadcast fusion

  • Make it type-stable

  • Reduce redundant calculations

  • Make use of BLAS calls

  • Optimize algorithm choice

We'll discuss these strategies in the context of small and large systems. Let's start with small systems.

Optimizing Small Systems (<100 DEs)

Let's take the classic Lorenz system from before. Let's start by naively writing the system in its out-of-place form:

function lorenz(u,p,t)
 dx = 10.0*(u[2]-u[1])
 dy = u[1]*(28.0-u[3]) - u[2]
 dz = u[1]*u[2] - (8/3)*u[3]
 [dx,dy,dz]
end
lorenz (generic function with 1 method)

Here, lorenz returns an object, [dx,dy,dz], which is created within the body of lorenz.

This is a common code pattern from high-level languages like MATLAB, SciPy, or R's deSolve. However, the issue with this form is that it allocates a vector, [dx,dy,dz], at each step. Let's benchmark the solution process with this choice of function:

using DifferentialEquations, BenchmarkTools
u0 = [1.0;0.0;0.0]
tspan = (0.0,100.0)
prob = ODEProblem(lorenz,u0,tspan)
@benchmark solve(prob,Tsit5())
BenchmarkTools.Trial: 
  memory estimate:  10.63 MiB
  allocs estimate:  99631
  --------------
  minimum time:     3.325 ms (0.00% GC)
  median time:      3.983 ms (0.00% GC)
  mean time:        8.005 ms (53.49% GC)
  maximum time:     16.798 ms (58.35% GC)
  --------------
  samples:          624
  evals/sample:     1

The BenchmarkTools package's @benchmark runs the code multiple times to get an accurate measurement. The minimum time is the time it takes when your OS and other background processes aren't getting in the way. Notice that in this case it takes about 5ms to solve and allocates around 11.11 MiB. However, if we were to use this inside of a real user code we'd see a lot of time spent doing garbage collection (GC) to clean up all of the arrays we made. Even if we turn off saving we have these allocations.

@benchmark solve(prob,Tsit5(),save_everystep=false)
BenchmarkTools.Trial: 
  memory estimate:  9.30 MiB
  allocs estimate:  87129
  --------------
  minimum time:     2.885 ms (0.00% GC)
  median time:      3.157 ms (0.00% GC)
  mean time:        6.832 ms (52.67% GC)
  maximum time:     16.506 ms (58.01% GC)
  --------------
  samples:          731
  evals/sample:     1

The problem of course is that arrays are created every time our derivative function is called. This function is called multiple times per step and is thus the main source of memory usage. To fix this, we can use the in-place form to ***make our code non-allocating***:

function lorenz!(du,u,p,t)
 du[1] = 10.0*(u[2]-u[1])
 du[2] = u[1]*(28.0-u[3]) - u[2]
 du[3] = u[1]*u[2] - (8/3)*u[3]
end
lorenz! (generic function with 1 method)

Here, instead of creating an array each time, we utilized the cache array du. When the inplace form is used, DifferentialEquations.jl takes a different internal route that minimizes the internal allocations as well. When we benchmark this function, we will see quite a difference.

u0 = [1.0;0.0;0.0]
tspan = (0.0,100.0)
prob = ODEProblem(lorenz!,u0,tspan)
@benchmark solve(prob,Tsit5())
BenchmarkTools.Trial: 
  memory estimate:  1.34 MiB
  allocs estimate:  12593
  --------------
  minimum time:     874.414 μs (0.00% GC)
  median time:      902.574 μs (0.00% GC)
  mean time:        1.414 ms (34.51% GC)
  maximum time:     10.179 ms (89.45% GC)
  --------------
  samples:          3526
  evals/sample:     1
@benchmark solve(prob,Tsit5(),save_everystep=false)
BenchmarkTools.Trial: 
  memory estimate:  6.86 KiB
  allocs estimate:  91
  --------------
  minimum time:     477.740 μs (0.00% GC)
  median time:      481.344 μs (0.00% GC)
  mean time:        485.000 μs (0.58% GC)
  maximum time:     9.914 ms (94.58% GC)
  --------------
  samples:          10000
  evals/sample:     1

There is a 4x time difference just from that change! Notice there are still some allocations and this is due to the construction of the integration cache. But this doesn't scale with the problem size:

tspan = (0.0,500.0) # 5x longer than before
prob = ODEProblem(lorenz!,u0,tspan)
@benchmark solve(prob,Tsit5(),save_everystep=false)
BenchmarkTools.Trial: 
  memory estimate:  6.86 KiB
  allocs estimate:  91
  --------------
  minimum time:     2.466 ms (0.00% GC)
  median time:      2.481 ms (0.00% GC)
  mean time:        2.485 ms (0.00% GC)
  maximum time:     2.909 ms (0.00% GC)
  --------------
  samples:          2010
  evals/sample:     1

since that's all just setup allocations.

But if the system is small we can optimize even more.

Allocations are only expensive if they are "heap allocations". For a more in-depth definition of heap allocations, there are a lot of sources online. But a good working definition is that heap allocations are variable-sized slabs of memory which have to be pointed to, and this pointer indirection costs time. Additionally, the heap has to be managed and the garbage controllers has to actively keep track of what's on the heap.

However, there's an alternative to heap allocations, known as stack allocations. The stack is statically-sized (known at compile time) and thus its accesses are quick. Additionally, the exact block of memory is known in advance by the compiler, and thus re-using the memory is cheap. This means that allocating on the stack has essentially no cost!

Arrays have to be heap allocated because their size (and thus the amount of memory they take up) is determined at runtime. But there are structures in Julia which are stack-allocated. structs for example are stack-allocated "value-type"s. Tuples are a stack-allocated collection. The most useful data structure for DiffEq though is the StaticArray from the package StaticArrays.jl. These arrays have their length determined at compile-time. They are created using macros attached to normal array expressions, for example:

using StaticArrays
A = @SVector [2.0,3.0,5.0]
3-element StaticArrays.SArray{Tuple{3},Float64,1,3}:
 2.0
 3.0
 5.0

Notice that the 3 after SVector gives the size of the SVector. It cannot be changed. Additionally, SVectors are immutable, so we have to create a new SVector to change values. But remember, we don't have to worry about allocations because this data structure is stack-allocated. SArrays have a lot of extra optimizations as well: they have fast matrix multiplication, fast QR factorizations, etc. which directly make use of the information about the size of the array. Thus, when possible they should be used.

Unfortunately static arrays can only be used for sufficiently small arrays. After a certain size, they are forced to heap allocate after some instructions and their compile time balloons. Thus static arrays shouldn't be used if your system has more than 100 variables. Additionally, only the native Julia algorithms can fully utilize static arrays.

Let's ***optimize lorenz using static arrays***. Note that in this case, we want to use the out-of-place allocating form, but this time we want to output a static array:

function lorenz_static(u,p,t)
 dx = 10.0*(u[2]-u[1])
 dy = u[1]*(28.0-u[3]) - u[2]
 dz = u[1]*u[2] - (8/3)*u[3]
 @SVector [dx,dy,dz]
end
lorenz_static (generic function with 1 method)

To make the solver internally use static arrays, we simply give it a static array as the initial condition:

u0 = @SVector [1.0,0.0,0.0]
tspan = (0.0,100.0)
prob = ODEProblem(lorenz_static,u0,tspan)
@benchmark solve(prob,Tsit5())
BenchmarkTools.Trial: 
  memory estimate:  461.59 KiB
  allocs estimate:  2583
  --------------
  minimum time:     498.084 μs (0.00% GC)
  median time:      504.044 μs (0.00% GC)
  mean time:        602.322 μs (15.86% GC)
  maximum time:     5.568 ms (89.09% GC)
  --------------
  samples:          8267
  evals/sample:     1
@benchmark solve(prob,Tsit5(),save_everystep=false)
BenchmarkTools.Trial: 
  memory estimate:  6.16 KiB
  allocs estimate:  73
  --------------
  minimum time:     399.249 μs (0.00% GC)
  median time:      404.952 μs (0.00% GC)
  mean time:        407.993 μs (0.43% GC)
  maximum time:     9.236 ms (95.17% GC)
  --------------
  samples:          10000
  evals/sample:     1

And that's pretty much all there is to it. With static arrays you don't have to worry about allocating, so use operations like * and don't worry about fusing operations (discussed in the next section). Do "the vectorized code" of R/MATLAB/Python and your code in this case will be fast, or directly use the numbers/values.

Exercise 1

Implement the out-of-place array, in-place array, and out-of-place static array forms for the Henon-Heiles System and time the results.

Optimizing Large Systems

Interlude: Managing Allocations with Broadcast Fusion

When your system is sufficiently large, or you have to make use of a non-native Julia algorithm, you have to make use of Arrays. In order to use arrays in the most efficient manner, you need to be careful about temporary allocations. Vectorized calculations naturally have plenty of temporary array allocations. This is because a vectorized calculation outputs a vector. Thus:

A = rand(1000,1000); B = rand(1000,1000); C = rand(1000,1000)
test(A,B,C) = A + B + C
@benchmark test(A,B,C)
BenchmarkTools.Trial: 
  memory estimate:  7.63 MiB
  allocs estimate:  3
  --------------
  minimum time:     3.327 ms (0.00% GC)
  median time:      3.390 ms (0.00% GC)
  mean time:        4.734 ms (29.21% GC)
  maximum time:     7.540 ms (56.46% GC)
  --------------
  samples:          1053
  evals/sample:     1

That expression A + B + C creates 2 arrays. It first creates one for the output of A + B, then uses that result array to + C to get the final result. 2 arrays! We don't want that! The first thing to do to fix this is to use broadcast fusion. Broadcast fusion puts expressions together. For example, instead of doing the + operations separately, if we were to add them all at the same time, then we would only have a single array that's created. For example:

test2(A,B,C) = map((a,b,c)->a+b+c,A,B,C)
@benchmark test2(A,B,C)
BenchmarkTools.Trial: 
  memory estimate:  7.63 MiB
  allocs estimate:  5
  --------------
  minimum time:     3.381 ms (0.00% GC)
  median time:      3.427 ms (0.00% GC)
  mean time:        4.754 ms (29.19% GC)
  maximum time:     7.835 ms (58.12% GC)
  --------------
  samples:          1050
  evals/sample:     1

Puts the whole expression into a single function call, and thus only one array is required to store output. This is the same as writing the loop:

function test3(A,B,C)
    D = similar(A)
    @inbounds for i in eachindex(A)
        D[i] = A[i] + B[i] + C[i]
    end
    D
end
@benchmark test3(A,B,C)
BenchmarkTools.Trial: 
  memory estimate:  7.63 MiB
  allocs estimate:  2
  --------------
  minimum time:     3.298 ms (0.00% GC)
  median time:      3.385 ms (0.00% GC)
  mean time:        4.731 ms (29.32% GC)
  maximum time:     7.617 ms (56.09% GC)
  --------------
  samples:          1054
  evals/sample:     1

However, Julia's broadcast is syntactic sugar for this. If multiple expressions have a ., then it will put those vectorized operations together. Thus:

test4(A,B,C) = A .+ B .+ C
@benchmark test4(A,B,C)
BenchmarkTools.Trial: 
  memory estimate:  7.63 MiB
  allocs estimate:  2
  --------------
  minimum time:     3.298 ms (0.00% GC)
  median time:      3.381 ms (0.00% GC)
  mean time:        4.731 ms (29.34% GC)
  maximum time:     7.653 ms (55.54% GC)
  --------------
  samples:          1054
  evals/sample:     1

is a version with only 1 array created (the output). Note that .s can be used with function calls as well:

sin.(A) .+ sin.(B)
1000×1000 Array{Float64,2}:
 0.883491  1.46941   1.31592   1.26031   …  0.278758   0.854864  0.546223
 0.441868  1.3608    1.0348    0.93745      0.963036   0.644012  0.79876 
 0.6326    1.3906    0.92889   0.95655      1.30921    1.60003   0.603888
 0.763656  0.89846   0.633754  0.751764     0.78233    0.976681  0.998619
 0.53243   1.14022   0.974028  0.331864     1.24247    1.31588   0.504788
 0.801383  0.463875  0.880555  0.717018  …  0.945005   1.13482   0.724712
 0.700827  1.6268    1.20322   0.770063     0.898224   1.00471   1.3217  
 1.02611   1.29816   0.898958  0.462772     0.0723639  1.14904   1.07084 
 1.17631   0.919471  1.498     0.919398     0.77762    0.499798  1.22517 
 1.15896   1.45735   1.36369   1.47575      1.15369    0.849005  0.949496
 ⋮                                       ⋱                               
 1.30411   1.19439   0.664451  0.543255     0.719017   0.475454  0.61456 
 1.27597   0.599419  0.303519  1.01371      1.25969    0.496569  1.33776 
 1.37308   0.150727  0.878723  1.15994      1.32613    0.707558  0.426822
 0.638035  0.735053  0.840038  0.814784     0.799952   1.19904   0.51726 
 1.24534   0.64412   1.27147   1.13847   …  1.47634    0.548757  1.53114 
 1.07713   0.951169  1.326     1.19763      0.998489   1.10742   0.669014
 0.603015  1.55926   0.976357  0.847166     1.09955    0.229294  1.08651 
 1.29188   1.14394   0.55431   0.88283      0.791295   0.914774  0.63556 
 1.05681   0.870731  0.820364  0.72312      1.39749    0.793006  1.0619

Also, the @. macro applys a dot to every operator:

test5(A,B,C) = @. A + B + C #only one array allocated
@benchmark test5(A,B,C)
BenchmarkTools.Trial: 
  memory estimate:  7.63 MiB
  allocs estimate:  3
  --------------
  minimum time:     3.350 ms (0.00% GC)
  median time:      3.395 ms (0.00% GC)
  mean time:        4.749 ms (29.40% GC)
  maximum time:     7.595 ms (56.66% GC)
  --------------
  samples:          1050
  evals/sample:     1

Using these tools we can get rid of our intermediate array allocations for many vectorized function calls. But we are still allocating the output array. To get rid of that allocation, we can instead use mutation. Mutating broadcast is done via .=. For example, if we pre-allocate the output:

D = zeros(1000,1000);

Then we can keep re-using this cache for subsequent calculations. The mutating broadcasting form is:

test6!(D,A,B,C) = D .= A .+ B .+ C #only one array allocated
@benchmark test6!(D,A,B,C)
BenchmarkTools.Trial: 
  memory estimate:  0 bytes
  allocs estimate:  0
  --------------
  minimum time:     3.274 ms (0.00% GC)
  median time:      3.294 ms (0.00% GC)
  mean time:        3.301 ms (0.00% GC)
  maximum time:     3.762 ms (0.00% GC)
  --------------
  samples:          1509
  evals/sample:     1

If we use @. before the =, then it will turn it into .=:

test7!(D,A,B,C) = @. D = A + B + C #only one array allocated
@benchmark test7!(D,A,B,C)
BenchmarkTools.Trial: 
  memory estimate:  0 bytes
  allocs estimate:  0
  --------------
  minimum time:     3.278 ms (0.00% GC)
  median time:      3.299 ms (0.00% GC)
  mean time:        3.307 ms (0.00% GC)
  maximum time:     3.622 ms (0.00% GC)
  --------------
  samples:          1506
  evals/sample:     1

Notice that in this case, there is no "output", and instead the values inside of D are what are changed (like with the DiffEq inplace function). Many Julia functions have a mutating form which is denoted with a !. For example, the mutating form of the map is map!:

test8!(D,A,B,C) = map!((a,b,c)->a+b+c,D,A,B,C)
@benchmark test8!(D,A,B,C)
BenchmarkTools.Trial: 
  memory estimate:  32 bytes
  allocs estimate:  1
  --------------
  minimum time:     3.286 ms (0.00% GC)
  median time:      3.387 ms (0.00% GC)
  mean time:        3.391 ms (0.00% GC)
  maximum time:     3.626 ms (0.00% GC)
  --------------
  samples:          1469
  evals/sample:     1

Some operations require using an alternate mutating form in order to be fast. For example, matrix multiplication via * allocates a temporary:

@benchmark A*B
BenchmarkTools.Trial: 
  memory estimate:  7.63 MiB
  allocs estimate:  2
  --------------
  minimum time:     21.865 ms (0.00% GC)
  median time:      21.989 ms (0.00% GC)
  mean time:        23.508 ms (6.28% GC)
  maximum time:     29.563 ms (15.02% GC)
  --------------
  samples:          213
  evals/sample:     1

Instead, we can use the mutating form mul! into a cache array to avoid allocating the output:

using LinearAlgebra
@benchmark mul!(D,A,B) # same as D = A * B
BenchmarkTools.Trial: 
  memory estimate:  0 bytes
  allocs estimate:  0
  --------------
  minimum time:     21.267 ms (0.00% GC)
  median time:      21.347 ms (0.00% GC)
  mean time:        21.362 ms (0.00% GC)
  maximum time:     23.280 ms (0.00% GC)
  --------------
  samples:          234
  evals/sample:     1

For repeated calculations this reduced allocation can stop GC cycles and thus lead to more efficient code. Additionally, ***we can fuse together higher level linear algebra operations using BLAS***. The package SugarBLAS.jl makes it easy to write higher level operations like alpha*B*A + beta*C as mutating BLAS calls.

Example Optimization: Gierer-Meinhardt Reaction-Diffusion PDE Discretization

Let's optimize the solution of a Reaction-Diffusion PDE's discretization. In its discretized form, this is the ODE:

\[ \begin{align} du &= D_1 (A_y u + u A_x) + \frac{au^2}{v} + \bar{u} - \alpha u\\ dv &= D_2 (A_y v + v A_x) + a u^2 + \beta v \end{align} \]

where $u$, $v$, and $A$ are matrices. Here, we will use the simplified version where $A$ is the tridiagonal stencil $[1,-2,1]$, i.e. it's the 2D discretization of the LaPlacian. The native code would be something along the lines of:

# Generate the constants
p = (1.0,1.0,1.0,10.0,0.001,100.0) # a,α,ubar,β,D1,D2
N = 100
Ax = Array(Tridiagonal([1.0 for i in 1:N-1],[-2.0 for i in 1:N],[1.0 for i in 1:N-1]))
Ay = copy(Ax)
Ax[2,1] = 2.0
Ax[end-1,end] = 2.0
Ay[1,2] = 2.0
Ay[end,end-1] = 2.0

function basic_version!(dr,r,p,t)
  a,α,ubar,β,D1,D2 = p
  u = r[:,:,1]
  v = r[:,:,2]
  Du = D1*(Ay*u + u*Ax)
  Dv = D2*(Ay*v + v*Ax)
  dr[:,:,1] = Du .+ a.*u.*u./v .+ ubar .- α*u
  dr[:,:,2] = Dv .+ a.*u.*u .- β*v
end

a,α,ubar,β,D1,D2 = p
uss = (ubar+β)/α
vss = (a/β)*uss^2
r0 = zeros(100,100,2)
r0[:,:,1] .= uss.+0.1.*rand.()
r0[:,:,2] .= vss

prob = ODEProblem(basic_version!,r0,(0.0,0.1),p)
ODEProblem with uType Array{Float64,3} and tType Float64. In-place: true
timespan: (0.0, 0.1)
u0: [11.0394 11.0028 … 11.0105 11.0494; 11.0489 11.0042 … 11.0641 11.0566; 
… ; 11.0191 11.0122 … 11.0218 11.0182; 11.0192 11.0377 … 11.0569 11.0811]

[12.1 12.1 … 12.1 12.1; 12.1 12.1 … 12.1 12.1; … ; 12.1 12.1 … 12.1 12.1; 1
2.1 12.1 … 12.1 12.1]

In this version we have encoded our initial condition to be a 3-dimensional array, with u[:,:,1] being the A part and u[:,:,2] being the B part.

@benchmark solve(prob,Tsit5())
BenchmarkTools.Trial: 
  memory estimate:  186.88 MiB
  allocs estimate:  8589
  --------------
  minimum time:     124.746 ms (31.02% GC)
  median time:      258.186 ms (66.57% GC)
  mean time:        227.409 ms (62.05% GC)
  maximum time:     313.404 ms (72.46% GC)
  --------------
  samples:          22
  evals/sample:     1

While this version isn't very efficient,

We recommend writing the "high-level" code first, and iteratively optimizing it!

The first thing that we can do is get rid of the slicing allocations. The operation r[:,:,1] creates a temporary array instead of a "view", i.e. a pointer to the already existing memory. To make it a view, add @view. Note that we have to be careful with views because they point to the same memory, and thus changing a view changes the original values:

A = rand(4)
@show A
A = [0.953358, 0.408393, 0.0122052, 0.277688]
B = @view A[1:3]
B[2] = 2
@show A
A = [0.953358, 2.0, 0.0122052, 0.277688]
4-element Array{Float64,1}:
 0.9533580491706126  
 2.0                 
 0.012205168934875665
 0.2776877187822635

Notice that changing B changed A. This is something to be careful of, but at the same time we want to use this since we want to modify the output dr. Additionally, the last statement is a purely element-wise operation, and thus we can make use of broadcast fusion there. Let's rewrite basic_version! to ***avoid slicing allocations*** and to ***use broadcast fusion***:

function gm2!(dr,r,p,t)
  a,α,ubar,β,D1,D2 = p
  u = @view r[:,:,1]
  v = @view r[:,:,2]
  du = @view dr[:,:,1]
  dv = @view dr[:,:,2]
  Du = D1*(Ay*u + u*Ax)
  Dv = D2*(Ay*v + v*Ax)
  @. du = Du + a.*u.*u./v + ubar - α*u
  @. dv = Dv + a.*u.*u - β*v
end
prob = ODEProblem(gm2!,r0,(0.0,0.1),p)
@benchmark solve(prob,Tsit5())
BenchmarkTools.Trial: 
  memory estimate:  119.55 MiB
  allocs estimate:  7119
  --------------
  minimum time:     98.919 ms (24.81% GC)
  median time:      168.334 ms (48.75% GC)
  mean time:        195.261 ms (56.70% GC)
  maximum time:     299.465 ms (71.34% GC)
  --------------
  samples:          26
  evals/sample:     1

Now, most of the allocations are taking place in Du = D1*(Ay*u + u*Ax) since those operations are vectorized and not mutating. We should instead replace the matrix multiplications with mul!. When doing so, we will need to have cache variables to write into. This looks like:

Ayu = zeros(N,N)
uAx = zeros(N,N)
Du = zeros(N,N)
Ayv = zeros(N,N)
vAx = zeros(N,N)
Dv = zeros(N,N)
function gm3!(dr,r,p,t)
  a,α,ubar,β,D1,D2 = p
  u = @view r[:,:,1]
  v = @view r[:,:,2]
  du = @view dr[:,:,1]
  dv = @view dr[:,:,2]
  mul!(Ayu,Ay,u)
  mul!(uAx,u,Ax)
  mul!(Ayv,Ay,v)
  mul!(vAx,v,Ax)
  @. Du = D1*(Ayu + uAx)
  @. Dv = D2*(Ayv + vAx)
  @. du = Du + a*u*u./v + ubar - α*u
  @. dv = Dv + a*u*u - β*v
end
prob = ODEProblem(gm3!,r0,(0.0,0.1),p)
@benchmark solve(prob,Tsit5())
BenchmarkTools.Trial: 
  memory estimate:  29.76 MiB
  allocs estimate:  5355
  --------------
  minimum time:     69.427 ms (6.83% GC)
  median time:      69.752 ms (6.91% GC)
  mean time:        71.224 ms (8.97% GC)
  maximum time:     74.430 ms (12.65% GC)
  --------------
  samples:          71
  evals/sample:     1

But our temporary variables are global variables. We need to either declare the caches as const or localize them. We can localize them by adding them to the parameters, p. It's easier for the compiler to reason about local variables than global variables. ***Localizing variables helps to ensure type stability***.

p = (1.0,1.0,1.0,10.0,0.001,100.0,Ayu,uAx,Du,Ayv,vAx,Dv) # a,α,ubar,β,D1,D2
function gm4!(dr,r,p,t)
  a,α,ubar,β,D1,D2,Ayu,uAx,Du,Ayv,vAx,Dv = p
  u = @view r[:,:,1]
  v = @view r[:,:,2]
  du = @view dr[:,:,1]
  dv = @view dr[:,:,2]
  mul!(Ayu,Ay,u)
  mul!(uAx,u,Ax)
  mul!(Ayv,Ay,v)
  mul!(vAx,v,Ax)
  @. Du = D1*(Ayu + uAx)
  @. Dv = D2*(Ayv + vAx)
  @. du = Du + a*u*u./v + ubar - α*u
  @. dv = Dv + a*u*u - β*v
end
prob = ODEProblem(gm4!,r0,(0.0,0.1),p)
@benchmark solve(prob,Tsit5())
BenchmarkTools.Trial: 
  memory estimate:  29.66 MiB
  allocs estimate:  1090
  --------------
  minimum time:     55.101 ms (8.56% GC)
  median time:      55.416 ms (8.62% GC)
  mean time:        56.881 ms (11.06% GC)
  maximum time:     60.062 ms (15.39% GC)
  --------------
  samples:          88
  evals/sample:     1

We could then use the BLAS gemmv to optimize the matrix multiplications some more, but instead let's devectorize the stencil.

p = (1.0,1.0,1.0,10.0,0.001,100.0,N)
function fast_gm!(du,u,p,t)
  a,α,ubar,β,D1,D2,N = p

  @inbounds for j in 2:N-1, i in 2:N-1
    du[i,j,1] = D1*(u[i-1,j,1] + u[i+1,j,1] + u[i,j+1,1] + u[i,j-1,1] - 4u[i,j,1]) +
              a*u[i,j,1]^2/u[i,j,2] + ubar - α*u[i,j,1]
  end

  @inbounds for j in 2:N-1, i in 2:N-1
    du[i,j,2] = D2*(u[i-1,j,2] + u[i+1,j,2] + u[i,j+1,2] + u[i,j-1,2] - 4u[i,j,2]) +
            a*u[i,j,1]^2 - β*u[i,j,2]
  end

  @inbounds for j in 2:N-1
    i = 1
    du[1,j,1] = D1*(2u[i+1,j,1] + u[i,j+1,1] + u[i,j-1,1] - 4u[i,j,1]) +
            a*u[i,j,1]^2/u[i,j,2] + ubar - α*u[i,j,1]
  end
  @inbounds for j in 2:N-1
    i = 1
    du[1,j,2] = D2*(2u[i+1,j,2] + u[i,j+1,2] + u[i,j-1,2] - 4u[i,j,2]) +
            a*u[i,j,1]^2 - β*u[i,j,2]
  end
  @inbounds for j in 2:N-1
    i = N
    du[end,j,1] = D1*(2u[i-1,j,1] + u[i,j+1,1] + u[i,j-1,1] - 4u[i,j,1]) +
           a*u[i,j,1]^2/u[i,j,2] + ubar - α*u[i,j,1]
  end
  @inbounds for j in 2:N-1
    i = N
    du[end,j,2] = D2*(2u[i-1,j,2] + u[i,j+1,2] + u[i,j-1,2] - 4u[i,j,2]) +
           a*u[i,j,1]^2 - β*u[i,j,2]
  end

  @inbounds for i in 2:N-1
    j = 1
    du[i,1,1] = D1*(u[i-1,j,1] + u[i+1,j,1] + 2u[i,j+1,1] - 4u[i,j,1]) +
              a*u[i,j,1]^2/u[i,j,2] + ubar - α*u[i,j,1]
  end
  @inbounds for i in 2:N-1
    j = 1
    du[i,1,2] = D2*(u[i-1,j,2] + u[i+1,j,2] + 2u[i,j+1,2] - 4u[i,j,2]) +
              a*u[i,j,1]^2 - β*u[i,j,2]
  end
  @inbounds for i in 2:N-1
    j = N
    du[i,end,1] = D1*(u[i-1,j,1] + u[i+1,j,1] + 2u[i,j-1,1] - 4u[i,j,1]) +
             a*u[i,j,1]^2/u[i,j,2] + ubar - α*u[i,j,1]
  end
  @inbounds for i in 2:N-1
    j = N
    du[i,end,2] = D2*(u[i-1,j,2] + u[i+1,j,2] + 2u[i,j-1,2] - 4u[i,j,2]) +
             a*u[i,j,1]^2 - β*u[i,j,2]
  end

  @inbounds begin
    i = 1; j = 1
    du[1,1,1] = D1*(2u[i+1,j,1] + 2u[i,j+1,1] - 4u[i,j,1]) +
              a*u[i,j,1]^2/u[i,j,2] + ubar - α*u[i,j,1]
    du[1,1,2] = D2*(2u[i+1,j,2] + 2u[i,j+1,2] - 4u[i,j,2]) +
              a*u[i,j,1]^2 - β*u[i,j,2]

    i = 1; j = N
    du[1,N,1] = D1*(2u[i+1,j,1] + 2u[i,j-1,1] - 4u[i,j,1]) +
             a*u[i,j,1]^2/u[i,j,2] + ubar - α*u[i,j,1]
    du[1,N,2] = D2*(2u[i+1,j,2] + 2u[i,j-1,2] - 4u[i,j,2]) +
             a*u[i,j,1]^2 - β*u[i,j,2]

    i = N; j = 1
    du[N,1,1] = D1*(2u[i-1,j,1] + 2u[i,j+1,1] - 4u[i,j,1]) +
             a*u[i,j,1]^2/u[i,j,2] + ubar - α*u[i,j,1]
    du[N,1,2] = D2*(2u[i-1,j,2] + 2u[i,j+1,2] - 4u[i,j,2]) +
             a*u[i,j,1]^2 - β*u[i,j,2]

    i = N; j = N
    du[end,end,1] = D1*(2u[i-1,j,1] + 2u[i,j-1,1] - 4u[i,j,1]) +
             a*u[i,j,1]^2/u[i,j,2] + ubar - α*u[i,j,1]
    du[end,end,2] = D2*(2u[i-1,j,2] + 2u[i,j-1,2] - 4u[i,j,2]) +
             a*u[i,j,1]^2 - β*u[i,j,2]
   end
end
prob = ODEProblem(fast_gm!,r0,(0.0,0.1),p)
@benchmark solve(prob,Tsit5())
BenchmarkTools.Trial: 
  memory estimate:  29.63 MiB
  allocs estimate:  504
  --------------
  minimum time:     17.318 ms (25.55% GC)
  median time:      17.547 ms (25.42% GC)
  mean time:        18.928 ms (31.18% GC)
  maximum time:     22.233 ms (39.53% GC)
  --------------
  samples:          265
  evals/sample:     1

Lastly, we can do other things like multithread the main loops, but these optimizations get the last 2x-3x out. The main optimizations which apply everywhere are the ones we just performed (though the last one only works if your matrix is a stencil. This is known as a matrix-free implementation of the PDE discretization).

This gets us to about 8x faster than our original MATLAB/SciPy/R vectorized style code!

The last thing to do is then ***optimize our algorithm choice***. We have been using Tsit5() as our test algorithm, but in reality this problem is a stiff PDE discretization and thus one recommendation is to use CVODE_BDF(). However, instead of using the default dense Jacobian, we should make use of the sparse Jacobian afforded by the problem. The Jacobian is the matrix $\frac{df_i}{dr_j}$, where $r$ is read by the linear index (i.e. down columns). But since the $u$ variables depend on the $v$, the band size here is large, and thus this will not do well with a Banded Jacobian solver. Instead, we utilize sparse Jacobian algorithms. CVODE_BDF allows us to use a sparse Newton-Krylov solver by setting linear_solver = :GMRES (see the solver documentation, and thus we can solve this problem efficiently. Let's see how this scales as we increase the integration time.

prob = ODEProblem(fast_gm!,r0,(0.0,10.0),p)
@benchmark solve(prob,Tsit5())
BenchmarkTools.Trial: 
  memory estimate:  2.76 GiB
  allocs estimate:  41670
  --------------
  minimum time:     2.705 s (44.96% GC)
  median time:      19.769 s (61.02% GC)
  mean time:        19.769 s (61.02% GC)
  maximum time:     36.833 s (62.20% GC)
  --------------
  samples:          2
  evals/sample:     1
using Sundials
@benchmark solve(prob,CVODE_BDF(linear_solver=:GMRES))
BenchmarkTools.Trial: 
  memory estimate:  117.33 MiB
  allocs estimate:  30688
  --------------
  minimum time:     715.584 ms (3.70% GC)
  median time:      897.363 ms (23.22% GC)
  mean time:        875.513 ms (21.20% GC)
  maximum time:     956.867 ms (27.77% GC)
  --------------
  samples:          6
  evals/sample:     1
prob = ODEProblem(fast_gm!,r0,(0.0,100.0),p)
# Will go out of memory if we don't turn off `save_everystep`!
@benchmark solve(prob,Tsit5(),save_everystep=false)
BenchmarkTools.Trial: 
  memory estimate:  2.91 MiB
  allocs estimate:  112
  --------------
  minimum time:     9.372 s (0.00% GC)
  median time:      9.372 s (0.00% GC)
  mean time:        9.372 s (0.00% GC)
  maximum time:     9.372 s (0.00% GC)
  --------------
  samples:          1
  evals/sample:     1
@benchmark solve(prob,CVODE_BDF(linear_solver=:GMRES))
BenchmarkTools.Trial: 
  memory estimate:  323.53 MiB
  allocs estimate:  84834
  --------------
  minimum time:     2.038 s (0.00% GC)
  median time:      2.232 s (8.59% GC)
  mean time:        2.260 s (9.74% GC)
  maximum time:     2.509 s (18.68% GC)
  --------------
  samples:          3
  evals/sample:     1

Now let's check the allocation growth.

@benchmark solve(prob,CVODE_BDF(linear_solver=:GMRES),save_everystep=false)
BenchmarkTools.Trial: 
  memory estimate:  4.36 MiB
  allocs estimate:  75057
  --------------
  minimum time:     1.988 s (0.00% GC)
  median time:      1.989 s (0.00% GC)
  mean time:        1.989 s (0.00% GC)
  maximum time:     1.990 s (0.00% GC)
  --------------
  samples:          3
  evals/sample:     1
prob = ODEProblem(fast_gm!,r0,(0.0,500.0),p)
@benchmark solve(prob,CVODE_BDF(linear_solver=:GMRES),save_everystep=false)
BenchmarkTools.Trial: 
  memory estimate:  5.99 MiB
  allocs estimate:  110232
  --------------
  minimum time:     2.918 s (0.00% GC)
  median time:      2.920 s (0.00% GC)
  mean time:        2.920 s (0.00% GC)
  maximum time:     2.922 s (0.00% GC)
  --------------
  samples:          2
  evals/sample:     1

Notice that we've elimated almost all allocations, allowing the code to grow without hitting garbage collection and slowing down.

Why is CVODE_BDF doing well? What's happening is that, because the problem is stiff, the number of steps required by the explicit Runge-Kutta method grows rapidly, whereas CVODE_BDF is taking large steps. Additionally, the GMRES linear solver form is quite an efficient way to solve the implicit system in this case. This is problem-dependent, and in many cases using a Krylov method effectively requires a preconditioner, so you need to play around with testing other algorithms and linear solvers to find out what works best with your problem.

Conclusion

Julia gives you the tools to optimize the solver "all the way", but you need to make use of it. The main thing to avoid is temporary allocations. For small systems, this is effectively done via static arrays. For large systems, this is done via in-place operations and cache arrays. Either way, the resulting solution can be immensely sped up over vectorized formulations by using these principles.

Appendix

This tutorial is part of the DiffEqTutorials.jl repository, found at: https://github.com/JuliaDiffEq/DiffEqTutorials.jl

To locally run this tutorial, do the following commands:

using DiffEqTutorials
DiffEqTutorials.weave_file("introduction","03-optimizing_diffeq_code.jmd")

Computer Information:

Julia Version 1.1.1
Commit 55e36cc308 (2019-05-16 04:10 UTC)
Platform Info:
  OS: Linux (x86_64-pc-linux-gnu)
  CPU: Intel(R) Core(TM) i7-3770 CPU @ 3.40GHz
  WORD_SIZE: 64
  LIBM: libopenlibm
  LLVM: libLLVM-6.0.1 (ORCJIT, ivybridge)

Package Information:

Status `~/.julia/environments/v1.1/Project.toml`
[7e558dbc-694d-5a72-987c-6f4ebed21442] ArbNumerics 0.5.4
[6e4b80f9-dd63-53aa-95a3-0cdb28fa8baf] BenchmarkTools 0.4.2
[be33ccc6-a3ff-5ff2-a52e-74243cff1e17] CUDAnative 2.2.0
[3a865a2d-5b23-5a0f-bc46-62713ec82fae] CuArrays 1.0.2
[55939f99-70c6-5e9b-8bb0-5071ed7d61fd] DecFP 0.4.8
[abce61dc-4473-55a0-ba07-351d65e31d42] Decimals 0.4.0
[ebbdde9d-f333-5424-9be2-dbf1e9acfb5e] DiffEqBayes 1.1.0
[eb300fae-53e8-50a0-950c-e21f52c2b7e0] DiffEqBiological 3.8.2
[459566f4-90b8-5000-8ac3-15dfb0a30def] DiffEqCallbacks 2.5.2
[f3b72e0c-5b89-59e1-b016-84e28bfd966d] DiffEqDevTools 2.9.0
[1130ab10-4a5a-5621-a13d-e4788d82bd4c] DiffEqParamEstim 1.6.0
[055956cb-9e8b-5191-98cc-73ae4a59e68a] DiffEqPhysics 3.1.0
[6d1b261a-3be8-11e9-3f2f-0b112a9a8436] DiffEqTutorials 0.1.0
[0c46a032-eb83-5123-abaf-570d42b7fbaa] DifferentialEquations 6.4.0
[31c24e10-a181-5473-b8eb-7969acd0382f] Distributions 0.20.0
[497a8b3b-efae-58df-a0af-a86822472b78] DoubleFloats 0.9.1
[f6369f11-7733-5829-9624-2563aa707210] ForwardDiff 0.10.3
[c91e804a-d5a3-530f-b6f0-dfbca275c004] Gadfly 1.0.1
[7073ff75-c697-5162-941a-fcdaad2a7d2a] IJulia 1.18.1
[4138dd39-2aa7-5051-a626-17a0bb65d9c8] JLD 0.9.1
[23fbe1c1-3f47-55db-b15f-69d7ec21a316] Latexify 0.8.2
[eff96d63-e80a-5855-80a2-b1b0885c5ab7] Measurements 2.0.0
[961ee093-0014-501f-94e3-6117800e7a78] ModelingToolkit 0.2.0
[76087f3c-5699-56af-9a33-bf431cd00edd] NLopt 0.5.1
[2774e3e8-f4cf-5e23-947b-6d7e65073b56] NLsolve 4.0.0
[429524aa-4258-5aef-a3af-852621145aeb] Optim 0.18.1
[1dea7af3-3e70-54e6-95c3-0bf5283fa5ed] OrdinaryDiffEq 5.8.1
[65888b18-ceab-5e60-b2b9-181511a3b968] ParameterizedFunctions 4.1.1
[91a5bcdd-55d7-5caf-9e0b-520d859cae80] Plots 0.25.1
[d330b81b-6aea-500a-939a-2ce795aea3ee] PyPlot 2.8.1
[731186ca-8d62-57ce-b412-fbd966d074cd] RecursiveArrayTools 0.20.0
[90137ffa-7385-5640-81b9-e52037218182] StaticArrays 0.11.0
[f3b207a7-027a-5e70-b257-86293d7955fd] StatsPlots 0.11.0
[c3572dad-4567-51f8-b174-8c6c989267f4] Sundials 3.6.1
[1986cc42-f94f-5a68-af5c-568840ba703d] Unitful 0.15.0
[44d3d7a6-8a23-5bf8-98c5-b353f8df5ec9] Weave 0.9.0
[b77e0a4c-d291-57a0-90e8-8db25a27a240] InteractiveUtils
[37e2e46d-f89d-539d-b4ee-838fcccc9c8e] LinearAlgebra
[44cfe95a-1eb2-52ea-b672-e2afdf69b78f] Pkg