nbodykit.algorithms

class nbodykit.algorithms.FFTPower(first, mode, Nmesh=None, BoxSize=None, second=None, los=[0, 0, 1], Nmu=5, dk=None, kmin=0.0, poles=[])[source]

Algorithm to compute the 1d or 2d power spectrum and/or multipoles in a periodic box, using a Fast Fourier Transform (FFT).

This computes the power spectrum as the square of the Fourier modes of the density field, which are computed via a FFT.

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

Note

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

Parameters:
  • first (CatalogSource, MeshSource) – the source for the first field; if a CatalogSource is provided, it is automatically converted to MeshSource using the default painting parameters (via to_mesh())
  • mode ({'1d', '2d'}) – compute either 1d or 2d power spectra
  • Nmesh (int, optional) – the number of cells per side in the particle mesh used to paint the source
  • BoxSize (int, 3-vector, optional) – the size of the box
  • second (CatalogSource, MeshSource, optional) – the second source for cross-correlations
  • los (array_like , optional) – the direction to use as the line-of-sight; must be a unit vector
  • Nmu (int, optional) – the number of mu bins to use from \(\mu=[0,1]\); if mode = 1d, then Nmu is set to 1
  • dk (float, optional) – the linear spacing of k bins to use; if not provided, the fundamental mode of the box is used
  • kmin (float, optional) – the lower edge of the first k bin to use
  • poles (list of int, optional) – a list of multipole numbers ell to compute \(P_\ell(k)\) from \(P(k,\mu)\)

Methods

load(output[, comm]) Load a saved FFTPower result.
run() Compute the power spectrum in a periodic box, using FFTs.
save(output) Save the FFTPower result to disk.
run()[source]

Compute the power spectrum in a periodic box, using FFTs. This function returns nothing, but attaches several attributes to the class:

edges

array_like – the edges of the wavenumber bins

power

BinnedStatistic – a BinnedStatistic object that holds the measured \(P(k)\) or \(P(k,\mu)\). It stores the following variables:

  • k :

    the mean value for each k bin

  • mu : mode=2d only

    the mean value for each mu bin

  • power :

    complex array storing the real and imaginary components of the power

  • modes :

    the number of Fourier modes averaged together in each bin

poles

BinnedStatistic or None – a BinnedStatistic object to hold the multipole results \(P_\ell(k)\); if no multipoles were requested by the user, this is None. It stores the following variables:

  • k :
    the mean value for each k bin
  • power_L :
    complex array storing the real and imaginary components for the \(\ell=L\) multipole
  • modes :
    the number of Fourier modes averaged together in each bin
attrs

dict – dictionary of meta-data; in addition to storing the input parameters, it includes the following fields computed during the algorithm execution:

  • shotnoise : float
    the power Poisson shot noise, equal to \(V/N\), where \(V\) is the volume of the box and N is the total number of objects; if a cross-correlation is computed, this will be equal to zero
  • N1 : int
    the total number of objects in the first source
  • N2 : int
    the total number of objects in the second source
class nbodykit.algorithms.ProjectedFFTPower(first, Nmesh=None, BoxSize=None, second=None, axes=(0, 1), dk=None, kmin=0.0)[source]

The power spectrum of a field in a periodic box, projected over certain axes.

This is not really always physically meaningful, but convenient for making sense of Lyman-Alpha forest or lensing maps.

This is usually called the 1d power spectrum or 2d power spectrum.

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

Parameters:
  • first (CatalogSource, MeshSource) – the source for the first field; if a CatalogSource is provided, it is automatically converted to MeshSource using the default painting parameters (via to_mesh())
  • Nmesh (int, optional) – the number of cells per side in the particle mesh used to paint the source
  • BoxSize (int, 3-vector, optional) – the size of the box
  • second (CatalogSource, MeshSource, optional) – the second source for cross-correlations
  • axes (tuple) – axes to measure the power on. The axes not in the list will be averaged out. For example: - (0, 1) : project to x,y and measure power - (0) : project to x and measure power.
  • dk (float, optional) – the linear spacing of k bins to use; if not provided, the fundamental mode of the box is used
  • kmin (float, optional) – the lower edge of the first k bin to use

