External API
This section documents the public functions intended for users.
Power Spectrum Computation
Effort.get_Pℓ — Function
get_Pℓ(cosmology::Array, D, bs::Array, cosmoemu::AbstractPℓEmulators; stoch_kwargs...)Compute the power spectrum multipole $P_\ell(k)$ given cosmological parameters, bias parameters, and growth factor.
Arguments
cosmology::Array: Array of cosmological parameters (format depends on the emulator training).D: Growth factor value at the redshift of interest.bs::Array: Array of bias parameters.cosmoemu::AbstractPℓEmulators: The multipole emulator containing P11, Ploop, Pct components.
Keyword Arguments
stoch_kwargs...: Additional keyword arguments passed to the stochastic model (e.g., shot noise parameters).
Returns
- Power spectrum multipole values evaluated on the emulator's k-grid.
Details
This function computes the power spectrum by:
- Evaluating the P11, Ploop, and Pct components using the trained neural network emulators.
- Computing the stochastic contribution.
- Combining all components via the bias combination function.
Effort.get_Pℓ_jacobian — Function
get_Pℓ_jacobian(cosmology::Array, D, bs::Array, cosmoemu::AbstractPℓEmulators; stoch_kwargs...)Compute both the power spectrum multipole $P_\ell(k)$ and its Jacobian with respect to bias parameters.
Arguments
cosmology::Array: Array of cosmological parameters (format depends on the emulator training).D: Growth factor value at the redshift of interest.bs::Array: Array of bias parameters.cosmoemu::AbstractPℓEmulators: The multipole emulator containing P11, Ploop, Pct components.
Keyword Arguments
stoch_kwargs...: Additional keyword arguments passed to the stochastic model (e.g., shot noise parameters).
Returns
A tuple (Pℓ, ∂Pℓ_∂b) where:
Pℓ: Power spectrum multipole values evaluated on the emulator's k-grid.∂Pℓ_∂b: Jacobian matrix of the power spectrum with respect to bias parameters.
Details
This function is optimized for inference workflows where both the power spectrum and its derivatives are needed (e.g., gradient-based MCMC, Fisher forecasts). It computes both quantities in a single pass, avoiding redundant neural network evaluations.
The Jacobian is computed using the analytical derivative of the bias combination function, which is significantly faster than automatic differentiation for this specific operation.
See Also
get_Pℓ: Compute only the power spectrum without Jacobian.
Alcock-Paczynski Corrections
Effort.q_par_perp — Method
q_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 — Function
apply_AP(k_input::AbstractVector, k_output::AbstractVector, mono::AbstractVector, quad::AbstractVector, hexa::AbstractVector, q_par, q_perp; n_GL_points=8, method::InterpolationMethod=Cubic())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.
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.method: Interpolation method to use (Akima or Cubic). Default:Cubic.
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.
apply_AP(k_input::AbstractVector, k_output::AbstractVector, mono::AbstractMatrix, quad::AbstractMatrix, hexa::AbstractMatrix, q_par, q_perp; n_GL_points=8, method::InterpolationMethod=Cubic())Batch version of apply_AP for processing multiple columns simultaneously using optimized matrix interpolation.
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.method: Interpolation method to use (Akima or Cubic). Default:Cubic.
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.
Window Convolution
Effort.window_convolution — Method
window_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 — Method
window_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