FlatLambdaCDM#

class astropy.cosmology.FlatLambdaCDM(H0, Om0, Tcmb0=<Quantity 0. K>, Neff=3.04, m_nu=<Quantity 0. eV>, Ob0=None, *, name=None, meta=None)[source]#

Bases: FlatFLRWMixin, LambdaCDM

FLRW cosmology with a cosmological constant and no curvature.

This has no additional attributes beyond those of FLRW.

Parameters:
H0float or scalar astropy:quantity-like [:ref: ‘frequency’]

Hubble constant at z = 0. If a float, must be in [km/sec/Mpc].

Om0float

Omega matter: density of non-relativistic matter in units of the critical density at z=0.

Tcmb0float or scalar astropy:quantity-like [:ref: ‘temperature’], optional

Temperature of the CMB z=0. If a float, must be in [K]. Default: 0 [K]. Setting this to zero will turn off both photons and neutrinos (even massive ones).

Nefffloat, optional

Effective number of Neutrino species. Default 3.04.

m_nuastropy:quantity-like [:ref: ‘energy’, :ref: ‘mass’] or array_like, optional

Mass of each neutrino species in [eV] (mass-energy equivalency enabled). If this is a scalar Quantity, then all neutrino species are assumed to have that mass. Otherwise, the mass of each species. The actual number of neutrino species (and hence the number of elements of m_nu if it is not scalar) must be the floor of Neff. Typically this means you should provide three neutrino masses unless you are considering something like a sterile neutrino.

Ob0float or None, optional

Omega baryons: density of baryonic matter in units of the critical density at z=0. If this is set to None (the default), any computation that requires its value will raise an exception.

namestr or None (optional, keyword-only)

Name for this cosmological object.

metamapping or None (optional, keyword-only)

Metadata for the cosmology, e.g., a reference.

Examples

>>> from astropy.cosmology import FlatLambdaCDM
>>> cosmo = FlatLambdaCDM(H0=70, Om0=0.3)

The comoving distance in Mpc at redshift z:

>>> z = 0.5
>>> dc = cosmo.comoving_distance(z)

To get an equivalent cosmology, but of type astropy.cosmology.LambdaCDM, use astropy.cosmology.FlatFLRWMixin.nonflat.