Methods

load(output[, comm]) Load a saved FFTPower result.
run() Run the algorithm.
save(output) Save the FFTPower result to disk.
run()[source]

Run the algorithm. This attaches the following attributes to the class:

edges

array_like – the edges of the wavenumber bins

power

BinnedStatistic – a BinnedStatistic object that holds the projected power. It stores the following variables:

  • k :
    the mean value for each k bin
  • power :
    complex array holding the real and imaginary components of the projected power
  • modes :
    the number of Fourier modes averaged together in each bin
class nbodykit.algorithms.FFTCorr(first, mode, Nmesh=None, BoxSize=None, second=None, los=[0, 0, 1], Nmu=5, dr=None, rmin=0.0, poles=[])[source]

Algorithm to compute the 1d or 2d correlation and/or multipoles in a periodic box, using a Fast Fourier Transform (FFT).

This computes the power spectrum as the square of the Fourier modes of the density field, which are computed via a FFT. Then it is transformed back to obtain the correlation function.

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

Note

This is very similar to FFTPower.

Parameters:
  • first (CatalogSource, MeshSource) – the source for the first field; if a CatalogSource is provided, it is automatically converted to MeshSource using the default painting parameters (via to_mesh())
  • mode ({'1d', '2d'}) – compute either 1d or 2d power spectra
  • Nmesh (int, optional) – the number of cells per side in the particle mesh used to paint the source
  • BoxSize (int, 3-vector, optional) – the size of the box
  • second (CatalogSource, MeshSource, optional) – the second source for cross-correlations
  • los (array_like , optional) – the direction to use as the line-of-sight; must be a unit vector
  • Nmu (int, optional) – the number of mu bins to use from \(\mu=[0,1]\); if mode = 1d, then Nmu is set to 1
  • dr (float, optional) – the linear spacing of r bins to use; if not provided, the fundamental mode of the box is used
  • rmin (float, optional) – the lower edge of the first r bin to use
  • poles (list of int, optional) – a list of multipole numbers ell to compute \(\xi_\ell(r)\) from \(\xi(r,\mu)\)

Methods

load(output[, comm]) Load a saved FFTPower result.
run() Compute the correlation function in a periodic box, using FFTs.
save(output) Save the FFTPower result to disk.
run()[source]

Compute the correlation function in a periodic box, using FFTs. This function returns nothing, but attaches several attributes to the class:

edges

array_like – the edges of the wavenumber bins

corr

BinnedStatistic – a BinnedStatistic object that holds the measured \(\xi(r)\) or \(\xi(r,\mu)\). It stores the following variables:

  • r :

    the mean value for each r bin

  • mu : mode=2d only

    the mean value for each mu bin

  • corr :

    real array storing the correlation function

  • modes :

    the number of modes averaged together in each bin

poles

BinnedStatistic or None – a BinnedStatistic object to hold the multipole results \(\xi_\ell(r)\); if no multipoles were requested by the user, this is None. It stores the following variables:

  • r :
    the mean value for each r bin
  • power_L :
    complex array storing the real and imaginary components for the \(\ell=L\) multipole
  • modes :
    the number of modes averaged together in each bin
attrs

dict – dictionary of meta-data; in addition to storing the input parameters, it includes the following fields computed during the algorithm execution:

  • shotnoise : float
    the power Poisson shot noise, equal to \(V/N\), where \(V\) is the volume of the box and N is the total number of objects; if a cross-correlation is computed, this will be equal to zero
  • N1 : int
    the total number of objects in the first source
  • N2 : int
    the total number of objects in the second source
class nbodykit.algorithms.KDDensity(source, margin=1.0)[source]

Estimate a proxy density based on the distance to the nearest neighbor. The result is proportional to the density but the scale is unspecified.

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

Parameters:
  • source (CatalogSource) – the input source of particles to compute the proxy density on; must specify the ‘Position’ column
  • margin (float, optional) – Padding region per parallel domain; relative to the mean seperation

Methods

run() Compute the density proxy.
run()[source]

Compute the density proxy. This attaches the following attribute:

density

