Reference

cosmotile

Cosmotile.

cosmotile.apply_rsds(field, los_displacement, distance, n_subcells=4)

Apply redshift-space distortions to a field.

Notes

To ensure that we cover all the slices in the field after the velocities have been applied, we extrapolate the densities and velocities on either end by the maximum velocity offset in the field. Then, to ensure we don’t pick up cells with zero particles (after displacement), we interpolate the slices onto a finer regular grid (in comoving distance) and then displace the field on that grid, before interpolating back onto the original slices.

Parameters:
  • field (ndarray) – The field to apply redshift-space distortions to, shape (nslices, ncoords).

  • los_displacement (ndarray) – The line-of-sight “apparent” displacement of the field, in pixel coordinates. Equal to v / H(z) / cell_size. Positive values are towards the observer, shape (nslices, ncoords).

  • distance (ndarray) – The comoving distance to each slice in the field, in units of the cell size. shape (nslices,).

  • n_subcells (int)

Return type:

ndarray

cosmotile.get_distance_to_shell_from_redshift(z, cell_size, cosmo=FlatLambdaCDM(name='Planck18', H0=<Quantity 67.66 km / (Mpc s)>, Om0=0.30966, Tcmb0=<Quantity 2.7255 K>, Neff=3.046, m_nu=<Quantity [0., 0., 0.06] eV>, Ob0=0.04897))

Get a distance to a shell, in units of cell size, from a given redshift.

Parameters:
  • z (float) – The redshift

  • cell_size (Quantity) – The resolution of the coeval simulation, in comoving units.

  • cosmo (FLRW) – The astropy cosmology.

Returns:

The distance, in units of pixels, to the shell.

Return type:

distance

cosmotile.make_healpix_lightcone_slice(nside, order='ring', **kwargs)

Create a healpix lightcone slice in angular coordinates.

This is a simple wrapper around make_lightcone_slice() that sets up angular co-ordinates from a healpix grid.

Parameters:
  • nside (int) – The Nside parameter of the healpix map.

  • order (Literal['ring', 'nested']) – The ordering of the pixels in the healpix map.

  • kwargs (Any)

Return type:

Generator

:param All other parameters are passed through to make_lightcone_slice().:

cosmotile.make_lightcone_slice(*, coevals, **kwargs)

Create a lightcone slice in angular coordinates from two coeval simulations.

Interpolates the input coeval box to angular coordinates.

Parameters:
  • coevals (Sequence[ndarray] | ndarray) – An iterable of rectangular coeval simulations to interpolate to the angular coordinates. Must have three dimensions (not necessarily the same size). Each box must have the same shape, and all are assumed to be at the same coordinates. Each coeval box can be a different simulated field.

  • kwargs (Any)

Return type:

Generator

:param All other parameters are passed to make_lightcone_slice_interpolator().:

Yields:

field – Each interpolated field on the angular coordinates.

Parameters:
  • coevals (Sequence[ndarray] | ndarray)

  • kwargs (Any)

Return type:

Generator

cosmotile.make_lightcone_slice_interpolator(*, latitude, longitude, distance_to_shell, interpolation_order=1, origin=None, rotation=None)

Create a callable interpolator for a lightcone slice.

Parameters:
  • latitude (ndarray) – An array of latitude coordinates onto which to tile the box. In radians from -pi/2 to pi/2

  • longitude (ndarray) – An array, same size as latitude, of longitude coordinates onto which to tile the box. In radians from 0 to 2pi.

  • distance_to_shell (float) – The distance to the spherical shell onto which to interpolate, in units of the cell-size of the coeval box(es) you wish to interpolate.

  • interpolation_order (int) – The order of interpolation. Must be in the range 0-5.

  • origin (ndarray | tuple[float, float, float] | None) – Define the location of the centre of the spherical shell, assuming that the (0,0,0) pixel of the coeval box is at (0,0,0) in cartesian coordinates.

  • rotation (Rotation | None) – The rotation by which to rotate the spherical coordinates before interpolation. This is done before shifting the origin, and is equivalent to rotating the coeval box beforing tiling it.

Returns:

A callable that takes a 3D array of coeval values and returns a 2D array of interpolated values on a redshift slice.

Return type:

interpolator

cosmotile.make_lightcone_slice_vector_field(coeval_vector_fields, interpolator)

Interpolate a 3D vector field to a lightcone slice as a line-of-sight component.

This takes a sequence of 3D vector fields, eg. the velocity field, and interpolates each component to the lightcone slice. It then computes the line-of-sight component of each interpolated vector field, where positive values are oriented towards the observer.

Parameters:
  • coeval_vector_fields (Sequence[Sequence[ndarray]]) – An iterable of 3D vector fields to interpolate to the lightcone slice. Each vector field must be an iterable of 3 3D arrays, each of the same shape.

  • interpolator (Callable[[ndarray], ndarray]) – A callable that takes a 3D array of coeval values and returns a 2D array of interpolated values on a redshift slice. This should be created by make_lightcone_slice_interpolator() using the properties of the coeval vector fields.

Yields:

los_component – The line-of-sight component of each interpolated vector field.

Return type:

Generator

cosmotile.transform_to_pixel_coords(*, comoving_radius, latitude, longitude, origin=None, rotation=None)

Transform input spherical coordinates to pixel coordinates wrt a coeval box.

Parameters:
  • comoving_radius (Quantity) – The radius of the spherical coordinates (in units of the cell size).

  • latitude (ndarray) – An array of latitude coordinates onto which to tile the box. In radians from -pi/2 to pi/2

  • longitude (ndarray) – An array, same size as latitude, of longitude coordinates onto which to tile the box. In radians from 0 to 2pi.

  • origin (Quantity | tuple[float, float, float] | None) – Define the location of the centre of the spherical shell, assuming that the (0,0,0) pixel of the coeval box is at (0,0,0) in cartesian coordinates. In units of the cell size.

  • rotation (Rotation | None) – The rotation by which to rotate the spherical coordinates before interpolation. This is done before shifting the origin, and is equivalent to rotating the coeval box beforing tiling it.

Return type:

ndarray