Using automatic differentiation

This tutorial shows how to compute the power spectrum $P(k; \theta)$ and its (logarithmic) derivatives

\[\frac{\partial \lg P}{\partial \lg \theta_i}\]

using automatic differentiation with ForwardDiff.jl. This technique can also differentiate any other quantity.

1. Wrap the evaluation

First, we must decide which parameters $\theta$ the power spectrum $P(k; \theta)$ should be considered a function of. To do so, let us write a small wrapper function that calculates the power spectrum as a function of the parameters $(T_{\gamma 0}, \Omega_{c0}, \Omega_{b0}, N_\textrm{eff}, h, Y_p)$, following the Getting started tutorial:

using SymBoltz
M = ΛCDM(K = nothing)
pars = [M.γ.T₀, M.c.Ω₀, M.b.Ω₀, M.ν.Neff, M.g.h, M.b.YHe, M.h.m_eV, M.I.ln_As1e10, M.I.ns]
prob0 = CosmologyProblem(M, Dict(pars .=> NaN))

probgen = parameter_updater(prob0, pars)
P(k, θ) = spectrum_matter(probgen(θ), k; verbose = true, ptopts = (reltol = 1e-3,))
P (generic function with 1 method)

It is now easy to evaluate the power spectrum:

using Unitful, UnitfulAstro
θ = [2.7, 0.27, 0.05, 3.0, 0.7, 0.25, 0.02, 3.0, 0.95]
ks = 10 .^ range(-3, 0, length=100) / u"Mpc"
Ps = P(ks, θ)
100-element Vector{Unitful.Quantity{Float64, 𝐋^3, Unitful.FreeUnits{(Mpc^3,), 𝐋^3, nothing}}}:
 15091.417709566385 Mpc^3
 16046.459687239376 Mpc^3
 17053.459245980255 Mpc^3
 18113.474192403635 Mpc^3
 19228.329564846626 Mpc^3
 20398.7495561304 Mpc^3
 21625.64662825441 Mpc^3
 22909.626555792373 Mpc^3
 24250.62729350334 Mpc^3
 25648.74261785118 Mpc^3
     ⋮
   234.8630248093804 Mpc^3
   200.6799788750559 Mpc^3
   171.18217336749046 Mpc^3
   145.98788436649153 Mpc^3
   124.34218351583318 Mpc^3
   105.6714380017265 Mpc^3
    89.79962705574447 Mpc^3
    76.2949236929751 Mpc^3
    64.6818460292661 Mpc^3

This can be plotted with

using Plots
plot(log10.(ks/u"1/Mpc"), log10.(Ps/u"Mpc^3"); xlabel = "lg(k/Mpc⁻¹)", ylabel = "lg(P/Mpc³)", label = nothing)
Example block output

2. Calculate the derivatives

To get $\partial \lg P / \partial \lg \theta$, we can simply pass the wrapper function P(k, θ) through ForwardDiff.jacobian:

using ForwardDiff
lgP(lgθ) = log10.(P(ks, 10 .^ lgθ) / u"Mpc^3") # in log-space
dlgP_dlgθs = ForwardDiff.jacobian(lgP, log10.(θ))
100×9 Matrix{Float64}:
 -0.148115  -1.23988  -0.229631  -0.00984893  …  -0.00456095  3.0  -3.71642
 -0.162171  -1.23403  -0.228547  -0.0108787      -0.00452996  3.0  -3.65014
 -0.177466  -1.22757  -0.227348  -0.0119977      -0.00452645  3.0  -3.58385
 -0.194103  -1.22058  -0.226052  -0.0132177      -0.00455518  3.0  -3.51756
 -0.212168  -1.21304  -0.224656  -0.0145427      -0.00461839  3.0  -3.45128
 -0.231778  -1.20501  -0.223167  -0.0159899   …  -0.00472015  3.0  -3.38499
 -0.252868  -1.19617  -0.22152   -0.0175526      -0.00486003  3.0  -3.3187
 -0.275677  -1.18662  -0.219765  -0.0192331      -0.0050385   3.0  -3.25242
 -0.300245  -1.17644  -0.21788   -0.0210689      -0.00525454  3.0  -3.18613
 -0.326617  -1.16552  -0.215866  -0.0230367      -0.00550654  3.0  -3.11984
  ⋮                                           ⋱                    
 -5.97639    1.88998  -0.368653  -0.382781       -0.0253154   3.0   2.31565
 -6.02217    1.90262  -0.366879  -0.385282       -0.0253177   3.0   2.38194
 -6.05807    1.91809  -0.364576  -0.386499       -0.0253177   3.0   2.44823
 -6.08032    1.93603  -0.367462  -0.386885       -0.0253198   3.0   2.51451
 -6.11367    1.94875  -0.365832  -0.388269    …  -0.0253207   3.0   2.5808
 -6.13775    1.96381  -0.367085  -0.389054       -0.0253215   3.0   2.64709
 -6.16419    1.9773   -0.368036  -0.39011        -0.0253213   3.0   2.71337
 -6.19072    1.99178  -0.367861  -0.39074        -0.0253216   3.0   2.77966
 -6.21756    2.0046   -0.367693  -0.391858       -0.0253223   3.0   2.84595

The matrix element dlgP_dlgθs[i, j] now contains $\partial \lg P(k_i) / \partial \lg \theta_j$. We can plot them all at once:

plot(
    log10.(ks/u"1/Mpc"), dlgP_dlgθs;
    xlabel = "lg(k/Mpc⁻¹)", ylabel = "∂ lg(P) / ∂ lg(θᵢ)",
    labels = "θᵢ=" .* ["Tγ0" "Ωc0" "Ωb0" "Neff" "h" "YHe" "mh" "ln(10¹⁰As)" "ns"]
)
Example block output

Get values and derivatives together

The above example showed how to calculate the power spectrum values and their derivatives through two separate calls. If you need both, it is faster to calculate them simultaneously with the package DiffResults.jl:

using DiffResults

# Following DiffResults documentation:
Pres = DiffResults.JacobianResult(ks/u"1/Mpc", θ) # allocate buffer for values+derivatives for a function with θ-sized input and ks-sized output
Pres = ForwardDiff.jacobian!(Pres, lgP, log10.(θ)) # evaluate values+derivatives of lgP(log10.(θ)) and store the results in Pres
lgPs = DiffResults.value(Pres) # extract values
dlgP_dlgθs = DiffResults.jacobian(Pres) # extract derivatives

p1 = plot(
   log10.(ks/u"1/Mpc"), lgPs;
   ylabel = "lg(P/Mpc³)", label = nothing
)
p2 = plot(
   log10.(ks/u"1/Mpc"), dlgP_dlgθs;
   xlabel = "lg(k/Mpc⁻¹)", ylabel = "∂ lg(P) / ∂ lg(θᵢ)",
   labels = "θᵢ=" .* ["Tγ0" "Ωc0" "Ωb0" "Neff" "h" "YHe" "mh" "ln(10¹⁰As)" "ns"]
)
plot(p1, p2, layout=(2, 1), size = (600, 600))
Example block output

General approach

The technique shown here can be used to calculate the derivative of any SymBoltz.jl output quantity:

  1. Write a wrapper function output(input) that calculates the desired output quantities from the desired input quantities.
  2. Use ForwardDiff.derivative(output, input) (scalar-to-scalar), ForwardDiff.gradient(output, input) (vector-to-scalar) or ForwardDiff.jacobian(output, input) (vector-to-vector) to evaluate the derivative of output at the values input. Or use the similar functions in DiffResults to calculate the value and derivatives simultaneously.