nbodykit.transform

Functions

ConcatenateSources(*sources, **kwargs) Concatenate CatalogSource objects together, optionally including only certain columns in the returned source.
ConstantArray(value, size[, chunks]) Return a dask array of the specified size holding a single value.
HaloConcentration(mass, cosmo, redshift[, mdef]) Return halo concentration from halo mass, based on the analytic fitting formulas presented in Dutton and Maccio 2014.
HaloRadius(mass, cosmo, redshift[, mdef]) Return halo radius from halo mass, based on the specified mass definition.
SkyToCartesian(ra, dec, redshift, cosmo[, …]) Convert sky coordinates (ra, dec, redshift) to a Cartesian Position column.
SkyToUnitSphere(ra, dec[, degrees]) Convert sky coordinates (ra, dec) to Cartesian coordinates on the unit sphere.
StackColumns(*cols) Stack the input dask arrays vertically, column by column.
nbodykit.transform.ConcatenateSources(*sources, **kwargs)[source]

Concatenate CatalogSource objects together, optionally including only certain columns in the returned source.

Note

The returned catalog object carries the meta-data from only the first catalog supplied to this function (in the attrs dict).

Parameters:
  • *sources (subclass of CatalogSource) – the catalog source objects to concatenate together
  • columns (str, list of str, optional) – the columns to include in the concatenated catalog
Returns:

the concatenated catalog source object

Return type:

CatalogSource

Examples

>>> from nbodykit.lab import *
>>> source1 = UniformCatalog(nbar=100, BoxSize=1.0)
>>> source2 = UniformCatalog(nbar=100, BoxSize=1.0)
>>> print(source1.csize, source2.csize)
>>> combined = transform.ConcatenateSources(source1, source2, columns=['Position', 'Velocity'])
>>> print(combined.csize)
nbodykit.transform.ConstantArray(value, size, chunks=100000)[source]

Return a dask array of the specified size holding a single value.

This uses numpy’s “stride tricks” to avoid replicating the data in memory for each element of the array.

Parameters:
  • value (float) – the scalar value to fill the array with
  • size (int) – the length of the returned dask array
  • chunks (int, optional) – the size of the dask array chunks
nbodykit.transform.HaloConcentration(mass, cosmo, redshift, mdef='vir')[source]

Return halo concentration from halo mass, based on the analytic fitting formulas presented in Dutton and Maccio 2014.

Note

The units of the input mass are assumed to be \(M_{\odot}/h\)

Parameters:
  • mass (array_like) – either a numpy or dask array specifying the halo mass; units assumed to be \(M_{\odot}/h\)
  • cosmo (Cosmology) – the cosmology instance used in the analytic formula
  • redshift (float) – compute the c(M) relation at this redshift
  • mdef (str, optional) – string specifying the halo mass definition to use; should be ‘vir’ or ‘XXXc’ or ‘XXXm’ where ‘XXX’ is an int specifying the overdensity
Returns:

concen – a dask array holding the analytic concentration values

Return type:

dask.array.Array

References

Dutton and Maccio, “Cold dark matter haloes in the Planck era: evolution of structural parameters for Einasto and NFW profiles”, 2014, arxiv:1402.7073

nbodykit.transform.HaloRadius(mass, cosmo, redshift, mdef='vir')[source]

Return halo radius from halo mass, based on the specified mass definition. This is independent of halo profile, and simply returns

\[R = \left [ 3 M /(4\pi\Delta) \right]^{1/3}\]

where \(\Delta\) is the density threshold, which depends on cosmology, redshift, and mass definition

Note

The units of the input mass are assumed to be \(M_{\odot}/h\)

Parameters:
  • mass (array_like) – either a numpy or dask array specifying the halo mass; units assumed to be \(M_{\odot}/h\)
  • cosmo (Cosmology) – the cosmology instance
  • redshift (float) – compute the density threshold which determines the R(M) relation at this redshift
  • mdef (str, optional) – string specifying the halo mass definition to use; should be ‘vir’ or ‘XXXc’ or ‘XXXm’ where ‘XXX’ is an int specifying the overdensity
Returns:

radius – a dask array holding the halo radius

Return type:

dask.array.Array

nbodykit.transform.SkyToCartesian(ra, dec, redshift, cosmo, degrees=True)[source]

Convert sky coordinates (ra, dec, redshift) to a Cartesian Position column.

Warning

The returned Cartesian position is in units of Mpc/h.

Parameters:
  • ra (dask.array.Array; shape: (N,)) – the right ascension angular coordinate
  • dec (dask.array.Array; shape: (N,)) – the declination angular coordinate
  • redshift (dask.array.Array; shape: (N,)) – the redshift coordinate
  • cosmo (Cosmology) – the cosmology used to meausre the comoving distance from redshift
  • degrees (bool, optional) – specifies whether ra and dec are in degrees
Returns:

pos – the cartesian position coordinates, where columns represent x, y, and z in units of Mpc/h

Return type:

dask.array.Array; shape: (N,3)

Raises:

TypeError – If the input columns are not dask arrays

nbodykit.transform.SkyToUnitSphere(ra, dec, degrees=True)[source]

Convert sky coordinates (ra, dec) to Cartesian coordinates on the unit sphere.

Parameters:
  • ra (dask.array.Array; shape: (N,)) – the right ascension angular coordinate
  • dec (dask.array.Array; ; shape: (N,)) – the declination angular coordinate
  • degrees (bool, optional) – specifies whether ra and dec are in degrees or radians
Returns:

pos – the cartesian position coordinates, where columns represent x, y, and z

Return type:

dask.array.Array; shape: (N,3)

Raises:

TypeError – If the input columns are not dask arrays

nbodykit.transform.StackColumns(*cols)[source]

Stack the input dask arrays vertically, column by column.

This uses dask.array.vstack().

Parameters:*cols (dask.array.Array) – the dask arrays to stack vertically together
Returns:the dask array where columns correspond to the input arrays
Return type:dask.array.Array
Raises:TypeError – If the input columns are not dask arrays