Shells¶
The following functions provide functionality for the definition of matter shells, i.e. the radial discretisation of the light cone.
Window functions¶
- class glass.RadialWindow(za, wa, zeff=0)[source]¶
A radial window, defined by a window function.
The radial window is defined by a window function in redshift, which is given by a pair of arrays
za
,wa
.The radial window also has an effective redshift, stored in the
zeff
attribute, which should be a representative redshift for the window function.To prevent accidental inconsistencies, instances of this type are immutable (however, the array entries may not be immutable; do not change them in place):
>>> from glass import RadialWindow >>> w1 = RadialWindow(..., ..., zeff=0.1) >>> w1.zeff = 0.15 Traceback (most recent call last): File "<stdin>", line 1, in <module> AttributeError: can't set attribute
To create a new instance with a changed attribute value, use the
._replace
method:>>> w1 = w1._replace(zeff=0.15) >>> w1 RadialWindow(za=..., wa=..., zeff=0.15)
- za¶
Redshift array; the abscissae of the window function.
- wa¶
Weight array; the values (ordinates) of the window function.
- zeff¶
Effective redshift of the window.
- _replace()¶
Create a new instance with changed attribute values.
- glass.tophat_windows(zbins, dz=0.001, weight=None)[source]¶
Tophat window functions from the given redshift bin edges.
Uses the N+1 given redshifts as bin edges to construct N tophat window functions. The redshifts of the windows have linear spacing approximately equal to
dz
.An optional weight function \(w(z)\) can be given using
weight
; it is applied to the tophat windows.The resulting windows functions are
RadialWindow
instances. Their effective redshifts are the mean redshifts of the (weighted) tophat bins.- Parameters:
zbins (
Sequence
[float
] |ndarray
[Any
,dtype
[float64
]]) – Redshift bin edges for the tophat window functions.dz (
float
) – Approximate spacing of the redshift grid.weight (
Callable
[[Sequence
[float
] |ndarray
[Any
,dtype
[float64
]]],ndarray
[Any
,dtype
[float64
]]] |None
) – If given, a weight function to be applied to the window functions.
- Returns:
A list of window functions.
- Raises:
ValueError – If the number of redshift bins is less than 2.
- Return type:
See also
- glass.linear_windows(zgrid, dz=0.001, weight=None)[source]¶
Linear interpolation window functions.
Uses the N+2 given redshifts as nodes to construct N triangular window functions between the first and last node. These correspond to linear interpolation of radial functions. The redshift spacing of the windows is approximately equal to
dz
.An optional weight function \(w(z)\) can be given using
weight
; it is applied to the triangular windows.The resulting windows functions are
RadialWindow
instances. Their effective redshifts correspond to the given nodes.- Parameters:
zgrid (
Sequence
[float
] |ndarray
[Any
,dtype
[float64
]]) – Redshift grid for the triangular window functions.dz (
float
) – Approximate spacing of the redshift grid.weight (
Callable
[[Sequence
[float
] |ndarray
[Any
,dtype
[float64
]]],ndarray
[Any
,dtype
[float64
]]] |None
) – If given, a weight function to be applied to the window functions.
- Returns:
A list of window functions.
- Raises:
ValueError – If the number of nodes is less than 3.
- Return type:
See also
- glass.cubic_windows(zgrid, dz=0.001, weight=None)[source]¶
Cubic interpolation window functions.
Uses the N+2 given redshifts as nodes to construct N cubic Hermite spline window functions between the first and last node. These correspond to cubic spline interpolation of radial functions. The redshift spacing of the windows is approximately equal to
dz
.An optional weight function \(w(z)\) can be given using
weight
; it is applied to the cubic spline windows.The resulting windows functions are
RadialWindow
instances. Their effective redshifts correspond to the given nodes.- Parameters:
zgrid (
Sequence
[float
] |ndarray
[Any
,dtype
[float64
]]) – Redshift grid for the cubic spline window functions.dz (
float
) – Approximate spacing of the redshift grid.weight (
Callable
[[Sequence
[float
] |ndarray
[Any
,dtype
[float64
]]],ndarray
[Any
,dtype
[float64
]]] |None
) – If given, a weight function to be applied to the window functions.
- Returns:
A list of window functions.
- Raises:
ValueError – If the number of nodes is less than 3.
- Return type:
See also
Window function tools¶
- glass.restrict(z, f, w)[source]¶
Restrict a function to a redshift window.
Multiply the function \(f(z)\) by a window function \(w(z)\) to produce \(w(z) f(z)\) over the support of \(w\).
The function \(f(z)\) is given by redshifts
z
of shape (N,) and function valuesf
of shape (…, N), with any number of leading axes allowed.The window function \(w(z)\) is given by
w
, which must be aRadialWindow
instance or compatible with it.The restriction has redshifts that are the union of the redshifts of the function and window over the support of the window. Intermediate function values are found by linear interpolation
- glass.partition(z, fz, shells, *, method='nnls')[source]¶
Partition a function by a sequence of windows.
Returns a vector of weights \(x_1, x_2, \ldots\) such that the weighted sum of normalised radial window functions \(x_1 \, w_1(z) + x_2 \, w_2(z) + \ldots\) approximates the given function \(f(z)\).
The function \(f(z)\) is given by redshifts z of shape (N,) and function values fz of shape (…, N), with any number of leading axes allowed.
The window functions are given by the sequence shells of
RadialWindow
or compatible entries.- Parameters:
z (
ndarray
[Any
,dtype
[float64
]]) – The function to be partitioned. If f is multi-dimensional, its last axis must agree with z.fz (
ndarray
[Any
,dtype
[float64
]]) – The function to be partitioned. If f is multi-dimensional, its last axis must agree with z.shells (
Sequence
[RadialWindow
]) – Ordered sequence of window functions for the partition.method (
str
) – Method for the partition. See notes for description. The options are “lstsq”, “nnls”, “restrict”.
- Returns:
The weights of the partition, where the leading axis corresponds to shells.
- Raises:
ValueError – If the method is not recognised.
- Return type:
Notes
Formally, if \(w_i\) are the normalised window functions, \(f\) is the target function, and \(z_i\) is a redshift grid with intervals \(\Delta z_i\), the partition problem seeks an approximate solution of
\[\begin{pmatrix} w_1(z_1) \Delta z_1 & w_2(z_1) \, \Delta z_1 & \cdots \\ w_1(z_2) \Delta z_2 & w_2(z_2) \, \Delta z_2 & \cdots \\ \vdots & \vdots & \ddots \end{pmatrix} \, \begin{pmatrix} x_1 \\ x_2 \\ \vdots \end{pmatrix} = \begin{pmatrix} f(z_1) \, \Delta z_1 \\ f(z_2) \, \Delta z_2 \\ \vdots \end{pmatrix} \;.\]The redshift grid is the union of the given array z and the redshift arrays of all window functions. Intermediate function values are found by linear interpolation.
When partitioning a density function, it is usually desirable to keep the normalisation fixed. In that case, the problem can be enhanced with the further constraint that the sum of the solution equals the integral of the target function,
\[\begin{pmatrix} w_1(z_1) \Delta z_1 & w_2(z_1) \, \Delta z_1 & \cdots \\ w_1(z_2) \Delta z_2 & w_2(z_2) \, \Delta z_2 & \cdots \\ \vdots & \vdots & \ddots \\ \hline \lambda & \lambda & \cdots \end{pmatrix} \, \begin{pmatrix} x_1 \\ x_2 \\ \vdots \end{pmatrix} = \begin{pmatrix} f(z_1) \, \Delta z_1 \\ f(z_2) \, \Delta z_2 \\ \vdots \\ \hline \lambda \int \! f(z) \, dz \end{pmatrix} \;,\]where \(\lambda\) is a multiplier to enforce the integral constraints.
The
partition()
function implements a number of methods to obtain a solution:If
method="nnls"
(the default), obtain a partition from a non-negative least-squares solution. This will usually match the shape of the input function closely. The contribution from each shell is a positive number, which is required to partition e.g. density functions.If
method="lstsq"
, obtain a partition from an unconstrained least-squares solution. This will more closely match the shape of the input function, but might lead to shells with negative contributions.If
method="restrict"
, obtain a partition by integrating the restriction (usingrestrict()
) of the function \(f\) to each window. For overlapping shells, this method might produce results which are far from the input function.
- glass.combine(z, weights, shells)[source]¶
Evaluate a linear combination of window functions.
Takes a vector of weights \(x_1, x_2, \ldots\) and computes the weighted sum of normalised radial window functions \(f(z) = x_1 \, w_1(z) + x_2 \, w_2(z) + \ldots\) in the given redshifts \(z\).
The window functions are given by the sequence shells of
RadialWindow
or compatible entries.- Parameters:
- Returns:
A linear combination of window functions, evaluated in z.
- Return type:
See also
partition
Find weights for a given function.
Redshift grids¶
- glass.redshift_grid(zmin, zmax, *, dz=None, num=None)[source]¶
Redshift grid with uniform spacing in redshift.
- glass.distance_grid(cosmo, zmin, zmax, *, dx=None, num=None)[source]¶
Redshift grid with uniform spacing in comoving distance.
- Parameters:
- Returns:
The redshift grid.
- Raises:
ValueError – If both
dx
andnum
are given.- Return type: