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.347613175154 Mpc^3
 16046.484622976406 Mpc^3
 17053.41010030566 Mpc^3
 18113.512044321735 Mpc^3
 19228.30177606352 Mpc^3
 20398.815914758343 Mpc^3
 21625.711290915657 Mpc^3
 22909.752769262435 Mpc^3
 24250.80018533999 Mpc^3
 25648.623134165347 Mpc^3
     ⋮
   234.89547125265833 Mpc^3
   200.59349180722376 Mpc^3
   171.25918850117833 Mpc^3
   146.0908699963315 Mpc^3
   124.32849526421477 Mpc^3
   105.71110166730908 Mpc^3
    89.78068654688566 Mpc^3
    76.25079842734091 Mpc^3
    64.68884756246591 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.148147  -1.23999  -0.22965   -0.0098548  …  -0.00456042  3.0  -3.71642
 -0.162166  -1.234    -0.228541  -0.010879      -0.00453003  3.0  -3.65014
 -0.177468  -1.22758  -0.227349  -0.0119988     -0.00452723  3.0  -3.58385
 -0.194107  -1.22059  -0.226055  -0.0132179     -0.00455534  3.0  -3.51756
 -0.212205  -1.21311  -0.22467   -0.0145479     -0.00461911  3.0  -3.45128
 -0.231761  -1.20492  -0.223151  -0.0159899  …  -0.0047199   3.0  -3.38499
 -0.252889  -1.19613  -0.22152   -0.0175504     -0.00485982  3.0  -3.3187
 -0.275644  -1.1866   -0.219749  -0.0192417     -0.00503861  3.0  -3.25242
 -0.300239  -1.17646  -0.217883  -0.0210705     -0.00525498  3.0  -3.18613
 -0.326597  -1.16545  -0.215849  -0.0230417     -0.005506    3.0  -3.11984
  ⋮                                          ⋱                    
 -5.97947    1.88708  -0.36891   -0.3837        -0.0253153   3.0   2.31565
 -6.01399    1.90433  -0.368799  -0.384556      -0.0253174   3.0   2.38194
 -6.0517     1.91753  -0.369855  -0.387364      -0.0253188   3.0   2.44823
 -6.07695    1.93555  -0.368038  -0.386918      -0.0253197   3.0   2.51451
 -6.11257    1.94947  -0.365679  -0.388009   …  -0.0253206   3.0   2.5808
 -6.13938    1.96304  -0.366437  -0.388973      -0.0253207   3.0   2.64709
 -6.16287    1.97828  -0.368643  -0.389958      -0.0253227   3.0   2.71337
 -6.19003    1.9912   -0.368325  -0.391005      -0.0253234   3.0   2.77966
 -6.21888    2.00428  -0.367995  -0.392301      -0.0253234   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.