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.06, 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}}}:
15200.828396694069 Mpc^3
16177.008593901686 Mpc^3
17208.627121615908 Mpc^3
18297.747676868334 Mpc^3
19445.685748071417 Mpc^3
20654.502223346502 Mpc^3
21925.18150897419 Mpc^3
23259.192083997845 Mpc^3
24657.268350056653 Mpc^3
26119.776392090313 Mpc^3
⋮
309.19751897582387 Mpc^3
264.36246009780007 Mpc^3
225.73890080477267 Mpc^3
192.54779287678053 Mpc^3
164.10453119805985 Mpc^3
139.6952572549726 Mpc^3
118.79163994609006 Mpc^3
100.92059285128487 Mpc^3
85.65865652784208 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)
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.126982 -1.24915 -0.231344 -0.010627 … -0.00509322 3.0 -3.71642
-0.138726 -1.24432 -0.230454 -0.0117581 -0.00499269 3.0 -3.65014
-0.151523 -1.23895 -0.229456 -0.0129933 -0.00488801 3.0 -3.58385
-0.165485 -1.233 -0.228348 -0.0143401 -0.00478057 3.0 -3.51756
-0.180772 -1.22664 -0.227171 -0.0158121 -0.00467357 3.0 -3.45128
-0.197395 -1.21963 -0.225872 -0.0174124 … -0.00456918 3.0 -3.38499
-0.215481 -1.21208 -0.224474 -0.0191534 -0.00447158 3.0 -3.3187
-0.235093 -1.20388 -0.222956 -0.0210431 -0.00438473 3.0 -3.25242
-0.256312 -1.19493 -0.221302 -0.0230897 -0.00431303 3.0 -3.18613
-0.279312 -1.18542 -0.219553 -0.0253069 -0.00426237 3.0 -3.11984
⋮ ⋱
-5.96369 1.906 -0.396841 -0.46792 -0.0298507 3.0 2.31565
-5.99511 1.92323 -0.397565 -0.468939 -0.0298623 3.0 2.38194
-6.02688 1.93958 -0.397973 -0.470059 -0.0298728 3.0 2.44823
-6.05978 1.95478 -0.397947 -0.471491 -0.029882 3.0 2.51451
-6.08764 1.97081 -0.398519 -0.472136 … -0.0298897 3.0 2.5808
-6.11799 1.98525 -0.398631 -0.473269 -0.0298969 3.0 2.64709
-6.14628 1.99947 -0.398858 -0.47426 -0.029903 3.0 2.71337
-6.17367 2.01352 -0.39901 -0.475044 -0.0299085 3.0 2.77966
-6.20088 2.027 -0.398889 -0.475909 -0.0299128 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"]
)
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))
General approach
The technique shown here can be used to calculate the derivative of any SymBoltz.jl output quantity:
- Write a wrapper function
output(input)
that calculates the desired output quantities from the desired input quantities. - Use
ForwardDiff.derivative(output, input)
(scalar-to-scalar),ForwardDiff.gradient(output, input)
(vector-to-scalar) orForwardDiff.jacobian(output, input)
(vector-to-vector) to evaluate the derivative ofoutput
at the valuesinput
. Or use the similar functions inDiffResults
to calculate the value and derivatives simultaneously.