nbodykit.algorithms.convpower

Functions

copy_meta(attrs, meta[, prefix])
get_compensation(mesh)
get_real_Ylm(l, m) Return a function that computes the real spherical
is_valid_crosscorr(first, second)

Classes

ConvolvedFFTPower(first, poles[, second, …]) Algorithm to compute power spectrum multipoles using FFTs for a data survey with non-trivial geometry.
class nbodykit.algorithms.convpower.ConvolvedFFTPower(first, poles, second=None, Nmesh=None, kmin=0.0, dk=None, use_fkp_weights=False, P0_FKP=None)[source]

Algorithm to compute power spectrum multipoles using FFTs for a data survey with non-trivial geometry.

Due to the geometry, the estimator computes the true power spectrum convolved with the window function (FFT of the geometry).

This estimator implemented in this class is described in detail in Hand et al. 2017 (arxiv:1704.02357). It uses the spherical harmonic addition theorem such that only \(2\ell+1\) FFTs are required to compute each multipole. This differs from the implementation in Bianchi et al. and Scoccimarro et al., which requires \((\ell+1)(\ell+2)/2\) FFTs.

Results are computed when the object is inititalized, and the result is stored in the poles attribute. Important meta-data computed during algorithm execution is stored in the attrs dict. See the documenation of run().

Note

A full tutorial on the class is available in the documentation here.

Note

Cross correlations are only supported when the FKP weight column differs between the two mesh objects, i.e., the underlying data and randoms must be the same. This allows users to compute the cross power spectrum of the same density field, weighted differently.

Parameters:
  • first (FKPCatalog, FKPCatalogMesh) – the first source to paint the data/randoms; FKPCatalog is automatically converted to a FKPCatalogMesh, using default painting parameters
  • poles (list of int) – a list of integer multipole numbers ell to compute
  • second (FKPCatalog, FKPCatalogMesh, optional) – the second source to paint the data/randoms; cross correlations are only supported when the weight column differs between the two mesh objects, i.e., the underlying data and randoms must be the same!
  • kmin (float, optional) – the edge of the first wavenumber bin; default is 0
  • dk (float, optional) – the spacing in wavenumber to use; if not provided; the fundamental mode of the box is used
  • use_fkp_weights (bool, optional) – if True, FKP weights will be added using P0_FKP such that the fkp weight is given by 1 / (1 + P0*NZ) where NZ is the number density as a function of redshift column
  • P0_FKP (float, optional) – the value of P0 to use when computing FKP weights; must not be None if use_fkp_weights=True

References

  • Hand, Nick et al. An optimal FFT-based anisotropic power spectrum estimator, 2017
  • Bianchi, Davide et al., Measuring line-of-sight-dependent Fourier-space clustering using FFTs, MNRAS, 2015
  • Scoccimarro, Roman, Fast estimators for redshift-space clustering, Phys. Review D, 2015

Methods

load(output[, comm]) Load a saved ConvolvedFFTPower result, which has been saved to disk with ConvolvedFFTPower.save().
normalization(name, alpha) Compute the power spectrum normalization, using either the data or randoms source.
run() Compute the power spectrum multipoles.
save(output) Save the ConvolvedFFTPower result to disk.
shotnoise(alpha) Compute the power spectrum shot noise, using either the data or randoms source.
to_pkmu(mu_edges, max_ell) Invert the measured multipoles \(P_\ell(k)\) into power spectrum wedges, \(P(k,\mu)\).
classmethod load(output, comm=None)[source]

Load a saved ConvolvedFFTPower result, which has been saved to disk with ConvolvedFFTPower.save().

The current MPI communicator is automatically used if the comm keyword is None

normalization(name, alpha)[source]

Compute the power spectrum normalization, using either the data or randoms source.

The normalization is given by:

\[A = \int d^3x \bar{n}'_1(x) \bar{n}'_2(x) w_{\mathrm{fkp},1} w_{\mathrm{fkp},2}.\]

The mean densities are assumed to be the same, so this can be converted to a summation over objects in the source, as

\[A = \sum w_{\mathrm{comp},1} \bar{n}_2 w_{\mathrm{fkp},1} w_{\mathrm{fkp},2}.\]

References

see Eqs. 13,14 of Beutler et al. 2014, “The clustering of galaxies in the SDSS-III Baryon Oscillation Spectroscopic Survey: testing gravity with redshift space distortions using the power spectrum multipoles”

run()[source]

Compute the power spectrum multipoles. This function does not return anything, but adds several attributes (see below).

edges

array_like – the edges of the wavenumber bins

poles

BinnedStatistic – a BinnedStatistic object that behaves similar to a structured array, with fancy slicing and re-indexing; it holds the measured multipole results, as well as the number of modes (modes) and average wavenumbers values in each bin (k)

attrs

dict – dictionary holding input parameters and several important quantites computed during execution:

  1. data.N, randoms.N :
    the unweighted number of data and randoms objects
  2. data.W, randoms.W :
    the weighted number of data and randoms objects, using the column specified as the completeness weights
  3. alpha :
    the ratio of data.W to randoms.W
  4. data.norm, randoms.norm :
    the normalization of the power spectrum, computed from either the “data” or “randoms” catalog (they should be similar). See equations 13 and 14 of arxiv:1312.4611.
  5. data.shotnoise, randoms.shotnoise :
    the shot noise values for the “data” and “random” catalogs; See equation 15 of arxiv:1312.4611.
  6. shotnoise :
    the total shot noise for the power spectrum, equal to data.shotnoise + randoms.shotnoise; this should be subtracted from the monopole.
  7. BoxSize :
    the size of the Cartesian box used to grid the data and randoms objects on a Cartesian mesh.

For further details on the meta-data, see the documentation.

save(output)[source]

Save the ConvolvedFFTPower result to disk.

The format is currently json.

Parameters:output (str) – the name of the file to dump the JSON results to
shotnoise(alpha)[source]

Compute the power spectrum shot noise, using either the data or randoms source.

This computes:

\[S = \sum (w_\mathrm{comp} w_\mathrm{fkp})^2\]

References

see Eq. 15 of Beutler et al. 2014, “The clustering of galaxies in the SDSS-III Baryon Oscillation Spectroscopic Survey: testing gravity with redshift space distortions using the power spectrum multipoles”

to_pkmu(mu_edges, max_ell)[source]

Invert the measured multipoles \(P_\ell(k)\) into power spectrum wedges, \(P(k,\mu)\).

Parameters:
  • mu_edges (array_like) – the edges of the \(\mu\) bins
  • max_ell (int) – the maximum multipole to use when computing the wedges; all even multipoles with \(ell\) less than or equal to this number are included
Returns:

pkmu – a data set holding the \(P(k,\mu)\) wedges

Return type:

BinnedStatistic

nbodykit.algorithms.convpower.get_real_Ylm(l, m)[source]

Return a function that computes the real spherical harmonic of order (l,m)

Parameters:
  • l (int) – the degree of the harmonic
  • m (int) – the order of the harmonic; abs(m) < l
Returns:

Ylm – a function that takes 4 arguments: (xhat, yhat, zhat) unit-normalized Cartesian coordinates and returns the specified Ylm

Return type:

callable

References

https://en.wikipedia.org/wiki/Spherical_harmonics#Real_form