array_like, length: size – a unit-less, proxy density value for each object on the local rank. This is computed as the inverse cube of the distance to the closest, nearest neighbor

class nbodykit.algorithms.FOF(source, linking_length, nmin, absolute=False, periodic=True)[source]

A friends-of-friends halo finder that computes the label for each particle, denoting which halo it belongs to.

Friends-of-friends was first used by Davis et al 1985 to define halos in hierachical structure formation of cosmological simulations. The algorithm is also known as DBSCAN in computer science. The subroutine here implements a parallel version of the FOF.

The underlying local FOF algorithm is from kdcount.cluster, which is an adaptation of the implementation in Volker Springel’s Gadget and Martin White’s PM.

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

For returning a CatalogSource of the FOF halos, see find_features() and for computing a halo catalog with added analytic information for a specific redshift and cosmology, see to_halos().

Parameters:
  • source (CatalogSource) – the source to run the FOF algorithm on; must support ‘Position’
  • linking_length (float) – the linking length, either in absolute units, or relative to the mean particle separation
  • nmin (int) – halo with fewer particles are ignored
  • absolute (bool, optional) – If True, the linking length is in absolute units, otherwise it is relative to the mean particle separation; default is False

Methods

find_features([peakcolumn]) Based on the particle labels, identify the groups, and return the center-of-mass CMPosition, CMVelocity, and Length of each feature.
run() Run the FOF algorithm.
to_halos(particle_mass, cosmo, redshift[, …]) Return a HaloCatalog, holding the center-of-mass position and velocity of each FOF halo, as well as the properly scaled mass, for a given cosmology and redshift.
find_features(peakcolumn=None)[source]

Based on the particle labels, identify the groups, and return the center-of-mass CMPosition, CMVelocity, and Length of each feature.

If a peakcolumn is given, the PeakPosition and PeakVelocity is also calculated for the particle at the peak value of the column.

Data is scattered evenly across all ranks.

Returns:a source holding the (‘CMPosition’, ‘CMVelocity’, ‘Length’) of each feature; optionaly, PeakPosition, PeakVelocity are also included if peakcolumn is not None
Return type:ArrayCatalog
run()[source]

Run the FOF algorithm. This function returns nothing, but does attach several attributes to the class instance:

  • attr:labels
  • max_labels

Note

The labels array is scattered evenly across all ranks.

labels

array_like, length: size – an array the label that specifies which FOF halo each particle belongs to

max_label

int – the maximum label across all ranks; this represents the total number of FOF halos found

to_halos(particle_mass, cosmo, redshift, mdef='vir', posdef='cm', peakcolumn='Density')[source]

Return a HaloCatalog, holding the center-of-mass position and velocity of each FOF halo, as well as the properly scaled mass, for a given cosmology and redshift.

The returned catalog also has default analytic prescriptions for halo radius and concentration.

The data is scattered evenly across all ranks.

Parameters:
  • particle_mass (float) – the particle mass, which is used to convert the number of particles in each halo to a total mass
  • cosmo (nbodykit.cosmology.core.Cosmology) – the cosmology of the catalog
  • redshift (float) – the redshift of the catalog
  • mdef (str, optional) – string specifying mass definition, used for computing default halo radii and concentration; should be ‘vir’ or ‘XXXc’ or ‘XXXm’ where ‘XXX’ is an int specifying the overdensity
  • posdef (str, optional) – position, can be cm (center of mass) or peak (particle with maximum value on a column)
  • peakcolumn (str , optional) – when posdef is ‘peak’, this is the column in source for identifying particles at the peak for the position and velocity.
Returns:

a HaloCatalog at the specified cosmology and redshift

Return type:

HaloCatalog

class nbodykit.algorithms.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

class nbodykit.algorithms.RedshiftHistogram(source, fsky, cosmo, bins=None, redshift='Redshift', weight=None)[source]

Compute the mean number density as a function of redshift \(n(z)\) from an input CatalogSource of particles.

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

Note

The units of the number density are \((\mathrm{Mpc}/h)^{-3}\)

