Basic tutorial
We present the main features of DifferentiationInterface.jl.
using DifferentiationInterface
Computing a gradient
A common use case of automatic differentiation (AD) is optimizing real-valued functions with first- or second-order methods. Let's define a simple objective (the squared norm) and a random input vector
f(x) = sum(abs2, x)
f (generic function with 1 method)
x = collect(1.0:5.0)
5-element Vector{Float64}:
1.0
2.0
3.0
4.0
5.0
To compute its gradient, we need to choose a "backend", i.e. an AD package to call under the hood. Most backend types are defined by ADTypes.jl and re-exported by DifferentiationInterface.jl.
ForwardDiff.jl is very generic and efficient for low-dimensional inputs, so it's a good starting point:
import ForwardDiff
backend = AutoForwardDiff()
AutoForwardDiff()
To avoid name conflicts, load AD packages with import
instead of using
. Indeed, most AD packages also export operators like gradient
and jacobian
, but you only want to use the ones from DifferentiationInterface.jl.
Now you can use the following syntax to compute the gradient:
gradient(f, backend, x)
5-element Vector{Float64}:
2.0
4.0
6.0
8.0
10.0
Was that fast? BenchmarkTools.jl helps you answer that question.
using BenchmarkTools
@benchmark gradient($f, $backend, $x)
BenchmarkTools.Trial: 10000 samples with 207 evaluations.
Range (min … max): 368.420 ns … 127.788 μs ┊ GC (min … max): 0.00% … 99.46%
Time (median): 568.341 ns ┊ GC (median): 0.00%
Time (mean ± σ): 557.542 ns ± 2.323 μs ┊ GC (mean ± σ): 10.14% ± 3.08%
▁▆█▇▄▂ ▁ ▂▄▅▆▆▆▅▅▆▅▅▄▂▁▁ ▂
████████▆▅▄▇███▇█▆▆▆▅▅▄▃▃▃▃▁▃▃▁▃▁▁▁▃▃▃▁▁▇█████████████████▇██ █
368 ns Histogram: log(frequency) by time 644 ns <
Memory estimate: 624 bytes, allocs estimate: 5.
Not bad, but you can do better.
Overwriting a gradient
Since you know how much space your gradient will occupy (the same as your input x
), you can pre-allocate that memory and offer it to AD. Some backends get a speed boost from this trick.
grad = similar(x)
gradient!(f, grad, backend, x)
grad # has been mutated
5-element Vector{Float64}:
2.0
4.0
6.0
8.0
10.0
The bang indicates that one of the arguments of gradient!
might be mutated. More precisely, our convention is that every positional argument between the function and the backend is mutated.
@benchmark gradient!($f, $grad, $backend, $x)
BenchmarkTools.Trial: 10000 samples with 208 evaluations.
Range (min … max): 357.452 ns … 112.649 μs ┊ GC (min … max): 0.00% … 99.41%
Time (median): 553.082 ns ┊ GC (median): 0.00%
Time (mean ± σ): 523.726 ns ± 1.931 μs ┊ GC (mean ± σ): 8.61% ± 3.17%
█ ▁▁
▅█▇▂▂▂▂▂▂▁▂▂▂▂▂▂▂▂▂▂▂▁▂▁▁▂▂▂▂▂▂▂▂▂▂▁▁▂▁▂▁▁▁▂▂▄▆███▅▃▂▂▂▂▂▂▂▂▂ ▃
357 ns Histogram: frequency by time 615 ns <
Memory estimate: 528 bytes, allocs estimate: 3.
For some reason the in-place version is not much better than your first attempt. However, it makes fewer allocations, thanks to the gradient vector you provided. Don't worry, you can get even more performance.
Preparing for multiple gradients
Internally, ForwardDiff.jl creates some data structures to keep track of things. These objects can be reused between gradient computations, even on different input values. We abstract away the preparation step behind a backend-agnostic syntax:
prep = prepare_gradient(f, backend, zero(x))
DifferentiationInterfaceForwardDiffExt.ForwardDiffGradientPrep{ForwardDiff.GradientConfig{ForwardDiff.Tag{typeof(Main.f), Float64}, Float64, 5, Vector{ForwardDiff.Dual{ForwardDiff.Tag{typeof(Main.f), Float64}, Float64, 5}}}}(ForwardDiff.GradientConfig{ForwardDiff.Tag{typeof(Main.f), Float64}, Float64, 5, Vector{ForwardDiff.Dual{ForwardDiff.Tag{typeof(Main.f), Float64}, Float64, 5}}}((Partials(1.0, 0.0, 0.0, 0.0, 0.0), Partials(0.0, 1.0, 0.0, 0.0, 0.0), Partials(0.0, 0.0, 1.0, 0.0, 0.0), Partials(0.0, 0.0, 0.0, 1.0, 0.0), Partials(0.0, 0.0, 0.0, 0.0, 1.0)), ForwardDiff.Dual{ForwardDiff.Tag{typeof(Main.f), Float64}, Float64, 5}[Dual{ForwardDiff.Tag{typeof(Main.f), Float64}}(0.0,0.0,0.0,0.0,1.0,0.0), Dual{ForwardDiff.Tag{typeof(Main.f), Float64}}(0.0,0.0,0.0,0.0,1.0,0.0), Dual{ForwardDiff.Tag{typeof(Main.f), Float64}}(0.0,0.0,0.0,0.0,1.0,0.0), Dual{ForwardDiff.Tag{typeof(Main.f), Float64}}(0.0,0.0,0.0,0.0,1.0,6.92716478046244e-310), Dual{ForwardDiff.Tag{typeof(Main.f), Float64}}(5.0,6.9271783602733e-310,1.0,0.0,0.0,0.0)]))
You don't need to know what this object is, you just need to pass it to the gradient operator. Note that preparation does not depend on the actual components of the vector x
, just on its type and size. You can thus reuse the prep
for different values of the input.
grad = similar(x)
gradient!(f, grad, prep, backend, x)
grad # has been mutated
5-element Vector{Float64}:
2.0
4.0
6.0
8.0
10.0
Preparation makes the gradient computation much faster, and (in this case) allocation-free.
@benchmark gradient!($f, $grad, $prep, $backend, $x)
BenchmarkTools.Trial: 10000 samples with 995 evaluations.
Range (min … max): 28.224 ns … 52.380 ns ┊ GC (min … max): 0.00% … 0.00%
Time (median): 28.869 ns ┊ GC (median): 0.00%
Time (mean ± σ): 28.978 ns ± 1.011 ns ┊ GC (mean ± σ): 0.00% ± 0.00%
▃▇▇█▅
▂▃▄███████▆▃▃▂▂▂▂▂▂▂▂▁▁▂▂▂▂▂▂▂▂▂▂▂▂▂▁▂▁▂▂▂▁▁▁▁▁▁▂▁▂▁▁▁▁▁▁▁▂ ▃
28.2 ns Histogram: frequency by time 33.9 ns <
Memory estimate: 0 bytes, allocs estimate: 0.
Beware that the prep
object is nearly always mutated by differentiation operators, even though it is given as the last positional argument.
Switching backends
The whole point of DifferentiationInterface.jl is that you can easily experiment with different AD solutions. Typically, for gradients, reverse mode AD might be a better fit, so let's try Zygote.jl!
import Zygote
backend2 = AutoZygote()
AutoZygote()
Once the backend is created, things run smoothly with exactly the same syntax as before:
gradient(f, backend2, x)
5-element Vector{Float64}:
2.0
4.0
6.0
8.0
10.0
And you can run the same benchmarks to see what you gained (although such a small input may not be realistic):
prep2 = prepare_gradient(f, backend2, zero(x))
@benchmark gradient!($f, $grad, $prep2, $backend2, $x)
BenchmarkTools.Trial: 10000 samples with 993 evaluations.
Range (min … max): 35.262 ns … 28.399 μs ┊ GC (min … max): 0.00% … 99.64%
Time (median): 60.880 ns ┊ GC (median): 0.00%
Time (mean ± σ): 66.248 ns ± 501.417 ns ┊ GC (mean ± σ): 14.92% ± 1.99%
█▇ ▁ ▂
▆██▆▂▂▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▃▇▆█▇▄▂▂▁▁▁▁▁▂▂▃▆▆███▆▄▂▁▁ ▂
35.3 ns Histogram: frequency by time 76.3 ns <
Memory estimate: 96 bytes, allocs estimate: 2.
In short, DifferentiationInterface.jl allows for easy testing and comparison of AD backends. If you want to go further, check out the documentation of DifferentiationInterfaceTest.jl. This related package provides benchmarking utilities to compare backends and help you select the one that is best suited for your problem.