z_at_value¶
- astropy.cosmology.z_at_value(func, fval, zmin=1e-08, zmax=1000, ztol=1e-08, maxfun=500, method='Brent', bracket=None, verbose=False)[source]¶
Find the redshift
z
at whichfunc(z) = fval
.This finds the redshift at which one of the cosmology functions or methods (for example Planck13.distmod) is equal to a known value.
Warning
Make sure you understand the behavior of the function that you are trying to invert! Depending on the cosmology, there may not be a unique solution. For example, in the standard Lambda CDM cosmology, there are two redshifts which give an angular diameter distance of 1500 Mpc, z ~ 0.7 and z ~ 3.8. To force
z_at_value
to find the solution you are interested in, use thezmin
andzmax
keywords to limit the search range (see the example below).- Parameters:
- funcpython:function or method
A function that takes a redshift as input.
- fval
Quantity
The (scalar or array) value of
func(z)
to recover.- zmin
python:float
or numpy:array_like[:ref: ‘dimensionless’] or astropy:quantity-like, optional The lower search limit for
z
. Beware of divergences in some cosmological functions, such as distance moduli, at z=0 (default 1e-8).- zmax
python:float
or numpy:array_like[:ref: ‘dimensionless’] or astropy:quantity-like, optional The upper search limit for
z
(default 1000).- ztol
python:float
or numpy:array_like[:ref: ‘dimensionless’], optional The relative error in
z
acceptable for convergence.- maxfun
python:int
or numpy:array_like, optional The maximum number of function evaluations allowed in the optimization routine (default 500).
- method
python:str
orpython:callable()
, optional Type of solver to pass to the minimizer. The built-in options provided by
minimize_scalar()
are ‘Brent’ (default), ‘Golden’ and ‘Bounded’ with names case insensitive - see documentation there for details. It also accepts a custom solver by passing any user-provided callable object that meets the requirements listed therein under the Notes on “Custom minimizers” - or in more detail in scipy:tutorial/optimize - although their use is currently untested.New in version 4.3.
- bracketpython:sequence or
object
array
[python:sequence], optional For methods ‘Brent’ and ‘Golden’,
bracket
defines the bracketing interval and can either have three items (z1, z2, z3) so that z1 < z2 < z3 andfunc(z2) < func (z1), func(z3)
or two items z1 and z3 which are assumed to be a starting interval for a downhill bracket search. For non-monotonic functions such as angular diameter distance this may be used to start the search on the desired side of the maximum, but see Examples below for usage notes.New in version 4.3.
- verbosebool, optional
Print diagnostic output from solver (default
False
).New in version 4.3.
- Returns:
- z
Quantity
[‘redshift’] The redshift
z
satisfyingzmin < z < zmax
andfunc(z) = fval
withinztol
. Has units of cosmological redshift.
- z
- Raises:
astropy.cosmology.CosmologyError
If the result is very close to either
zmin
orzmax
.ValueError
If
bracket
is not an array nor a 2 (or 3) element sequence.TypeError
If
bracket
is not an object array. 2 (or 3) element sequences will be turned into object arrays, so this error should only occur if a non-object array is used forbracket
.
- Warns:
AstropyUserWarning
If
fval
is not bracketed byfunc(zmin)=fval(zmin)
andfunc(zmax)=fval(zmax)
.If the solver was not successful.
Notes
This works for any arbitrary input cosmology, but is inefficient if you want to invert a large number of values for the same cosmology. In this case, it is faster to instead generate an array of values at many closely-spaced redshifts that cover the relevant redshift range, and then use interpolation to find the redshift at each value you are interested in. For example, to efficiently find the redshifts corresponding to 10^6 values of the distance modulus in a Planck13 cosmology, you could do the following:
>>> import astropy.units as u >>> from astropy.cosmology import Planck13, z_at_value
Generate 10^6 distance moduli between 24 and 44 for which we want to find the corresponding redshifts:
>>> Dvals = (24 + np.random.rand(1000000) * 20) * u.mag
Make a grid of distance moduli covering the redshift range we need using 50 equally log-spaced values between zmin and zmax. We use log spacing to adequately sample the steep part of the curve at low distance moduli:
>>> zmin = z_at_value(Planck13.distmod, Dvals.min()) >>> zmax = z_at_value(Planck13.distmod, Dvals.max()) >>> zgrid = np.geomspace(zmin, zmax, 50) >>> Dgrid = Planck13.distmod(zgrid)
Finally interpolate to find the redshift at each distance modulus:
>>> zvals = np.interp(Dvals.value, Dgrid.value, zgrid)
Examples
>>> import astropy.units as u >>> from astropy.cosmology import Planck13, Planck18, z_at_value
The age and lookback time are monotonic with redshift, and so a unique solution can be found:
>>> z_at_value(Planck13.age, 2 * u.Gyr) <Quantity 3.19812268 redshift>
The angular diameter is not monotonic however, and there are two redshifts that give a value of 1500 Mpc. You can use the zmin and zmax keywords to find the one you are interested in:
>>> z_at_value(Planck18.angular_diameter_distance, ... 1500 * u.Mpc, zmax=1.5) <Quantity 0.68044452 redshift> >>> z_at_value(Planck18.angular_diameter_distance, ... 1500 * u.Mpc, zmin=2.5) <Quantity 3.7823268 redshift>
Alternatively the
bracket
option may be used to initialize the function solver on a desired region, but one should be aware that this does not guarantee it will remain close to this starting bracket. For the example of angular diameter distance, which has a maximum near a redshift of 1.6 in this cosmology, defining a bracket on either side of this maximum will often return a solution on the same side:>>> z_at_value(Planck18.angular_diameter_distance, 1500 * u.Mpc, ... method="Brent", bracket=(1.0, 1.2)) <Quantity 0.68044452 redshift>
But this is not ascertained especially if the bracket is chosen too wide and/or too close to the turning point:
>>> z_at_value(Planck18.angular_diameter_distance, ... 1500 * u.Mpc, bracket=(0.1, 1.5)) <Quantity 3.7823268 redshift>
Likewise, even for the same minimizer and same starting conditions different results can be found depending on architecture or library versions:
>>> z_at_value(Planck18.angular_diameter_distance, ... 1500 * u.Mpc, bracket=(2.0, 2.5)) <Quantity 3.7823268 redshift>
>>> z_at_value(Planck18.angular_diameter_distance, ... 1500 * u.Mpc, bracket=(2.0, 2.5)) <Quantity 0.68044452 redshift>
It is therefore generally safer to use the 3-parameter variant to ensure the solution stays within the bracketing limits:
>>> z_at_value(Planck18.angular_diameter_distance, 1500 * u.Mpc, method="Brent", ... bracket=(0.1, 1.0, 1.5)) <Quantity 0.68044452 redshift>
Also note that the luminosity distance and distance modulus (two other commonly inverted quantities) are monotonic in flat and open universes, but not in closed universes.
All the arguments except
func
,method
andverbose
accept array inputs. This does NOT use interpolation tables or any method to speed up evaluations, rather providing a convenient means to broadcast arguments over an element-wise scalar evaluation.The most common use case for non-scalar input is to evaluate ‘func’ for an array of
fval
:>>> z_at_value(Planck13.age, [2, 7] * u.Gyr) <Quantity [3.19812061, 0.75620443] redshift>
fval
can be any shape:>>> z_at_value(Planck13.age, [[2, 7], [1, 3]]*u.Gyr) <Quantity [[3.19812061, 0.75620443], [5.67661227, 2.19131955]] redshift>
Other arguments can be arrays. For non-monotic functions – for example, the angular diameter distance – this can be useful to find all solutions.
>>> z_at_value(Planck13.angular_diameter_distance, 1500 * u.Mpc, ... zmin=[0, 2.5], zmax=[2, 4]) <Quantity [0.68127747, 3.79149062] redshift>
The
bracket
argument can likewise be be an array. However, since bracket must already be a sequence (or None), it MUST be given as an objectnumpy.ndarray
. Importantly, the depth of the array must be such that each bracket subsequence is an object. Errors or unexpected results will happen otherwise. A convenient means to ensure the right depth is by including a length-0 tuple as a bracket and then truncating the object array to remove the placeholder. This can be seen in the following example:>>> bracket=np.array([(1.0, 1.2),(2.0, 2.5), ()], dtype=object)[:-1] >>> z_at_value(Planck18.angular_diameter_distance, 1500 * u.Mpc, ... bracket=bracket) <Quantity [0.68044452, 3.7823268] redshift>