External API
This section documents the public functions intended for users.
Effort.q_par_perp — Methodq_par_perp(z, cosmo_mcmc::AbstractCosmology, cosmo_ref::AbstractCosmology)Calculates the parallel (q_par) and perpendicular (q_perp) Alcock-Paczynski (AP) parameters at a given redshift z, comparing a varying cosmology to a reference cosmology.
The AP parameters quantify the distortion of observed clustering due to assuming a different cosmology than the true one when converting redshifts and angles to distances.
Arguments
z: The redshift at which to calculate the AP parameters.cosmo_mcmc: AnAbstractCosmologystruct representing the varying cosmology (e.g., from an MCMC chain).cosmo_ref: AnAbstractCosmologystruct representing the reference cosmology used for measurements.
Returns
A tuple (q_par, q_perp) containing the calculated parallel and perpendicular AP parameters at redshift z.
Details
The parallel AP parameter q_par is the ratio of the Hubble parameter in the reference cosmology to that in the varying cosmology. The perpendicular AP parameter q_perp is the ratio of the conformal angular diameter distance in the varying cosmology to that in the reference cosmology.
Formula
The formulas for the Alcock-Paczynski parameters are:
\[q_\parallel(z) = \frac{E_{\text{ref}}(z)}{E_{\text{mcmc}}(z)}\]
\[q_\perp(z) = \frac{\tilde{d}_{A,\text{mcmc}}(z)}{\tilde{d}_{A,\text{ref}}(z)}\]
where $E(z)$ is the normalized Hubble parameter and $\tilde{d}_A(z)$ is the conformal angular diameter distance.
Effort.apply_AP — Functionapply_AP(k_input::AbstractVector, k_output::AbstractVector, mono::AbstractVector, quad::AbstractVector, hexa::AbstractVector, q_par, q_perp; n_GL_points=8)Calculates the observed power spectrum multipole moments (monopole, quadrupole, hexadecapole) on a given observed wavenumber grid k_output, using arrays of true multipole moments provided on an input wavenumber grid k_input, and employing Gauss-Lobatto quadrature.
This is the standard, faster implementation for applying the Alcock-Paczynski (AP) effect to the power spectrum multipoles, designed for performance compared to the check version using generic numerical integration.
Arguments
k_input: A vector of wavenumber values on which the input true multipole moments (mono,quad,hexa) are defined.k_output: A vector of observed wavenumber values at which to calculate the output observed multipoles.mono: A vector containing the values of the true monopole moment $I_0(k)$ on thek_inputgrid.quad: A vector containing the values of the true quadrupole moment $I_2(k)$ on thek_inputgrid.hexa: A vector containing the values of the true hexadecapole moment $I_4(k)$ on thek_inputgrid.q_par: A parameter related to parallel anisotropic scaling.q_perp: A parameter related to perpendicular anisotropic scaling.
Keyword Arguments
n_GL_points: The number of Gauss-Lobatto points to use for the integration overμ. The actual number of nodes used corresponds to2 * n_GL_points. Defaults to 8.
Returns
A tuple (P0_obs, P2_obs, P4_obs), where each element is a vector containing the calculated observed monopole, quadrupole, and hexadecapole moments respectively, evaluated at the observed wavenumbers in k_output.
Details
The function applies the AP and RSD effects by integrating the observed anisotropic power spectrum $P_{\text{obs}}(k_o, \mu_o)$ over the observed cosine of the angle to the line-of-sight $\mu_o \in [0, 1]$ (assuming symmetry for even multipoles), weighted by the corresponding Legendre polynomial $\mathcal{L}_\ell(\mu_o)$.
The process involves:
- Determine Gauss-Lobatto nodes and weights for the interval
[0, 1]. - For each observed wavenumber
k_oin the inputk_outputarray and eachμ_onode: a. Calculate the true wavenumber $k_t(k_o, \mu_o)$ using_k_true. b. Calculate the true angle cosine $\mu_t(\mu_o)$ using_μ_true. c. Interpolate the true multipole moments $I_\ell(k_t)$ using_akima_interpolation, interpolating from thek_inputgrid to the newk_tvalues. d. Calculate the true Legendre polynomials $\mathcal{L}_\ell(\mu_t)$ using_Legendre_0,_Legendre_2,_Legendre_4. e. Reconstruct the true power spectrum $P(k_t, \mu_t)$ using_Pk_recon. f. Calculate the observed power spectrum $P_{\text{obs}}(k_o, \mu_o) = P(k_t, \mu_t) / (q_\parallel q_\perp^2)$. - Perform the weighted sum (quadrature) over the
μ_onodes to get the observed multipoles $P_\ell(k_o)$ on thek_outputgrid.
This function is the standard, performant implementation for applying AP compared to the slower apply_AP_check.
Formula
The observed multipole moments are calculated using the formula:
\[P_\ell(k_o) = (2\ell + 1) \int_{0}^1 P_{\text{obs}}(k_o, \mu_o) \mathcal{L}_\ell(\mu_o) d\mu_o\]
for $\ell \in \{0, 2, 4\}$. The integral is approximated using Gauss-Lobatto quadrature.
See Also
apply_AP_check: The slower, check version using generic numerical integration._k_true: Transforms observed wavenumber to true wavenumber._μ_true: Transforms observed angle cosine to true angle cosine._Legendre_0,_Legendre_2,_Legendre_4: Calculate the Legendre polynomials._akima_interpolation: Interpolates the true multipole moments._Pk_recon: Reconstructs the true power spectrum on a grid.gausslobatto: Function used to get quadrature nodes and weights.
apply_AP(k_input::AbstractVector, k_output::AbstractVector, mono::AbstractMatrix, quad::AbstractMatrix, hexa::AbstractMatrix, q_par, q_perp; n_GL_points=8)Batch version of apply_AP for processing multiple columns simultaneously using optimized matrix Akima interpolation.
This method applies the Alcock-Paczynski effect to multiple sets of multipole moments (e.g., multiple Jacobian columns or parameter variations) in a single call. It leverages the optimized matrix Akima spline implementation to interpolate all columns at once, providing significant performance improvements over column-by-column processing.
Arguments
k_input::AbstractVector: Input wavenumber grid.k_output::AbstractVector: Output wavenumber grid.mono::AbstractMatrix: Monopole moments with shape(n_k, n_cols).quad::AbstractMatrix: Quadrupole moments with shape(n_k, n_cols).hexa::AbstractMatrix: Hexadecapole moments with shape(n_k, n_cols).q_par: Parallel AP parameter.q_perp: Perpendicular AP parameter.
Keyword Arguments
n_GL_points::Int: Number of Gauss-Lobatto points. Default: 8.
Returns
A tuple (mono_AP, quad_AP, hexa_AP) where each is a matrix of shape (n_k_output, n_cols) containing the AP-corrected multipoles for all input columns.
Details
This optimized implementation uses matrix Akima interpolation to process all columns simultaneously rather than iterating column-by-column. For N columns, this reduces 3×N Akima calls to just 3 matrix Akima calls, providing a ~2-3× speedup.
This is particularly useful for computing Jacobians where each column represents the derivative with respect to a different parameter (typically 11 bias parameters).
Performance
For typical DESI-like scenarios (50 input k-points, 100 output k-points, 11 columns):
- Old implementation: ~6 ms (33 scalar Akima calls)
- New implementation: ~2 ms (3 matrix Akima calls)
- Speedup: ~2.5-3×
See Also
apply_AP(k_input::AbstractVector, k_output::AbstractVector, mono::AbstractVector, quad::AbstractVector, hexa::AbstractVector, q_par, q_perp): Single-column version._akima_interpolation(u::AbstractMatrix, t, t_new): Optimized matrix Akima interpolation.
Effort.window_convolution — Methodwindow_convolution(W::Array{T, 4}, v::Matrix) where {T}Applies a 4-dimensional window function or kernel W to a 2-dimensional input matrix v.
This operation performs a transformation or generalized convolution, summing over the j and l indices of the inputs to produce a 2D result indexed by i and k. This is commonly used in analyses where a 4D kernel relates input data in two dimensions to output data in another two dimensions.
Arguments
W: A 4-dimensional array representing the window function or kernel.v: A 2-dimensional matrix representing the input data.
Returns
A 2-dimensional matrix representing the result of the convolution or transformation.
Details
The function implements the summation using the @tullio macro, which provides an efficient way to express tensor contractions and generalized convolutions. The operation can be thought of as applying a 4D kernel to a 2D input, resulting in a 2D output.
Formula
The operation is defined as:
\[C_{ik} = \sum_{j,l} W_{ijkl} v_{jl}\]
See Also
window_convolution(W::AbstractMatrix, v::AbstractVector): Method for a matrix kernel and vector input.
References
- The methodology for this type of window measurement is discussed in: arXiv:1810.05051
Effort.window_convolution — Methodwindow_convolution(W::AbstractMatrix, v::AbstractVector)Performs matrix-vector multiplication, where the matrix W acts as a linear transformation or window applied to the vector input v.
Arguments
W: An abstract matrix representing the linear transformation or window.v: An abstract vector representing the input data.
Returns
An abstract vector representing the result of the matrix-vector multiplication.
Details
This method is a direct implementation of standard matrix-vector multiplication. It applies the linear transformation defined by matrix W to the vector v.
Formula
The operation is defined as:
\[\mathbf{c} = \mathbf{W} \mathbf{v}\]
or element-wise:
\[c_i = \sum_j W_{ij} v_j\]
See Also
window_convolution(W::Array{T, 4}, v::Matrix) where {T}: Method for a 4D kernel and matrix input.
References
- The methodology for this type of window measurement is discussed in: arXiv:1810.05051