Parameters:
  • source (CatalogSource) – the source of particles holding the redshift column to histogram
  • fsky (float) – the sky area fraction, which is used in the volume calculation when normalizing \(n(z)\)
  • cosmo (nbodykit.cosmology.core.Cosmology) – the cosmological parameters, which are used to compute the volume from redshift shells when normalizing \(n(z)\)
  • bins (int or sequence of scalars, optional) – If bins is an int, it defines the number of equal-width bins in the given range. If bins is a sequence, it defines the bin edges, including the rightmost edge, allowing for non-uniform bin widths. If not provided, Scott’s rule is used to estimate the optimal bin width from the input data (default)
  • redshift (str, optional) – the name of the column specifying the redshift data
  • weight (str, optional) – the name of the column specifying weights to use when histogramming the data

Methods

load(output[, comm]) Load a saved RedshiftHistogram result.
run() Run the algorithm, which computes the histogram.
save(output) Save the RedshiftHistogram result to disk.
classmethod load(output, comm=None)[source]

Load a saved RedshiftHistogram result.

The result has been saved to disk with RedshiftHistogram.save().

run()[source]

Run the algorithm, which computes the histogram. This function does not return anything, but adds the following attributes to the class:

Note

All ranks store the same result attributes.

bin_edges

array_like – the edges of the redshift bins

bin_centers

array_like – the center values of each redshift bin

dV

array_like – the volume of each redshift shell in units of \((\mathrm{Mpc}/h)^3\)

nbar

array_like – the values of the redshift histogram, normalized to number density (in units of \((\mathrm{Mpc}/h)^{-3}\))

save(output)[source]

Save the RedshiftHistogram result to disk.

The format is JSON.

class nbodykit.algorithms.FiberCollisions(ra, dec, collision_radius=0.017222222222222226, seed=None, degrees=True, comm=None)[source]

Run an angular FOF algorithm to determine fiber collision groups from an input catalog, and then assign fibers such that the maximum amount of object receive a fiber.

This amounts to determining the following population of objects:

  • population 1:
    the maximal “clean” sample of objects in which each object is not angularly collided with any other object in this subsample
  • population 2:
    the potentially-collided objects; these objects are those that are fiber collided + those that have been “resolved” due to multiple coverage in tile overlap regions

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

Parameters:
  • ra (array_like) – the right ascension coordinate column
  • dec (array_like) – the declination coordinate column
  • collision_radius (float, optional) – the size of the angular collision radius (in degrees); default is 62 arcseconds
  • seed (int, optional) – the random seed to use when determining which objects get fibers
  • degrees (bool, optional) – set to True if the units of ra and dec are degrees

References

Methods

run() Run the fiber assignment algorithm.
run()[source]

Run the fiber assignment algorithm. This attaches the following attribute to the object:

Note

The labels attribute has a 1-to-1 correspondence with the rows in the input source.

labels

ArrayCatalog; size: size – a CatalogSource that has the following columns:

  • Label :
    the group labels for each object in the input CatalogSource; label == 0 objects are not in a group
  • Collided :
    a flag array specifying which objects are collided, i.e., do not receive a fiber
  • NeighborID :
    for those objects that are collided, this gives the (global) index of the nearest neighbor on the sky (0-indexed) in the input catalog source, else it is set to -1
class nbodykit.algorithms.Multipoles3PCF(source, poles, edges, BoxSize=None, periodic=True, weight='Weight', selection='Selection')[source]

Compute the multipoles of the isotropic, three-point correlation function in configuration space.

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.

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
  • selection (str, optional) – if not None, the name of the column to use to select certain objects; should be the name of a boolean column

References

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

Methods

load(filename[, comm]) Load a Multipoles3PCF result from filename that has been saved to disk with save().
run() Compute the three-point CF multipoles.
save(output) Save the poles result to a JSON file with name output
classmethod load(filename, comm=None)[source]

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

run()[source]

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

poles

BinnedStatistic – a BinnedStatistic object to hold the multipole results; the binned statistics stores the multipoles as variables zeta_0, zeta_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

save(output)[source]

Save the poles result to a JSON file with name output

class nbodykit.algorithms.SimulationBoxPairCount(mode, first, redges, BoxSize=None, periodic=True, second=None, Nmu=5, los='z', weight='Weight', show_progress=True, **config)[source]

