nbodykit.algorithms.threeptcf

Classes

Base3PCF(source, poles, edges, required_cols)

Base class for implementing common 3PCF calculations.

SimulationBox3PCF(source, poles, edges[, …])

Compute the multipoles of the isotropic, three-point correlation function in configuration space for data in a simulation box.

SurveyData3PCF(source, poles, edges, cosmo)

Compute the multipoles of the isotropic, three-point correlation function in configuration space for observational survey data.

YlmCache(ells, comm)

A class to compute spherical harmonics \(Y_{lm}\) up to a specified maximum \(\ell\).

class nbodykit.algorithms.threeptcf.Base3PCF(source, poles, edges, required_cols, BoxSize=None, periodic=None)[source]

Base class for implementing common 3PCF calculations.

Users should use SimulationBox3PCF or SurveyData3PCF.

Methods

load(filename[, comm])

Load a result from filename that has been saved to disk with save().

save(self, output)

Save the poles result to a JSON file with name output.

classmethod load(filename, comm=None)[source]

Load a result from filename that has been saved to disk with save().

save(self, output)[source]

Save the poles result to a JSON file with name output.

class nbodykit.algorithms.threeptcf.SimulationBox3PCF(source, poles, edges, BoxSize=None, periodic=True, weight='Weight', position='Position')[source]

Compute the multipoles of the isotropic, three-point correlation function in configuration space for data in a simulation box.

This uses the algorithm of Slepian and Eisenstein, 2015 which scales as \(\mathcal{O}(N^2)\), where \(N\) is the number of objects.

Results are computed when the object is inititalized. See the documenation of run() for the attributes storing the results.

Note

The algorithm expects the positions of objects in a simulation box to be the Cartesian x, y, and z vectors. For survey data, in the form of right ascension, declination, and redshift, see SurveyData3PCF.

Parameters
  • source (CatalogSource) – the input source of particles providing the ‘Position’ column

  • poles (list of int) – the list of multipole numbers to compute

  • edges (array_like) – the edges of the bins of separation to use; length of nbins+1

  • BoxSize (float, 3-vector, optional) – the size of the box; if periodic boundary conditions used, and ‘BoxSize’ not provided in the source attrs, it must be provided here

  • periodic (bool, optional) – whether to use periodic boundary conditions when computing separations between objects

  • weight (str, optional) – the name of the column in the source specifying the particle weights

References

Slepian and Eisenstein, MNRAS 454, 4142-4158 (2015)

Methods

load(filename[, comm])

Load a result from filename that has been saved to disk with save().

run(self[, pedantic])

Compute the three-point CF multipoles.

save(self, output)

Save the poles result to a JSON file with name output.

classmethod load(filename, comm=None)

Load a result from filename that has been saved to disk with save().

run(self, pedantic=False)[source]

Compute the three-point CF multipoles. This attaches the following the attributes to the class:

poles

a BinnedStatistic object to hold the multipole results; the binned statistics stores the multipoles as variables corr_0, corr_1, etc for \(\ell=0,1,\) etc. The coordinates of the binned statistic are r1 and r2, which give the separations between the three objects in CF.

Type

BinnedStatistic

save(self, output)

Save the poles result to a JSON file with name output.

class nbodykit.algorithms.threeptcf.SurveyData3PCF(source, poles, edges, cosmo, domain_factor=4, ra='RA', dec='DEC', redshift='Redshift', weight='Weight')[source]

Compute the multipoles of the isotropic, three-point correlation function in configuration space for observational survey data.

This uses the algorithm of Slepian and Eisenstein, 2015 which scales as \(\mathcal{O}(N^2)\), where \(N\) is the number of objects.

Results are computed when the object is inititalized. See the documenation of run() for the attributes storing the results.

Note

The algorithm expects the positions of objects from a survey catalog be the sky coordinates, right ascension and declination, and redshift. For simulation box data in Cartesian coordinates, see SimulationBox3PCF.

Warning

The right ascension and declination columns should be specified in degrees.

Parameters
  • source (CatalogSource) – the input source of particles providing the ‘Position’ column

  • poles (list of int) – the list of multipole numbers to compute

  • edges (array_like) – the edges of the bins of separation to use; length of nbins+1

  • cosmo (Cosmology) – the cosmology instance used to convert redshifts into comoving distances

  • ra (str, optional) – the name of the column in the source specifying the right ascension coordinates in units of degrees; default is ‘RA’

  • dec (str, optional) – the name of the column in the source specifying the declination coordinates; default is ‘DEC’

  • redshift (str, optional) – the name of the column in the source specifying the redshift coordinates; default is ‘Redshift’

  • weight (str, optional) – the name of the column in the source specifying the object weights

  • domain_factor (int, optional) – the integer value by which to oversubscribe the domain decomposition mesh before balancing loads; this number can affect the distribution of loads on the ranks – an optimal value will lead to balanced loads

References

Slepian and Eisenstein, MNRAS 454, 4142-4158 (2015)

Methods

load(filename[, comm])

Load a result from filename that has been saved to disk with save().

run(self)

Compute the three-point CF multipoles.

save(self, output)

Save the poles result to a JSON file with name output.

classmethod load(filename, comm=None)

Load a result from filename that has been saved to disk with save().

run(self)[source]

Compute the three-point CF multipoles. This attaches the following the attributes to the class:

poles

a BinnedStatistic object to hold the multipole results; the binned statistics stores the multipoles as variables corr_0, corr_1, etc for \(\ell=0,1,\) etc. The coordinates of the binned statistic are r1 and r2, which give the separations between the three objects in CF.

Type

BinnedStatistic

save(self, output)

Save the poles result to a JSON file with name output.

class nbodykit.algorithms.threeptcf.YlmCache(ells, comm)[source]

A class to compute spherical harmonics \(Y_{lm}\) up to a specified maximum \(\ell\).

During calculation, the necessary power of Cartesian unit vectors are cached in memory to avoid repeated calculations for separate harmonics.

Methods

__call__(self, xpyhat, zhat)

Return a dictionary holding Ylm for each (l,m) combination required

__call__(self, xpyhat, zhat)[source]

Return a dictionary holding Ylm for each (l,m) combination required

Parameters
  • xpyhat (array_like) – a complex array holding xhat + i * yhat, where xhat and yhat are the two cartesian unit vectors

  • zhat (array_like) – the third cartesian unit vector