>>> print(cosmo.nonflat)
LambdaCDM(H0=70.0 km / (Mpc s), Om0=0.3, Ode0=0.7, ...

Attributes Summary

H0

Hubble constant as an Quantity at z=0.

Neff

Number of effective neutrino species.

Ob0

Omega baryon; baryonic matter density/critical density at z=0.

Ode0

Omega dark energy; dark energy density/critical density at z=0.

Odm0

Omega dark matter; dark matter density/critical density at z=0.

Ogamma0

Omega gamma; the density/critical density of photons at z=0.

Ok0

Omega curvature; the effective curvature density/critical density at z=0.

Om0

Omega matter; matter density/critical density at z=0.

Onu0

Omega nu; the density/critical density of neutrinos at z=0.

Otot0

Omega total; the total density/critical density at z=0.

Tcmb0

Temperature of the CMB as Quantity at z=0.

Tnu0

Temperature of the neutrino background as Quantity at z=0.

critical_density0

Critical density as Quantity at z=0.

from_format

Transform object to a Cosmology.

h

Dimensionless Hubble constant: h = H_0 / 100 [km/sec/Mpc].

has_massive_nu

Does this cosmology have at least one massive neutrino species?

hubble_distance

Hubble distance as Quantity.

hubble_time

Hubble time as Quantity.

is_flat

Return True, the cosmology is flat.

m_nu

Mass of neutrino species.

meta

name

The name of the Cosmology instance.

nonflat

Return the equivalent non-flat-class instance of this cosmology.

parameters

Immutable mapping of the Parameters.

read

Read and parse data to a Cosmology.

scale_factor0

Scale factor at redshift 0.

to_format

Transform this Cosmology to another format.

write

Write this Cosmology object out in the specified format.

Methods Summary

H(z)

Hubble parameter (km/s/Mpc) at redshift z.

Ob(z)

Return the density parameter for baryonic matter at redshift z.

Ode(z)

Return the density parameter for dark energy at redshift z.

Odm(z)

Return the density parameter for dark matter at redshift z.

Ogamma(z)

Return the density parameter for photons at redshift z.

Ok(z)

Return the equivalent density parameter for curvature at redshift z.

Om(z)

Return the density parameter for non-relativistic matter at redshift z.

Onu(z)

Return the density parameter for neutrinos at redshift z.

Otot(z)

The total density parameter at redshift z.

Tcmb(z)

Return the CMB temperature at redshift z.

Tnu(z)

Return the neutrino temperature at redshift z.

abs_distance_integrand(z)

Integrand of the absorption distance [1].

absorption_distance(z, /)

Absorption distance at redshift z.

age(z)

Age of the universe in Gyr at redshift z.

angular_diameter_distance(z)

Angular diameter distance in Mpc at a given redshift.

angular_diameter_distance_z1z2(z1, z2)

Angular diameter distance between objects at 2 redshifts.

arcsec_per_kpc_comoving(z)

Angular separation in arcsec equal to a comoving kpc at redshift z.

arcsec_per_kpc_proper(z)

Angular separation in arcsec corresponding to a proper kpc at redshift z.

clone(*[, meta, to_nonflat])

Returns a copy of this object with updated parameters, as specified.

comoving_distance(z)

Comoving line-of-sight distance in Mpc at a given redshift.

comoving_transverse_distance(z)

Comoving transverse distance in Mpc at a given redshift.

comoving_volume(z)

Comoving volume in cubic Mpc at redshift z.

critical_density(z)

Critical density in grams per cubic cm at redshift z.

de_density_scale(z)

Evaluates the redshift dependence of the dark energy density.

differential_comoving_volume(z)

Differential comoving volume at redshift z.

distmod(z)

Distance modulus at redshift z.

efunc(z)

Function used to calculate H(z), the Hubble parameter.

inv_efunc(z)

Function used to calculate \(\frac{1}{H_z}\).

is_equivalent(other, /, *[, format])

Check equivalence between Cosmologies.

kpc_comoving_per_arcmin(z)

Separation in transverse comoving kpc equal to an arcmin at redshift z.

kpc_proper_per_arcmin(z)

Separation in transverse proper kpc equal to an arcminute at redshift z.

lookback_distance(z)

The lookback distance is the light travel time distance to a given redshift.

lookback_time(z)

Lookback time in Gyr to redshift z.

lookback_time_integrand(z)

Integrand of the lookback time (equation 30 of [1]).

luminosity_distance(z)

Luminosity distance in Mpc at redshift z.

nu_relative_density(z)

Neutrino density function relative to the energy density in photons.

scale_factor(z)

Scale factor at redshift z.

w(z)

Returns dark energy equation of state at redshift z.

Attributes Documentation

H0#

Hubble constant as an Quantity at z=0.

Neff#

Number of effective neutrino species.

Ob0#

Omega baryon; baryonic matter density/critical density at z=0.

Ode0#

Omega dark energy; dark energy density/critical density at z=0.

Odm0#

Omega dark matter; dark matter density/critical density at z=0.

Ogamma0#

Omega gamma; the density/critical density of photons at z=0.

Ok0#

Omega curvature; the effective curvature density/critical density at z=0.

Om0#

Omega matter; matter density/critical density at z=0.

Onu0#

Omega nu; the density/critical density of neutrinos at z=0.

Otot0#

Omega total; the total density/critical density at z=0.

Tcmb0#

Temperature of the CMB as Quantity at z=0.

Tnu0#

Temperature of the neutrino background as Quantity at z=0.

critical_density0#

Critical density as Quantity at z=0.

from_format#

Transform object to a Cosmology.

This function provides the Cosmology interface to the Astropy unified I/O layer. This allows easily parsing supported data formats using syntax such as:

>>> from astropy.cosmology import Cosmology
>>> cosmo1 = Cosmology.from_format(cosmo_mapping, format='mapping')

When the from_format method is called from a subclass the subclass will provide a keyword argument cosmology=<class> to the registered parser. The method uses this cosmology class, regardless of the class indicated in the data, and sets parameters’ default values from the class’ signature.

Get help on the available readers using the help() method:

>>> Cosmology.from_format.help()  # Get help and list supported formats
>>> Cosmology.from_format.help('<format>')  # Get detailed help on a format
>>> Cosmology.from_format.list_formats()  # Print list of available formats

See also: https://docs.astropy.org/en/stable/io/unified.html

Parameters:
objobject

The object to parse according to ‘format’

*args

Positional arguments passed through to data parser.

formatstr or None, optional keyword-only

Object format specifier. For None (default) CosmologyFromFormat tries to identify the correct format.

**kwargs

Keyword arguments passed through to data parser. Parsers should accept the following keyword arguments:

  • cosmologythe class (or string name thereof) to use / check when

    constructing the cosmology instance.

Returns:
outCosmology subclass instance

Cosmology corresponding to obj contents.

h#

Dimensionless Hubble constant: h = H_0 / 100 [km/sec/Mpc].

has_massive_nu#

Does this cosmology have at least one massive neutrino species?

hubble_distance#

Hubble distance as Quantity.

hubble_time#

Hubble time as Quantity.

is_flat#

Return True, the cosmology is flat.

m_nu#

Mass of neutrino species.

meta = None#
name#

The name of the Cosmology instance.

nonflat#
parameters = mappingproxy({'H0': Parameter(derived=False, unit=Unit("km / (Mpc s)"), equivalencies=[], fvalidate='scalar', doc='Hubble constant as an `~astropy.units.Quantity` at z=0.'), 'Om0': Parameter(derived=False, unit=None, equivalencies=[], fvalidate='non-negative', doc='Omega matter; matter density/critical density at z=0.'), 'Tcmb0': Parameter(default=<Quantity 0. K>, derived=False, unit=Unit("K"), equivalencies=[], fvalidate='scalar', doc='Temperature of the CMB as `~astropy.units.Quantity` at z=0.'), 'Neff': Parameter(default=3.04, derived=False, unit=None, equivalencies=[], fvalidate='non-negative', doc='Number of effective neutrino species.'), 'm_nu': Parameter(default=<Quantity 0. eV>, derived=False, unit=Unit("eV"), equivalencies=[(Unit("kg"), Unit("J"), <function mass_energy.<locals>.<lambda>>, <function mass_energy.<locals>.<lambda>>), (Unit("kg / m2"), Unit("J / m2"), <function mass_energy.<locals>.<lambda>>, <function mass_energy.<locals>.<lambda>>), (Unit("kg / m3"), Unit("J / m3"), <function mass_energy.<locals>.<lambda>>, <function mass_energy.<locals>.<lambda>>), (Unit("kg / s"), Unit("J / s"), <function mass_energy.<locals>.<lambda>>, <function mass_energy.<locals>.<lambda>>)], fvalidate=<function FLRW.m_nu>, doc='Mass of neutrino species.'), 'Ob0': Parameter(default=None, derived=False, unit=None, equivalencies=[], fvalidate=<function FLRW.Ob0>, doc='Omega baryon; baryonic matter density/critical density at z=0.')})#

Immutable mapping of the Parameters.

If accessed from the class, this returns a mapping of the Parameter objects themselves. If accessed from an instance, this returns a mapping of the values of the Parameters.

read#

Read and parse data to a Cosmology.

This function provides the Cosmology interface to the Astropy unified I/O layer. This allows easily reading a file in supported data formats using syntax such as:

>>> from astropy.cosmology import Cosmology
>>> cosmo1 = Cosmology.read('<file name>')

When the read method is called from a subclass the subclass will provide a keyword argument cosmology=<class> to the registered read method. The method uses this cosmology class, regardless of the class indicated in the file, and sets parameters’ default values from the class’ signature.

Get help on the available readers using the help() method:

>>> Cosmology.read.help()  # Get help reading and list supported formats
>>> Cosmology.read.help(format='<format>')  # Get detailed help on a format
>>> Cosmology.read.list_formats()  # Print list of available formats

See also: https://docs.astropy.org/en/stable/io/unified.html

Parameters:
*args

Positional arguments passed through to data reader. If supplied the first argument is typically the input filename.

formatstr (optional, keyword-only)

File format specifier.

**kwargs

Keyword arguments passed through to data reader.

Returns:
outCosmology subclass instance

Cosmology corresponding to file contents.

scale_factor0#

Scale factor at redshift 0.

The scale factor is defined as \(a = \frac{a_0}{1 + z}\). The common convention is to set \(a_0 = 1\). However, in some cases, e.g. in some old CMB papers, \(a_0\) is used to normalize a to be a convenient number at the redshift of interest for that paper. Explicitly using \(a_0\) in both calculation and code avoids ambiguity.

to_format#

Transform this Cosmology to another format.

This function provides the Cosmology interface to the astropy unified I/O layer. This allows easily transforming to supported data formats using syntax such as:

>>> from astropy.cosmology import Planck18
>>> Planck18.to_format("mapping")
{'cosmology': astropy.cosmology.core.FlatLambdaCDM,
 'name': 'Planck18',
 'H0': <Quantity 67.66 km / (Mpc s)>,
 'Om0': 0.30966,
 ...

Get help on the available representations for Cosmology using the help() method:

>>> Cosmology.to_format.help()  # Get help and list supported formats
>>> Cosmology.to_format.help('<format>')  # Get detailed help on format
>>> Cosmology.to_format.list_formats()  # Print list of available formats
Parameters:
formatstr

Format specifier.

*args

Positional arguments passed through to data writer. If supplied the first argument is the output filename.

**kwargs

Keyword arguments passed through to data writer.

write#

Write this Cosmology object out in the specified format.

This function provides the Cosmology interface to the astropy unified I/O layer. This allows easily writing a file in supported data formats using syntax such as:

>>> from astropy.cosmology import Planck18
>>> Planck18.write('<file name>')

Get help on the available writers for Cosmology using the help() method:

>>> Cosmology.write.help()  # Get help writing and list supported formats
>>> Cosmology.write.help(format='<format>')  # Get detailed help on format
>>> Cosmology.write.list_formats()  # Print list of available formats
Parameters:
*args

Positional arguments passed through to data writer. If supplied the first argument is the output filename.

formatstr (optional, keyword-only)

File format specifier.

**kwargs

Keyword arguments passed through to data writer.

Methods Documentation

H(z)#

Hubble parameter (km/s/Mpc) at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
HQuantity [:ref: ‘frequency’]

Hubble parameter at each input redshift.

Ob(z)#

Return the density parameter for baryonic matter at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Obndarray or float

The density of baryonic matter relative to the critical density at each redshift. Returns float if the input is scalar.

Raises:
ValueError

If Ob0 is None.

Ode(z)#

Return the density parameter for dark energy at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Odendarray or float

The density of non-relativistic matter relative to the critical density at each redshift. Returns float if the input is scalar.

Odm(z)#

Return the density parameter for dark matter at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Odmndarray or float

The density of non-relativistic dark matter relative to the critical density at each redshift. Returns float if the input is scalar.

Raises:
ValueError

If Ob0 is None.

Notes

This does not include neutrinos, even if non-relativistic at the redshift of interest.

Ogamma(z)#

Return the density parameter for photons at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Ogammandarray or float

The energy density of photons relative to the critical density at each redshift. Returns float if the input is scalar.

Ok(z)#

Return the equivalent density parameter for curvature at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Okndarray or float

The equivalent density parameter for curvature at each redshift. Returns float if the input is scalar.

Om(z)#

Return the density parameter for non-relativistic matter at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Omndarray or float

The density of non-relativistic matter relative to the critical density at each redshift. Returns float if the input is scalar.

Notes

This does not include neutrinos, even if non-relativistic at the redshift of interest; see Onu.

Onu(z)#

Return the density parameter for neutrinos at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Onundarray or float

The energy density of neutrinos relative to the critical density at each redshift. Note that this includes their kinetic energy (if they have mass), so it is not equal to the commonly used \(\sum \frac{m_{\nu}}{94 eV}\), which does not include kinetic energy. Returns float if the input is scalar.

Otot(z)#

The total density parameter at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshifts.

Returns:
Ototndarray or float

Returns float if input scalar. Value of 1.

Tcmb(z)#

Return the CMB temperature at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
TcmbQuantity [:ref: ‘temperature’]

The temperature of the CMB in K.

Tnu(z)#

Return the neutrino temperature at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
TnuQuantity [:ref: ‘temperature’]

The temperature of the cosmic neutrino background in K.

abs_distance_integrand(z)#

Integrand of the absorption distance [1].

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Xfloat or array

The integrand for the absorption distance.

References

[1] (1,2)

Hogg, D. (1999). Distance measures in cosmology, section 11. arXiv e-prints, astro-ph/9905116.

absorption_distance(z, /)#

Absorption distance at redshift z.

This is used to calculate the number of objects with some cross section of absorption and number density intersecting a sightline per unit redshift path ([1], [2]).

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
dfloat or ndarray

Absorption distance (dimensionless) at each input redshift. Returns float if input scalar, ndarray otherwise.

References

[1]

Hogg, D. (1999). Distance measures in cosmology, section 11. arXiv e-prints, astro-ph/9905116.

[2]

Bahcall, John N. and Peebles, P.J.E. 1969, ApJ, 156L, 7B

age(z)#

Age of the universe in Gyr at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
tQuantity [:ref: ‘time’]

The age of the universe in Gyr at each input redshift.

See also

z_at_value

Find the redshift corresponding to an age.

angular_diameter_distance(z)#

Angular diameter distance in Mpc at a given redshift.

This gives the proper (sometimes called ‘physical’) transverse distance corresponding to an angle of 1 radian for an object at redshift z ([1], [2], [3]).

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
dQuantity [:ref: ‘length’]

Angular diameter distance in Mpc at each input redshift.

References

[1]

Weinberg, 1972, pp 420-424; Weedman, 1986, pp 421-424.

[2]

Weedman, D. (1986). Quasar astronomy, pp 65-67.

[3]

Peebles, P. (1993). Principles of Physical Cosmology, pp 325-327.

angular_diameter_distance_z1z2(z1, z2)#

Angular diameter distance between objects at 2 redshifts.

Useful for gravitational lensing, for example computing the angular diameter distance between a lensed galaxy and the foreground lens.

Parameters:
z1, z2Quantity-like [‘redshift’], array_like, or Number

Input redshifts. For most practical applications such as gravitational lensing, z2 should be larger than z1. The method will work for z2 < z1; however, this will return negative distances.

Returns:
dQuantity

The angular diameter distance between each input redshift pair. Returns scalar if input is scalar, array else-wise.

arcsec_per_kpc_comoving(z)#

Angular separation in arcsec equal to a comoving kpc at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
thetaQuantity [:ref: ‘angle’]

The angular separation in arcsec corresponding to a comoving kpc at each input redshift.

arcsec_per_kpc_proper(z)#

Angular separation in arcsec corresponding to a proper kpc at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
thetaQuantity [:ref: ‘angle’]

The angular separation in arcsec corresponding to a proper kpc at each input redshift.

clone(*, meta: Mapping | None = None, to_nonflat: bool = False, **kwargs)#

Returns a copy of this object with updated parameters, as specified.

This cannot be used to change the type of the cosmology, except for changing to the non-flat version of this cosmology.

Parameters:
metamapping or None (optional, keyword-only)

Metadata that will update the current metadata.

to_nonflatbool, optional keyword-only

Whether to change to the non-flat version of this cosmology.

**kwargs

Cosmology parameter (and name) modifications. If any parameter is changed and a new name is not given, the name will be set to “[old name] (modified)”.

Returns:
newcosmoCosmology subclass instance

A new instance of this class with updated parameters as specified. If no arguments are given, then a reference to this object is returned instead of copy.

Examples

To make a copy of the Planck13 cosmology with a different matter density (Om0), and a new name:

>>> from astropy.cosmology import Planck13
>>> Planck13.clone(name="Modified Planck 2013", Om0=0.35)
FlatLambdaCDM(name='Modified Planck 2013', H0=<Quantity 67.77 km / (Mpc s)>,
              Om0=0.35, ...

If no name is specified, the new name will note the modification.

>>> Planck13.clone(Om0=0.35).name
'Planck13 (modified)'

The keyword ‘to_nonflat’ can be used to clone on the non-flat equivalent cosmology. For FLRW cosmologies this means Ode0 can be modified:

>>> Planck13.clone(to_nonflat=True, Ode0=1)
LambdaCDM(name='Planck13 (modified)', H0=<Quantity 67.77 km / (Mpc s)>,
          Om0=0.30712, Ode0=1.0, ...
comoving_distance(z)#

Comoving line-of-sight distance in Mpc at a given redshift.

The comoving distance along the line-of-sight between two objects remains constant with time for objects in the Hubble flow.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
dQuantity [:ref: ‘length’]

Comoving distance in Mpc to each input redshift.

comoving_transverse_distance(z)#

Comoving transverse distance in Mpc at a given redshift.

This value is the transverse comoving distance at redshift z corresponding to an angular separation of 1 radian. This is the same as the comoving distance if \(\Omega_k\) is zero (as in the current concordance Lambda-CDM model).

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
dQuantity [:ref: ‘length’]

Comoving transverse distance in Mpc at each input redshift.

Notes

This quantity is also called the ‘proper motion distance’ in some texts.

comoving_volume(z)#

Comoving volume in cubic Mpc at redshift z.

This is the volume of the universe encompassed by redshifts less than z. For the case of \(\Omega_k = 0\) it is a sphere of radius comoving_distance but it is less intuitive if \(\Omega_k\) is not.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
VQuantity

Comoving volume in \(Mpc^3\) at each input redshift.

critical_density(z)#

Critical density in grams per cubic cm at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
rhoQuantity

Critical density in g/cm^3 at each input redshift.

de_density_scale(z)#

Evaluates the redshift dependence of the dark energy density.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Indarray or float

The scaling of the energy density of dark energy with redshift. Returns float if the input is scalar.

Notes

The scaling factor, I, is defined by \(\rho(z) = \rho_0 I\), and in this case is given by \(I = 1\).

differential_comoving_volume(z)#

Differential comoving volume at redshift z.

Useful for calculating the effective comoving volume. For example, allows for integration over a comoving volume that has a sensitivity function that changes with redshift. The total comoving volume is given by integrating differential_comoving_volume to redshift z and multiplying by a solid angle.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
dVQuantity

Differential comoving volume per redshift per steradian at each input redshift.

distmod(z)#

Distance modulus at redshift z.

The distance modulus is defined as the (apparent magnitude - absolute magnitude) for an object at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
distmodQuantity [:ref: ‘length’]

Distance modulus at each input redshift, in magnitudes.

See also

z_at_value

Find the redshift corresponding to a distance modulus.

efunc(z)[source]#

Function used to calculate H(z), the Hubble parameter.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Endarray or float

The redshift scaling of the Hubble constant. Returns float if the input is scalar. Defined such that \(H(z) = H_0 E(z)\).

inv_efunc(z)[source]#

Function used to calculate \(\frac{1}{H_z}\).

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Endarray or float

The inverse redshift scaling of the Hubble constant. Returns float if the input is scalar. Defined such that \(H_z = H_0 / E\).

is_equivalent(other: Any, /, *, format: _FormatType = False) bool#

Check equivalence between Cosmologies.

Two cosmologies may be equivalent even if not the same class. For example, an instance of LambdaCDM might have \(\Omega_0=1\) and \(\Omega_k=0\) and therefore be flat, like FlatLambdaCDM.

Parameters:
otherCosmology subclass instance, positional-only

The object to which to compare.

formatbool or None or str, optional keyword-only

Whether to allow, before equivalence is checked, the object to be converted to a Cosmology. This allows, e.g. a Table to be equivalent to a Cosmology. False (default) will not allow conversion. True or None will, and will use the auto-identification to try to infer the correct format. A str is assumed to be the correct format to use when converting. format is broadcast to match the shape of other. Note that the cosmology arguments are not broadcast against format, so it cannot determine the output shape.

Returns:
bool

True if cosmologies are equivalent, False otherwise.

Examples

Two cosmologies may be equivalent even if not of the same class. In this examples the LambdaCDM has Ode0 set to the same value calculated in FlatLambdaCDM.

>>> import astropy.units as u
>>> from astropy.cosmology import LambdaCDM, FlatLambdaCDM
>>> cosmo1 = LambdaCDM(70 * (u.km/u.s/u.Mpc), 0.3, 0.7)
>>> cosmo2 = FlatLambdaCDM(70 * (u.km/u.s/u.Mpc), 0.3)
>>> cosmo1.is_equivalent(cosmo2)
True

While in this example, the cosmologies are not equivalent.

>>> cosmo3 = FlatLambdaCDM(70 * (u.km/u.s/u.Mpc), 0.3, Tcmb0=3 * u.K)
>>> cosmo3.is_equivalent(cosmo2)
False

Also, using the keyword argument, the notion of equivalence is extended to any Python object that can be converted to a Cosmology.

>>> from astropy.cosmology import Planck18
>>> tbl = Planck18.to_format("astropy.table")
>>> Planck18.is_equivalent(tbl, format=True)
True

The list of valid formats, e.g. the Table in this example, may be checked with Cosmology.from_format.list_formats().

As can be seen in the list of formats, not all formats can be auto-identified by Cosmology.from_format.registry. Objects of these kinds can still be checked for equivalence, but the correct format string must be used.

>>> tbl = Planck18.to_format("yaml")
>>> Planck18.is_equivalent(tbl, format="yaml")
True
kpc_comoving_per_arcmin(z)#

Separation in transverse comoving kpc equal to an arcmin at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
dQuantity [:ref: ‘length’]

The distance in comoving kpc corresponding to an arcmin at each input redshift.

kpc_proper_per_arcmin(z)#

Separation in transverse proper kpc equal to an arcminute at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
dQuantity [:ref: ‘length’]

The distance in proper kpc corresponding to an arcmin at each input redshift.

lookback_distance(z)#

The lookback distance is the light travel time distance to a given redshift.

It is simply c * lookback_time. It may be used to calculate the proper distance between two redshifts, e.g. for the mean free path to ionizing radiation.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
dQuantity [:ref: ‘length’]

Lookback distance in Mpc

lookback_time(z)#

Lookback time in Gyr to redshift z.

The lookback time is the difference between the age of the Universe now and the age at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
tQuantity [:ref: ‘time’]

Lookback time in Gyr to each input redshift.

See also

z_at_value

Find the redshift corresponding to a lookback time.

lookback_time_integrand(z)#

Integrand of the lookback time (equation 30 of [1]).

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
Ifloat or array

The integrand for the lookback time.

References

[1] (1,2)

Hogg, D. (1999). Distance measures in cosmology, section 11. arXiv e-prints, astro-ph/9905116.

luminosity_distance(z)#

Luminosity distance in Mpc at redshift z.

This is the distance to use when converting between the bolometric flux from an object at redshift z and its bolometric luminosity [1].

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
dQuantity [:ref: ‘length’]

Luminosity distance in Mpc at each input redshift.

See also

z_at_value

Find the redshift corresponding to a luminosity distance.

References

[1]

Weinberg, 1972, pp 420-424; Weedman, 1986, pp 60-62.

nu_relative_density(z)#

Neutrino density function relative to the energy density in photons.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
fndarray or float

The neutrino density scaling factor relative to the density in photons at each redshift. Only returns float if z is scalar.

Notes

The density in neutrinos is given by

\[\rho_{\nu} \left(a\right) = 0.2271 \, N_{eff} \, f\left(m_{\nu} a / T_{\nu 0} \right) \, \rho_{\gamma} \left( a \right)\]

where

\[f \left(y\right) = \frac{120}{7 \pi^4} \int_0^{\infty} \, dx \frac{x^2 \sqrt{x^2 + y^2}} {e^x + 1}\]

assuming that all neutrino species have the same mass. If they have different masses, a similar term is calculated for each one. Note that f has the asymptotic behavior \(f(0) = 1\). This method returns \(0.2271 f\) using an analytical fitting formula given in Komatsu et al. 2011, ApJS 192, 18.

scale_factor(z)#

Scale factor at redshift z.

The scale factor is defined as \(a = 1 / (1 + z)\).

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
andarray or float

Scale factor at each input redshift. Returns float if the input is scalar.

w(z)#

Returns dark energy equation of state at redshift z.

Parameters:
zQuantity-like [‘redshift’], array_like, or Number

Input redshift.

Returns:
wndarray or float

The dark energy equation of state. Returns float if the input is scalar.

Notes

The dark energy equation of state is defined as \(w(z) = P(z)/\rho(z)\), where \(P(z)\) is the pressure at redshift z and \(\rho(z)\) is the density at redshift z, both in units where c=1. Here this is \(w(z) = -1\).