Count (weighted) pairs of objects in a simulation box using the Corrfunc package.

This uses the Corrfunc.theory.DD.DD() and Corrfunc.theory.DDsmu.DDsmu() functions to count pairs.

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 particles in a simulation box to be the Cartesian x, y, and z vectors. To compute pair counts on survey data, using right ascension, declination, and redshift, see SurveyDataPairCount.

Parameters:
  • mode ({'1d', '2d'}) – compute pair counts as a function of r and mu or just r
  • first (CatalogSource) – the first source of particles, providing the ‘Position’ column
  • redges (array_like) – the radius bin edges; length of nbins+1
  • BoxSize (float, 3-vector, optional) – the size of the box; if ‘BoxSize’ is not provided in the source ‘attrs’, it must be provided here
  • second (CatalogSource, optional) – the second source of particles to cross-correlate
  • Nmu (int, optional) – the number of mu bins, ranging from 0 to 1
  • los ({'x', 'y', 'z'}, int, optional) – the axis of the simulation box to treat as the line-of-sight direction; this can be provided as string identifying one of ‘x’, ‘y’, ‘z’ or the equivalent integer number of the axis
  • periodic (bool, optional) – whether to use periodic boundary conditions
  • weight (str, optional) – the name of the column in the source specifying the particle weights
  • show_progress (bool, optional) – if True, perform the pair counting calculation in 10 iterations, logging the progress after each iteration; this is useful for understanding the scaling of the code
  • **config (key/value pairs) – additional keywords to pass to the Corrfunc function

Methods

load(output[, comm]) Load a result has been saved to disk with save().
run() Calculate the 3D pair-counts in a simulation box as a function of separation r or separation and angle to line-of-sight (r, mu).
save(output) Save result as a JSON file with name output
run()[source]

Calculate the 3D pair-counts in a simulation box as a function of separation r or separation and angle to line-of-sight (r, mu). This adds the following attributes to the class:

result

BinnedStatistic – a BinnedStatistic object holding the pair count and correlation function results. The coordinate grid is either r or r and mu. It stores the following variables:

  • r: the mean separation value in the bin
  • xi: the mean correlation function value in the bin, computed as \(DD/RR - 1\), where \(RR\) is the number of random pairs in the bin
  • npairs: the number of pairs in the bin
  • weightavg: the average weight value in the bin; each pair contributes the product of the individual weight values
class nbodykit.algorithms.SurveyDataPairCount(mode, first, redges, cosmo, second=None, Nmu=5, ra='RA', dec='DEC', redshift='Redshift', weight='Weight', show_progress=True, **config)[source]

Count (weighted) pairs of objects from a survey data catalog using the Corrfunc package.

This uses the:func:Corrfunc.mocks.DDsmu_mocks.DDsmu_mocks function to count pairs.

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 particles from a survey catalog be the sky coordinates, right ascension and declination, and redshift. To compute pair counts in a simulation box, using the Cartesian coordinate vectors, see SimulationBoxPairCount.

Warning

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

Parameters:
  • mode ({'1d', '2d'}) – compute paircounts as a function of r and mu or just r
  • first (CatalogSource) – the first source of particles, providing the ‘Position’ column
  • redges (array_like) – the radius bin edges; length of nbins+1
  • cosmo (Cosmology) – the cosmology instance used to convert redshift into comoving distance
  • second (CatalogSource, optional) – the second source of particles to cross-correlate
  • Nmu (int, optional) – the number of mu bins, ranging from 0 to 1
  • 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
  • show_progress (bool, optional) – if True, perform the pair counting calculation in 10 iterations, logging the progress after each iteration; this is useful for understanding the scaling of the code
  • **config (key/value pairs) – additional keywords to pass to the Corrfunc function

Methods

load(output[, comm]) Load a result has been saved to disk with save().
run() Calculate the 3D pair-counts of a survey data catalog as a function of separation r or separation and angle to line-of-sight (r, mu).
save(output) Save result as a JSON file with name output
run()[source]

Calculate the 3D pair-counts of a survey data catalog as a function of separation r or separation and angle to line-of-sight (r, mu). This adds the following attribute:

result

