In [1]:
import pymaster as nmt
import numpy as np
import healpy as hp
import matplotlib.pyplot as plt

These notebooks aim to provide a comprehensive review of the functionality implemented in NaMaster. No rigorous mathematical derivations of any statements made will be provided, and users are referred to the references below for further details:

  • A19. The original NaMaster paper (Alonso et al. 2019 https://arxiv.org/abs/1809.09603).
  • G19. Garcia-Garcia et al. 2019 (https://arxiv.org/abs/1906.11765), which introduces the basic approximation used by NaMaster to estimate covariance matrix.
  • N20. Nicola et al. 2020 (https://arxiv.org/abs/2010.09717), which refined these approximations and described in detail the procedure to estimate cosmic shear power spectra.
  • W24. Wolz et al. 2024 (TBD), which introduced the formalism for catalog-based $C_\ell$s.

1 Maps, masks, pseudo-$C_\ell$s, bandpowers ¶

This first tutorial presents the main basic concepts behind power spectrum estimation: maps, masks, pseudo-$C_\ell$s, mode-coupling matrices, bandpowers, and bandpower window functions.

In most common situations, we aim to calculate the angular power spectrum of a map $a(\hat{\bf n})$, of which we can only access a masked version $\tilde{a}(\hat{\bf n})\equiv w(\hat{\bf n})\,a(\hat{\bf n})$, where $w(\hat{\bf n})$ is the mask.

Naively, we could ignore the masked nature of our map and use the following simple estimator to calculate its power spectrum: $$ \tilde{C}_\ell=\frac{1}{2\ell+1}\sum_{\ell=-m}^m\,|\tilde{a}_{\ell m}|^2, $$ where $a_{\ell m}$ are the spherical harmonic coefficients of the field. This is commonly called the "pseudo-$C_\ell$" of the map (or in general, of any pair of maps), and gives its name to the estimator used by NaMaster. Below, we may abbreviate this as PCL.

As per the convolution theorem, the spherical harmonic transform of the masked map is a convolution of the harmonic coefficients of the mask and the true map. By the same token, the power spectrum of the masked map (or, more precisely, the ensemble average of the pseudo-$C_\ell$ as defined above), is a convolution of the power spectrum of the true map, $C_\ell$, with a mode-coupling matrix (MCM) $M_{\ell \ell'}$ that depends exclusively on the power spectrum of the mask: $$ \langle \tilde{C}_\ell\rangle = \sum_{\ell'}\,M_{\ell\ell'}\,C_{\ell'}. $$ Analytical expressions for $M_{\ell\ell'}$ can be found in Section 2.1.2 of A19, and constitute a key component of NaMaster.

Since the mode-coupling matrix is not necessarily invertible, a common procedure to obtain a close-to-unbiased estimator for $C_\ell$ is to consider bins of $\ell$ (as an effective way to regularise $M_{\ell \ell'}$). These bins of $\ell$ are often called bandpowers. The resulting binned mode-coupling matrix is thus invertible and can be used to estimate mode-decoupled bandpowers (the NaMaster documentation often refers to these as the "decoupled $C_\ell$").

This decoupled bandpowers $C_b$ (where $b$ labels the bandpower), are related to the true $C_\ell$ through the bandpower window functions $W_{b\ell}$ (labelled ${\cal F}_{b\ell}$ in Eq. 19 of A19). Under ideal circumstances (no mode coupling), $W_{b\ell}$ would be a set of top-hat functions representing the averaging of the original $C_\ell$ into bandpowers. In most practical situations, however, the bandpower window functions encode any residual mode-coupling effects, which may be significant and should therefore be taken into account.

The code below is a quick demonstration of the PCL estimator as implemented in NaMaster. It generates 100 simulations with a known true power spectrum and, for a non-trivial mask, shows that the resulting estimated power spectra are unbiased once the bandpower window functions are taken into account. In order to get there, we will introduce the different elements of the calculation one by one.

Table of contents¶

  • 1.1 Maps and masks
  • 1.2 Fields
  • 1.3 Bins
  • 1.4 Pseudo-$C_\ell$s
  • 1.5 Mode-coupling matrices
  • 1.6 Mode-decoupled power spectra
  • 1.7 Unbiased estimator

1.1 Maps and masks ¶

The code below generates Gaussian realisations of a power spectrum of the form $C_\ell = 1/(\ell+10)$. We will assume a relatively simple mask: a circular footprint with 20-degree radius.

It is often useful (although not strictly necessary) to apodise your mask. Although this can lead to a reduction in the uncertainties of the resulting estimator, if the true underlying power spectrum is steep, or if the mask is very patchy, it can also do the opposite: apodising means downweighting parts of the sky that are close to the edges of your mask, which effectively translates into a loss of area or a suppression of certain modes. Depending on the case, this can therefore lead to an increase in the final uncertainties.

In [2]:
# HEALPix resolution parameter
nside = 256
npix = hp.nside2npix(nside)

# Range of ells to explore and input power spectrum
ls = np.arange(3*nside)
cl_true = 1/(ls+10)

# Simulator
def gen_sim():
    return hp.synfast(cl_true, nside)

# Binary mask
mask_binary = np.zeros(npix)
mask_binary[hp.query_disc(nside, [1, 0, 0], np.radians(20))] = 1

# Apodised mask
# We use a 5-degree "C1" apodisation (see Grain et al. 2009).
mask = nmt.mask_apodization(mask_binary, 5.0, "C1")

fig, axes = plt.subplots(2, 2, figsize=(14, 10))
plt.axes(axes[0][0])
hp.mollview(mask_binary, title='Binary mask', hold=True)
plt.axes(axes[0][1])
hp.mollview(mask, title='Apodized mask', hold=True)
mp = gen_sim()
plt.axes(axes[1][0])
hp.mollview(mp, title='Full-sky map', hold=True)
plt.axes(axes[1][1])
hp.mollview(mp*mask, title='Masked map', hold=True)
No description has been provided for this image

1.2 Fields ¶

Once you have a map and its mask, the first step to interact with NaMaster is to create a NmtField object from them. NmtField objects take care of storing all relevant information about any given map we want to correlate, and of, internally, performing operations on them, including spherical harmonic transforms, deprojection, purification, etc. (which we will cover below).

Note that, when constructing an NmtField object, the spherical harmonic transforms of the map and the associated mask will be computed automatically. Depending on the resolution of your map, this can be an expensive operation! You can limit the time taken to do this by using the n_iter argument, which limits the number of Jacobi iterations used to increase the accuracy of the recovered harmonic coefficients. The code below shows that, once initialised, you can use a field to access the mask and masked maps that were passed to create it, as well as their spherical harmonic coefficients.

In [3]:
f = nmt.NmtField(mask, [mp], n_iter=0)
fig, axes = plt.subplots(2, 2, figsize=(14, 10))
plt.axes(axes[0][0])
hp.mollview(mask, title='Mask', hold=True)
plt.axes(axes[0][1])
hp.mollview(f.get_mask(), title='Mask from NmtField', hold=True)
plt.axes(axes[1][0])
hp.mollview(mp*mask, title='Masked map', hold=True)
plt.axes(axes[1][1])
hp.mollview(f.get_maps()[0], title='Masked map from NmtField', hold=True)

print("The first 5 harmonic coefficients of the masked map are:")
print(f.get_alms()[0][:5])
print("The first 5 harmonic coefficients of the mask are:")
print(f.get_mask_alms()[:5])
The first 5 harmonic coefficients of the masked map are:
[-0.01108506+0.j  0.00632131+0.j  0.0126594 +0.j -0.01295433+0.j
 -0.01352294+0.j]
The first 5 harmonic coefficients of the mask are:
[ 8.27262120e-02+0.j -6.89504011e-20+0.j -8.60449099e-02+0.j
  1.49830650e-19+0.j  7.26256197e-02+0.j]
No description has been provided for this image

1.3 Bins ¶

Integral to NaMaster's estimator is the binning scheme you want to use for your bandpowers. This is controlled by NmtBin objects. The code below showcases their basic functionality.

In [4]:
# The default constructor for NmtBins allows you to fine tune
# each of your bandpowers. Here we show only some of the more
# basic constructors (provided as class methods).

# This constructor builds equi-space bandpowers from ell=2
# up to the maximum ell that maps with resolution nside should
# be able to reproduce (3*nside-1)
delta_ell = 30
b1 = nmt.NmtBin.from_nside_linear(nside, nlb=delta_ell)

# This one is very similar. Instead of passing nside, you
# pass the maximum multipole.
b2 = nmt.NmtBin.from_lmax_linear(3*nside-1, nlb=delta_ell)

# Finally, this one constructs bandpowers from their
# corresponding edges.
edges = np.arange(2, 3+delta_ell*((3*nside-1)//delta_ell), delta_ell)
b3 = nmt.NmtBin.from_edges(edges[:-1], edges[1:])

# All the bins above should be equivalent, so let's just
# pick one.
b = b1

# NmtBin contain all information you need about the
# binning scheme you chose. Some basic functionality:
# - effective ells for each bandpower
leff = b1.get_effective_ells()
# - simple binning function
cl_binned = b.bin_cell(cl_true)
# - unbinning (assuming stepwise)
cl_unbinned = b.unbin_cell(cl_binned)
# - binning weights for each bandpower
n_bins = b.get_n_bands()
Fls = np.zeros([n_bins, 3*nside])
for i in range(n_bins):
    Fls[i, b.get_ell_list(i)] = b.get_weight_list(i)

plt.figure(figsize=(7, 5))
plt.plot(ls, cl_true, label=r'Input $C_\ell$')
plt.plot(leff, cl_binned, label=r'Binned $C_\ell$')
plt.plot(ls, cl_unbinned, label=r'Binned-unbinned $C_\ell$')
plt.yscale('log')
plt.xlabel(r'$\ell$', fontsize=15)
plt.ylabel(r'$C_\ell$', fontsize=15)
plt.legend(fontsize=12)

plt.figure(figsize=(7, 5))
plt.title('Bin weights')
for i, fl in enumerate(Fls):
    plt.plot(ls, fl, 'k-', alpha=(i+1)/n_bins)
plt.xlabel(r'$\ell$', fontsize=15)
Out[4]:
Text(0.5, 0, '$\\ell$')
No description has been provided for this image
No description has been provided for this image

1.4 Pseudo-$C_\ell$s ¶

The code below shows how to calculate the PCL of a given map, and how the result compares with the input power spectrum.

It is interesting that, in the case of perfecly white noise maps, in which the power spectrum is constant, the PCL is exactly related to the true power spectrum via a rescaling by the sky fraction, given by the sky average of the mask product $f_{\rm sky}^{ab}\equiv\langle w_a w_b\rangle$, where $w_a/w_b$ are the masks of the two fields under study (see N20): $$ \left\langle \frac{1}{2\ell+1}\sum_{m=-\ell}^\ell \tilde{a}_{\ell m}\tilde{b}^*_{\ell m}\right\rangle=\langle w_a w_b\rangle C_\ell^{a,b}\hspace{12pt}{\rm for\,\,white\,\,noise}. $$ The code below shows that this approximation can be used as a first-order correction to the PCL, although in general it is not enough to fully account for the effects of mode-coupling (you can in fact see how the different $\ell$s display correlated oscillations by eye). NaMaster thus does not rely on this approximation to compute power spectra, but we will be able to take advantage of it to estimate accurate but cheap covariance matrices.

In [5]:
pcl = nmt.compute_coupled_cell(f, f)
fsky = np.mean(f.get_mask()**2)
pcl_corrected = pcl / fsky

plt.figure(figsize=(7, 5))
plt.plot(ls, cl_true, label=r'Input $C_\ell$')
plt.plot(ls, pcl[0], label=r'PCL')
plt.plot(ls, pcl_corrected[0], label=r'PCL/$f_{\rm sky}$')
plt.yscale('log')
plt.xlabel(r'$\ell$', fontsize=15)
plt.ylabel(r'$C_\ell$', fontsize=15)
plt.legend(fontsize=12)
Out[5]:
<matplotlib.legend.Legend at 0x70a588478d40>
No description has been provided for this image

1.5 Mode-coupling matrices ¶

The code below shows how to compute the mode-coupling matrix of a given pair of maps. We will estimate power spectra in the same bandpowers displayed in Section 1.1.3, and the code shows how to calculate the associated bandpower window functions, and their differences with respect to the naive top-hat binning they represent (a result of the residual mode-coupling caused by inverting the mode-coupling matrix only after binning).

In [6]:
# Create a NaMaster workspace
w = nmt.NmtWorkspace.from_fields(f, f, b)

# Extract the mode-coupling matrix
mcm = w.get_coupling_matrix()

# Extract bandpower window functions
Bbl = w.get_bandpower_windows().squeeze()

# Plot MCM
plt.figure(figsize=(5, 5))
plt.title('Mode-coupling matrix')
plt.imshow(mcm, cmap='bone')
plt.ylabel("$\\ell$", fontsize=15)
plt.xlabel("$\\ell'$", fontsize=15)
plt.colorbar()

plt.figure(figsize=(7, 5))
for ll in [10, 50, 100, 150, 200]:
    plt.plot(ls, mcm[ll]/fsky, label=f'$\\ell={ll}$')
plt.xlabel("$\\ell'$", fontsize=15)
plt.ylabel("$M_{\\ell \\ell'}/f_{\\rm sky}$", fontsize=15)
plt.legend(fontsize=12)
plt.xlim([0, 300])

plt.figure(figsize=(7, 5))
plt.title('Bins vs. bandpower windows')
for i, fl in enumerate(Fls):
    plt.plot(ls, fl, 'k-', alpha=(i+1)/n_bins)
    plt.plot(ls, Bbl[i], 'r-', alpha=(i+1)/n_bins)
plt.xlabel(r'$\ell$', fontsize=15);
No description has been provided for this image
No description has been provided for this image
No description has been provided for this image

Note that you can also normalise your bandpower window functions according to the "FKP" convention, where instead of inverting the binned mode-coupling matrix, you simply divide it by a quantity (the effective sky fraction) to guarantee that the amplitude of any flat power spectrum is conserved. For more details see Section 2.1 of Baleato & White 2024. The code below just show the differences between both normalisations (FKP in blue, standard in red).

In [7]:
# Create a NaMaster workspace with FKP normalization
w = nmt.NmtWorkspace.from_fields(f, f, b, normalization='FKP')

# Extract bandpower window functions
Bbl_f = w.get_bandpower_windows().squeeze()

plt.figure(figsize=(7, 5))
plt.title('FKP vs. standard normalization')
for i, fl in enumerate(Fls):
    plt.plot(ls, Bbl[i], 'r-', alpha=(i+1)/n_bins)
    plt.plot(ls, Bbl_f[i], 'b-', alpha=(i+1)/n_bins)
plt.xlabel(r'$\ell$', fontsize=15);
No description has been provided for this image

Now that we have introduced mode-coupling matrices, we can use them to understand the $f_{\rm sky}$ approximation we saw in the last section. It turns out this approximation is valid due to a neat result regarding mode-coupling matrices: the sum over $\ell'$ of the MCM is exactly $f_{\rm sky}^{ab}$ as defined above (i.e. $\sum_{\ell'} M_{\ell\ell'}=\langle w_a w_b\rangle$). Let's show this!

In [7]:
print("Sum of MCM over ell': %.5lf" % (np.dot(mcm, np.ones(len(mcm)))[0]))
print("fsky: %.5lf" % fsky)
Sum of MCM over ell': 0.02196
fsky: 0.02196

1.6 Mode-decoupled power spectra ¶

Once the mode-coupling matrix has been estimated, we can then bin it, invert it, and use it to compute the mode-decoupled bandpowers. The result should then be compared with the theoretical prediction convolved with the bandpower window functions. The code below shows how to do this with inbuilt NaMaster methods.

In [8]:
# First, estimated decoupled bandpowers. The explicit way of doing this would be to compute
# the binned PCL and the binned MCM. Then, invert the binned MCM, and multiply it by the PCL.
# All this can be done with a simple NaMaster call:
cl_dec = w.decouple_cell(pcl)

# Explicitly, let's do this step by step (although there is no need to!)
mcm_binned_rows = np.array([b.bin_cell(row) for row in mcm])
mcm_binned = np.array([b.bin_cell(col) for col in mcm_binned_rows.T]).T*b.get_nell_list()
mcm_binned_inv = np.linalg.inv(mcm_binned)
pcl_binned = b.bin_cell(pcl)
cl_dec_b = np.einsum('jk,ik', mcm_binned_inv, pcl_binned)
assert np.allclose(cl_dec, cl_dec_b)

# Now let's compute the bandpower-window-convolved theoretical prediction. To do this, we
# should subject the input C_ell to the same steps that the real maps have undergone, namely
# mode-coupling, binning, and multiplication by the inverse binned MCM. Again, there are
# different ways to do this with NaMaster (some shorter than others). We show these below
# in order of convenience.

# 1. Just convolve theory with bandpower window functions:
cl_dec_th = np.dot(Bbl, cl_true)

# 2. Mode-couple theory and then decouple it
cl_dec_th_b = w.decouple_cell(w.couple_cell([cl_true]))[0]
assert np.allclose(cl_dec_th, cl_dec_th_b)

# 3. The artisan way
# Mode-couple and bin
pcl_th_binned = b.bin_cell(np.dot(mcm, cl_true))
# Multiply by inverse binned MCM
cl_dec_th_c = np.dot(mcm_binned_inv, pcl_th_binned)
assert np.allclose(cl_dec_th, cl_dec_th_c)


plt.figure(figsize=(7, 5))
plt.plot(ls, cl_true, 'r-', label='Input, unbinned')
plt.plot(leff, cl_dec_th, 'b-', label='Input, bpw-convolved')
plt.plot(leff, cl_dec[0], 'k.', label=r'Decoupled $C_\ell$')
plt.yscale('log')
plt.xlabel(r'$\ell$', fontsize=15)
plt.ylabel(r'$C_\ell$', fontsize=15)
plt.legend(fontsize=12)
Out[8]:
<matplotlib.legend.Legend at 0x71e72ed7c860>
No description has been provided for this image

1.7 Unbiased estimator ¶

We now show, by averaging over 100 simulations, that the power spectrum estimator implemented in NaMaster is unbiased once the residual impact of mode coupling has been taken into account. The gray band in the final plot shows the range of $\ell$s ($\ell>2 N_{\rm side}$) over which the result should not be trusted without testint it first (although in this case we see the estimator actually does a very good job! This is likely due to the simplicity of our mask.).

In [9]:
nsims = 100
cls_sims = []
for i in range(nsims):
    f = nmt.NmtField(mask, [gen_sim()], n_iter=0)
    cls_sims.append(w.decouple_cell(nmt.compute_coupled_cell(f, f)))
cls_sims = np.array(cls_sims).squeeze()

cls_mean = np.mean(cls_sims, axis=0)
cls_scatter = np.std(cls_sims, axis=0)
cls_error = cls_scatter/np.sqrt(nsims)

fig = plt.figure(figsize=(8, 5))
ax = fig.add_axes((.1,.3,.8,.6))
ax.plot(leff, cl_dec_th, 'b-', label='Input, bpw-convolved')
ax.plot(leff, b.bin_cell(cl_true), 'r--', label='Input, naive binning')
ax.errorbar(leff, cls_mean, yerr=cls_error, fmt='k.', label=r'$C_\ell$ from sims')
ax.set_yscale('log')
ax.set_ylabel(r'$C_\ell$', fontsize=15)
ax.set_xticklabels([])
ax.legend()
ax.set_xlim([-0.5, 3*nside])
ax.axvspan(2*nside, 3*nside, facecolor='#AAAAAA', alpha=0.5)

ax = fig.add_axes((.1,.1,.8,.2))
ax.errorbar(leff, (cls_mean-b.bin_cell(cl_true))/cls_error, yerr=np.ones_like(leff), fmt='r.--')
ax.errorbar(leff, (cls_mean-cl_dec_th)/cls_error, yerr=np.ones_like(leff), fmt='b.')
ax.axhline(0, c='k', ls='--')
ax.set_ylim([-2.9, 2.9])
ax.set_xlim([-0.5, 3*nside])
ax.axvspan(2*nside, 3*nside, facecolor='#AAAAAA', alpha=0.5)
ax.set_ylabel(r'$\Delta C_\ell/\sigma(C_\ell)$', fontsize=16)
ax.set_xlabel(r'$\ell$', fontsize=16);
No description has been provided for this image

Before we conclude, it is worth addressing a few common questions most users have when confronted with pseudo-$C_\ell$s for the first time.

  • Pixel window functions. It is worth noting that, above, we have not corrected for the impact of the pixel window function in any of our power spectrum calculations. This is because all our simulated realisations were generated simply using healpy's synfast method, which effectively evaluates the spherical harmonic functions at the pixel positions, rather than averaging them over the pixel area (which would induce a pixel window function). In fact, it is important not to correct for pixel window functions by default, since this correction is not always appropriate. I.e. before correcting for any smoothing effects, it is worth thinking about how the map you are analysing came to be, and how fluctuations below your pixel scale were treated. The usual pixel window function is appropriate if your map values at a given pixel is the average of an infinite-resolution map over the pixel area (which admittedly is a very common scenario). This notebook has not actually discussed the treatment of beam window functions in NaMaster's fields. Hopefully this will be remedied soon, but in the meantime do take a look at the documentation of the NmtField class.
  • Pseudo-$C_\ell$ bandpowers can be negative, even for auto-spectra. Indeed, even though the auto-spectrum of a given field is, by definition, positive, the pseudo-$C_\ell$ estimator of that power spectrum for a given realisation need not be so. The only guarantee is that the estimator will be unbiased in the mean over many realisations. Some negative bandpowers, particularly at low-$\ell$ where the variance is higher, are far from a rare occurrence, especially if your mask is a bit complex.