yt.utilities.cosmology module

class yt.utilities.cosmology.Cosmology(hubble_constant=0.71, omega_matter=0.27, omega_lambda=0.73, omega_radiation=0.0, omega_curvature=0.0, unit_registry=None, unit_system='cgs', use_dark_factor=False, w_0=-1.0, w_a=0.0)[source]

Bases: object

Create a cosmology calculator to compute cosmological distances and times.

For an explanation of the various cosmological measures, see, for example Hogg (1999, https://arxiv.org/abs/astro-ph/9905116).

WARNING: Cosmological distance calculations return values that are either in the comoving or proper frame, depending on the specific quantity. For simplicity, the proper and comoving frames are set equal to each other within the cosmology calculator. This means that for some distance value, x, x.to(“Mpc”) and x.to(“Mpccm”) will be the same. The user should take care to understand which reference frame is correct for the given calculation.

Parameters:
  • hubble_constant (float) – The Hubble parameter at redshift zero in units of 100 km/s/Mpc. Default: 0.71.

  • omega_matter (the fraction of the energy density of the Universe in) – matter at redshift zero. Default: 0.27.

  • omega_lambda (the fraction of the energy density of the Universe in) – a cosmological constant. Default: 0.73.

  • omega_radiation (the fraction of the energy density of the Universe in) – relativistic matter at redshift zero.

  • omega_curvature (the fraction of the energy density of the Universe in) – curvature. Default: 0.0.

  • unit_system (yt.units.unit_systems.UnitSystem, optional) – The units system to use when making calculations. If not specified, cgs units are assumed.

  • use_dark_factor (Bool, optional) – The flag to either use the cosmological constant (False, default) or to use the parameterization of w(a) as given in Linder 2002. This, along with w_0 and w_a, only matters in the function expansion_factor.

  • w_0 (float, optional) – The Linder 2002 parameterization of w(a) is: w(a) = w_0 + w_a(1 - a). w_0 is w(a = 1). Only matters if use_dark_factor = True. Default is None. Cosmological constant case corresponds to w_0 = -1.

  • w_a (float, optional) – See w_0. w_a is the derivative of w(a) evaluated at a = 1. Cosmological constant case corresponds to w_a = 0. Default is None.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.t_from_z(0.0).in_units("Gyr"))
a_from_t(t)[source]

Compute the scale factor for a given age of the Universe.

Parameters:

t (YTQuantity or float) – Time since the Big Bang. If a float is given, units are assumed to be seconds.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.a_from_t(4.0e17))
age_integrand(z)[source]
angular_diameter_distance(z_i, z_f)[source]

Following Hogg (1999), the angular diameter distance is ‘the ratio of an object’s physical transverse size to its angular size in radians.’

Parameters:
  • z_i (float) – The redshift of the observer.

  • z_f (float) – The redshift of the observed object.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.angular_diameter_distance(0.0, 1.0).in_units("Mpc"))
angular_scale(z_i, z_f)[source]

The proper transverse distance between two points at redshift z_f observed at redshift z_i per unit of angular separation.

Parameters:
  • z_i (float) – The redshift of the observer.

  • z_f (float) – The redshift of the observed object.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.angular_scale(0.0, 1.0).in_units("kpc / arcsec"))
property arr
comoving_radial_distance(z_i, z_f)[source]

The comoving distance along the line of sight to on object at redshift, z_f, viewed at a redshift, z_i.

Parameters:
  • z_i (float) – The redshift of the observer.

  • z_f (float) – The redshift of the observed object.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.comoving_radial_distance(0.0, 1.0).in_units("Mpccm"))
comoving_transverse_distance(z_i, z_f)[source]

When multiplied by some angle, the distance between two objects observed at redshift, z_f, with an angular separation given by that angle, viewed by an observer at redshift, z_i (Hogg 1999).

Parameters:
  • z_i (float) – The redshift of the observer.

  • z_f (float) – The redshift of the observed object.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.comoving_transverse_distance(0.0, 1.0).in_units("Mpccm"))
comoving_volume(z_i, z_f)[source]

“The comoving volume is the volume measure in which number densities of non-evolving objects locked into Hubble flow are constant with redshift.” – Hogg (1999)

Parameters:
  • z_i (float) – The lower redshift of the interval.

  • z_f (float) – The higher redshift of the interval.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.comoving_volume(0.0, 1.0).in_units("Gpccm**3"))
critical_density(z)[source]

The density required for closure of the Universe at a given redshift in the proper frame.

Parameters:

z (float) – Redshift.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.critical_density(0.0).in_units("g/cm**3"))
>>> print(co.critical_density(0).in_units("Msun/Mpc**3"))
expansion_factor(z)[source]

The ratio between the Hubble parameter at a given redshift and redshift zero.

This is also the primary function integrated to calculate the cosmological distances.

get_dark_factor(z)[source]

This function computes the additional term that enters the expansion factor when using non-standard dark energy. See Dolag et al 2004 eq. 7 for ref (but note that there’s a typo in his eq. There should be no negative sign).

At the moment, this only works using the parameterization given in Linder 2002 eq. 7: w(a) = w0 + wa(1 - a) = w0 + wa * z / (1+z). This gives rise to an analytic expression. It is also only functional for Gadget simulations, at the moment.

Parameters:

z (float) – Redshift

hubble_distance()[source]

The distance corresponding to c / h, where c is the speed of light and h is the Hubble parameter in units of 1 / time.

hubble_parameter(z)[source]

The value of the Hubble parameter at a given redshift.

Parameters:

z (float) – Redshift.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.hubble_parameter(1.0).in_units("km/s/Mpc"))
inverse_expansion_factor(z)[source]
lookback_time(z_i, z_f)[source]

The difference in the age of the Universe between the redshift interval z_i to z_f.

Parameters:
  • z_i (float) – The lower redshift of the interval.

  • z_f (float) – The higher redshift of the interval.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.lookback_time(0.0, 1.0).in_units("Gyr"))
luminosity_distance(z_i, z_f)[source]

The distance that would be inferred from the inverse-square law of light and the measured flux and luminosity of the observed object.

Parameters:
  • z_i (float) – The redshift of the observer.

  • z_f (float) – The redshift of the observed object.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.luminosity_distance(0.0, 1.0).in_units("Mpc"))
path_length(z_i, z_f)[source]
path_length_function(z)[source]
property quan
t_from_a(a)[source]

Compute the age of the Universe for a given scale factor.

Parameters:

a (float) – Scale factor.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.t_from_a(1.0).in_units("Gyr"))
t_from_z(z)[source]

Compute the age of the Universe for a given redshift.

Parameters:

z (float) – Redshift.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.t_from_z(0.0).in_units("Gyr"))
z_from_t(t)[source]

Compute the redshift for a given age of the Universe.

Parameters:

t (YTQuantity or float) – Time since the Big Bang. If a float is given, units are assumed to be seconds.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.z_from_t(4.0e17))
class yt.utilities.cosmology.InterpTable(x, y)[source]

Bases: object

Generate a function to linearly interpolate from provided arrays.

yt.utilities.cosmology.trapezoid_cumulative_integral(f, x)[source]

Perform cumulative integration using the trapezoid rule.

yt.utilities.cosmology.trapezoid_int(f, a, b, bins=10000)[source]
yt.utilities.cosmology.trapzint(f, a, b, bins=10000)[source]