BinnedStatistic – a BinnedStatistic object holding the pair count results. The coordinate grid is either r or r and mu. It stores the following variables:

  • r: the mean separation value in the bin
  • npairs: the number of pairs in the bin
  • weightavg: the average weight value in the bin; each pair contributes the product of the individual weight values
class nbodykit.algorithms.AngularPairCount(first, edges, second=None, ra='RA', dec='DEC', weight='Weight', show_progress=True, **config)[source]

Count (weighted) angular pairs of objects from a survey data catalog using the Corrfunc package.

This uses the:func:Corrfunc.mocks.DDtheta_mocks.DDtheta_mocks function to count pairs.

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 particles from the survey catalog be the right ascension and declination angular coordinates.

Warning

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

Parameters:
  • first (CatalogSource) – the first source of particles, providing the ‘Position’ column
  • edges (array_like) – the angular separation bin edges (in degrees); length of nbins+1
  • second (CatalogSource, optional) – the second source of particles to cross-correlate
  • 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’
  • weight (str, optional) – the name of the column in the source specifying the object weights
  • show_progress (bool, optional) – if True, perform the pair counting calculation in 10 iterations, logging the progress after each iteration; this is useful for understanding the scaling of the code
  • **config (key/value pairs) – additional keywords to pass to the Corrfunc function

Methods

load(output[, comm]) Load a result has been saved to disk with save().
run() Calculate the angular pair counts of a survey data catalog as a function of separation angle on the sky.
save(output) Save result as a JSON file with name output
run()[source]

Calculate the angular pair counts of a survey data catalog as a function of separation angle on the sky. This adds the following attribute:

result

BinnedStatistic – a BinnedStatistic object holding the pair count results. The coordinate grid is theta, the angular separation bins. It stores the following variables:

  • theta: the mean separation value in the bin
  • npairs: the number of pairs in the bin
  • weightavg: the average weight value in the bin; each pair contributes the product of the individual weight values
class nbodykit.algorithms.CylindricalGroups(source, rankby, rperp, rpar, flat_sky_los=None, periodic=False, BoxSize=None)[source]

Compute groups of objects using a cylindrical grouping method. We identify all satellites within a given cylindrical volume around a central object.

Results are computed when the object is inititalized, and the result is stored in the groups attribute; see the documenation of run().

Input parameters are stored in the attrs attribute dictionary.

Parameters:
  • source (subclass of CatalogSource) – the input source of particles providing the ‘Position’ column; the grouping algorithm is run on this catalog
  • rperp (float) – the radius of the cylinder in the sky plane (i.e., perpendicular to the line-of-sight)
  • rpar (float) – the radius along the line-of-sight direction; this is 1/2 the height of the cylinder
  • rankby (str, list, None) – a single or list of column names to rank order the input source by before computing the cylindrical groups, such that objects ranked first are marked as CGM centrals; if None is supplied, no sorting will be done
  • flat_sky_los (bool, optional) – a unit vector of length 3 providing the line-of-sight direction, assuming a fixed line-of-sight across the box, e.g., [0,0,1] to use the z-axis. If None, the observer at (0,0,0) is used to compute the line-of-sight for each pair
  • periodic (bool, optional) – whether to use periodic boundary conditions
  • BoxSize (float, 3-vector, optional) – the size of the box of the input data; must be provided as a keyword or in source.attrs if periodic=True

References

Okumura, Teppei, et al. “Reconstruction of halo power spectrum from redshift-space galaxy distribution: cylinder-grouping method and halo exclusion effect”, arXiv:1611.04165, 2016.

Methods

run() Compute the cylindrical groups, saving the results to the
run()[source]

Compute the cylindrical groups, saving the results to the groups attribute

groups

ArrayCatalog – a catalog holding the result of the grouping. The length of the catalog is equal to the length of the input size, i.e., the length is equal to the size attribute. The relevant fields are:

  1. cgm_type :
    a flag specifying the type for each object, with 0 specifying CGM central and 1 denoting CGM satellite
  2. cgm_haloid :
    The index of the CGM object this object belongs to; an integer between 0 and the total number of CGM halos
  3. num_cgm_sats :
    The number of satellites in the CGM halo