py21cmfast.AngularLightconer

class py21cmfast.AngularLightconer

Bases: Lightconer

Angular lightcone slices constructed from rectlinear coevals.

__attrs_post_init__()

Post-init.

Return type:

None

classmethod __init_subclass__()

Enabe plugin-style behaviour.

Return type:

None

classmethod between_redshifts(min_redshift, max_redshift, resolution, cosmo=Planck18, **kw)

Construct a Lightconer with regular comoving dist. slices between two z’s.

Parameters:
  • min_redshift (float)

  • max_redshift (float)

  • resolution (astropy.units.Quantity[_LENGTH])

coeval_subselect(lcd, coeval, coeval_res)

Sub-Select the coeval box required for interpolation at one slice.

Parameters:
  • lcd (float)

  • coeval (numpy.ndarray)

  • coeval_res (astropy.units.Quantity[_LENGTH])

Return type:

numpy.ndarray

construct_lightcone(lcd, box)

Construct the lightcone slices from bracketing coevals.

Parameters:
  • lcd (astropy.units.Quantity[astropy.units.pixel])

  • box (numpy.ndarray)

Return type:

tuple[numpy.ndarray, numpy.ndarray]

construct_los_velocity_lightcone(lcd, velocities)

Construct the LoS velocity lightcone from 3D velocities.

Parameters:
  • lcd (astropy.units.Quantity[astropy.units.pixel])

  • velocities (collections.abc.Sequence[numpy.ndarray])

find_required_lightcone_limits(classy_output, inputs)

Obtain the redshift limits required for the lightcone to include RSDs.

This is a crude estimation of the maximum/minimum lightcone limits that are required in order to simulate all the “mass” that enters the requested ligthcone (due to RSD shift). We use the rms of the velocity field from linear perturbation theory in order to determine the required lightcone limits. If no limit is found, it means that the limits of node_redshifts are not sufficient and an error is raised.

Parameters:
  • classy_output (classy.Class) – An object containing all the information from the CLASS calculation.

  • inputs (InputParameters) – The input parameters corresponding to the box.

Returns:

lcd_limits_rsd (list) – List that contains the limits of the required lightcone distances.

Parameters:
Return type:

list[astropy.units.Quantity]

get_lc_distances_in_pixels(resolution)

Get the lightcone distances in pixels, given a resolution.

Parameters:

resolution (astropy.units.Quantity[_LENGTH])

get_shape(simulation_options)

Get the shape of the lightcone slices.

Parameters:

simulation_options (py21cmfast.wrapper.inputs.SimulationOptions)

Return type:

tuple[int, int]

classmethod like_rectilinear(simulation_options, match_at_z, cosmo=Planck18, **kw)

Create an angular lightconer with the same pixel size as a rectilinear one.

This is useful for comparing the two lightconer types.

Parameters:
  • simulation_options – The user parameters.

  • match_at_z – The redshift at which the angular lightconer should match the rectilinear one.

  • cosmo – The cosmology to use.

Other Parameters:

All other parameters passed through to the constructor.

Returns:

AngularLightconer – The angular lightconer.

Parameters:
make_lightcone_slices(c1, c2)

Make lightcone slices out of two coeval objects.

Parameters:

c1, c2 (Coeval) – The coeval boxes to interpolate.

Returns:

  • quantity – The field names of the quantities required by the lightcone.

  • lcidx – The indices of the lightcone to which these slices belong.

  • scalar_field_slices – The scalar fields evaluated on the “lightcone” slices that exist within the redshift range spanned by c1 and c2.

Parameters:
Return type:

tuple[dict[str, numpy.ndarray], numpy.ndarray]

redshift_interpolation(dc, coeval_a, coeval_b, dc_a, dc_b, kind='mean')

Perform redshift interpolation to a new box given two bracketing coevals.

Parameters:
  • dc (float)

  • coeval_a (numpy.ndarray)

  • coeval_b (numpy.ndarray)

  • dc_a (float)

  • dc_b (float)

  • kind (str)

Return type:

numpy.ndarray

validate_options(inputs, include_dvdr_in_tau21, apply_rsds, classy_output=None)

Validate 21cmFAST options.

Parameters:
Return type:

Lightconer

classmethod with_equal_cdist_slices(min_redshift, max_redshift, resolution, cosmo=Planck18, **kw)

Create a lightconer with equally spaced slices in comoving distance.

This method is deprecated and will be removed in future versions. Instead, use between_redshifts to create a lightconer with equally spaced slices in comoving distance.

Parameters:
  • min_redshift (float)

  • max_redshift (float)

  • resolution (astropy.units.Quantity[_LENGTH])

__slots__ = ()
cosmo: astropy.cosmology.FLRW
interp_kinds: dict[str, str]
interpolation_order: int
latitude: numpy.ndarray
lc_distances: astropy.units.Quantity[_LENGTH]
property lc_redshifts: numpy.ndarray

The redshifts of the lightcone slices.

Return type:

numpy.ndarray

longitude: numpy.ndarray
origin: astropy.units.Quantity[astropy.units.pixel, 3, float]
quantities: tuple[str]
rotation: scipy.spatial.transform.Rotation