Home | Trees | Indices | Help |
---|
|
Cosmological distance measures.
Mostly follows David Hogg's pedagogical paper arXiv:astro-ph/9905116v4 .
Distance units are Mpc, time units are seconds.
Functions | |||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
Variables | |
__package__ =
|
Function Details |
'Spatial curvature density' omega_k_0 for a cosmology (if needed). If omega_k_0 is specified, return it. Otherwise return: 1.0 - omega_M_0 - omega_lambda_0 |
Returns the cosmo dictionary with omega_k_0 set. See get_omega_k_0. Note that cosmo is not passed as **cosmo for once. This function modifies the dictionary in place and returns the result. |
The unitless Hubble expansion rate at redshift z. In David Hogg's (arXiv:astro-ph/9905116v4) formalism, this is equivalent to E(z), defined in his eq. 14. Modified (JBJ, 29-Feb-2012) to include scalar w parameter |
The value of the Hubble constant at redshift z. Units are s^-1 In David Hogg's (arXiv:astro-ph/9905116v4) formalism, this is equivalent to H_0 * E(z) (see his eq. 14). |
The value of the Hubble distance at redshift z. Units are Mpc. In David Hogg's (arXiv:astro-ph/9905116v4) formalism, this is equivalent to D_H / E(z) = c / (H_0 E(z)) [see his eq. 14], which appears in the definitions of many other distance measures. |
The derivative of the comoving distance with redshift: dd_c/dz. See equation 15 of David Hogg's arXiv:astro-ph/9905116v4 Units are Mpc. |
The line-of-sight comoving distance (in Mpc) to redshift z. See equation 15 of David Hogg's arXiv:astro-ph/9905116v4 Units are Mpc. Optionally calculate the integral from z0 to z. Returns
Examples>>> import cosmolopy.distance as cd >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> d_co = cd.comoving_distance(6., **cosmo) >>> print "Comoving distance to z=6 is %.1f Mpc" % (d_co) Comoving distance to z=6 is 8017.8 Mpc |
Returns comoving_distance_transverse. Units are Mpc. Examples>>> import cosmolopy.distance as cd >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> d_M = cd.proper_motion_distance(6., **cosmo) >>> print "Transverse comoving distance to z=6 is %.1f Mpc" % (d_M) Transverse comoving distance to z=6 is 8017.8 Mpc |
The transverse comoving distance (in Mpc) to redshift z. This is also called the proper motion distance, D_M. See equation 16 of David Hogg's arXiv:astro-ph/9905116v4 Units are Mpc. This is the distance d_m, such that the comoving distance between two events at the same redshift, but separated on the sky by some angle delta_theta is d_m * delta_theta. Examples>>> import cosmolopy.distance as cd >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> d_M = cd.comoving_distance_transverse(6., **cosmo) >>> print "Transverse comoving distance to z=6 is %.1f Mpc" % (d_M) Transverse comoving distance to z=6 is 8017.8 Mpc |
The angular-diameter distance (Mpc) to redshift z. Optionally find the angular diameter distance between objects at z0 and z (only implemented for omega_k_0 >= 0). See equations 18-19 of David Hogg's arXiv:astro-ph/9905116v4 Units are Mpc. Examples>>> import cosmolopy.distance as cd >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> d_a = cd.angular_diameter_distance(6., **cosmo) >>> print "Angular diameter distance = %.1f Mpc" % (d_a) Angular diameter distance = 1145.4 Mpc |
The luminosity distance to redshift z. Units are Mpc. See, for example, David Hogg's arXiv:astro-ph/9905116v4 |
The differential comoving volume element dV_c/dz/dSolidAngle. Dimensions are volume per unit redshift per unit solid angle. Units are Mpc**3 Steradians^-1. See David Hogg's arXiv:astro-ph/9905116v4, equation 28. Examples>>> import cosmolopy.distance as cd >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> dVc = cd.diff_comoving_volume(6.0, **cosmo) >>> print "dV/dz/dSolidAngle at z=6 is %.3g Mpc**3" % (dVc) dV/dz/dSolidAngle at z=6 is 2.63e+10 Mpc**3 |
The comoving volume out to redshift z. See David Hogg's arXiv:astro-ph/9905116v4, equation 29. Examples>>> import cosmolopy.distance as cd >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> Vc = cd.comoving_volume(6.0, **cosmo) >>> print "Vc = %.3g Mpc**3" % (Vc) Vc = 2.16e+12 Mpc**3 >>> import cosmolopy.distance as cd >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.0, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> Vc = cd.comoving_volume(6.0, **cosmo) >>> print "Vc = %.3g Mpc**3" % (Vc) Vc = 1.68e+12 Mpc**3 |
The derivative of the lookback time with redshift: dt_L/dz. See equation 30 of David Hogg's arXiv:astro-ph/9905116v4 Units are seconds. |
The lookback time (in s) to redshift z. See equation 30 of David Hogg's arXiv:astro-ph/9905116v4 Units are s. Optionally calculate the integral from z0 to z. Returns
|
The age of the universe as seen at redshift z. Age at z is lookback time at z'->Infinity minus lookback time at z. See also: lookback_time. Units are s. Examples>>> import cosmolopy.distance as cd >>> import cosmolopy.constants as cc >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> t = cd.age(6.0, **cosmo) >>> print "age at z=6.0 = %.3g Gyr" % (t/cc.Gyr_s) age at z=6.0 = 0.892 Gyr |
The age of the universe assuming a flat cosmology. Units are s. Analytical formula from Peebles, p. 317, eq. 13.2. Examples>>> import cosmolopy.distance as cd >>> import cosmolopy.constants as cc >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> t = cd.age_flat(6.0, **cosmo) >>> print "age at z=6.0 is %.3g Gyr" % (t/cc.Gyr_s) age at z=6.0 is 0.892 Gyr |
Return an interpolation function that will give distance as a funtion of z If return_inverse is True, will also return a function giving z as a function of distance. Inputsfunction -- the distance function to interpolate (can be any callable that takes a redshift argument plus cosmology keywords). k -- spline order ( Returnsdistfunc or distfunc, zfunc Examples>>> import cosmolopy.distance as cd >>> import cosmolopy.constants as cc >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> distfunc, redfunc = cd.quick_distance_function(cd.luminosity_distance, return_inverse=True, **cosmo) >>> d = distfunc(6.3333) >>> z = redfunc(d) >>> "%.1g" % (distfunc(6.3333)/cd.luminosity_distance(6.3333, **cosmo) - 1.0) '-2e-16' >>> "%.1g" % (z/6.3333 - 1.0) '0' |
Return an interpolation function that will give age as a funtion of z Units are s. If return_inverse is True, will also return a function giving z as a function of age. Returnsagefunc or agefunc, redfunc Examples>>> import cosmolopy.distance as cd >>> import cosmolopy.constants as cc >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> agefunc = cd.quick_age_function(**cosmo) >>> t = agefunc(6.0) >>> print "age at z=6.0 is %.3g Gyr" % (t/cc.Gyr_s) age at z=6.0 is 0.892 Gyr |
Return an interpolation function giving z as a funtion of age of the universe. Units of time are s. Returnsredfunc Examples>>> import cosmolopy.distance as cd >>> import cosmolopy.constants as cc >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> redfunc = cd.quick_redshift_age_function(**cosmo) >>> z = redfunc(1.0 * cc.Gyr_s) >>> print "When age=1.0Gyr z=%.2f" % (z) When age=1.0Gyr z=5.49 |
The light travel distance to redshift z. Units are Mpc. Examples>>> import cosmolopy.distance as cd >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> dlookback = cd.light_travel_distance(3.0, 2.0, **cosmo) >>> print "Lookback distance from z=2 to 3 is %.2g Mpc" % (dlookback) Lookback distance from z=2 to 3 is 3.3e+02 Mpc |
The redshift corresponding to a given light travel distance. Units are the same as light_travel_distance (Mpc). Examples>>> import cosmolopy.distance as cd >>> import cosmolopy.constants as cc >>> cosmo = {'omega_M_0' : 0.3, 'omega_lambda_0' : 0.7, 'h' : 0.72} >>> cosmo = cd.set_omega_k_0(cosmo) >>> z = cd.redshift_d_light(10. * cc.c_light_Mpc_Gyr, **cosmo) Optimization terminated successfully. Current function value: 0.000112 Iterations: 26 Function evaluations: 52 >>> print "Redshift for a lookback time of 10Gyr is z=%.3f" % (z) Redshift for a lookback time of 10Gyr is z=2.025 |
Home | Trees | Indices | Help |
---|
Generated by Epydoc 3.0.1 on Tue Mar 27 20:49:14 2012 | http://epydoc.sourceforge.net |