External API

This section documents the public functions intended for users.

Effort.w0waCDMCosmologyType
w0waCDMCosmology(ln10Aₛ::Number, nₛ::Number, h::Number, ωb::Number, ωc::Number, mν::Number=0., w0::Number=-1., wa::Number=0.)

This struct contains the value of the cosmological parameters for $w_0 w_a$CDM cosmologies.

Keyword arguments

  • ln10Aₛ and nₛ, the amplitude and the tilt of the primordial power spectrum fluctuations
  • h, the value of the reduced Hubble paramater
  • ωb and ωc, the physical energy densities of baryons and cold dark matter
  • , the sum of the neutrino masses in eV
  • w₀ and wₐ, the Dark Energy equation of state parameters in the CPL parameterization
source
Effort.q_par_perpMethod
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: An AbstractCosmology struct representing the varying cosmology (e.g., from an MCMC chain).
  • cosmo_ref: An AbstractCosmology struct 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.

The Hubble parameter E(z) is calculated using _E_z, and the conformal angular diameter distance d̃_A(z) is calculated using _d̃A_z.

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.

See Also

  • _E_z: Calculates the normalized Hubble parameter.
  • _d̃A_z: Calculates the conformal angular diameter distance.
  • AbstractCosmology: The abstract type for cosmology structs.
source
Effort.apply_APFunction
apply_AP(k_input::Array, k_output::Array, mono::Array, quad::Array, hexa::Array, 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 and redshift-space distortions (RSD) to the power spectrum multipoles, designed for performance compared to the check version using generic numerical integration.

Arguments

  • k_input: An array of wavenumber values on which the input true multipole moments (mono, quad, hexa) are defined.
  • k_output: An array of observed wavenumber values at which to calculate the output observed multipoles.
  • mono: An array containing the values of the true monopole moment $I_0(k)$ on the k_input grid.
  • quad: An array containing the values of the true quadrupole moment $I_2(k)$ on the k_input grid.
  • hexa: An array containing the values of the true hexadecapole moment $I_4(k)$ on the k_input grid.
  • 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 to 2 * n_GL_points. Defaults to 8.

Returns

A tuple (P0_obs, P2_obs, P4_obs), where each element is an array 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:

  1. Determine Gauss-Lobatto nodes and weights for the interval [0, 1].
  2. For each observed wavenumber k_o in the input k_output array and each μ_o node: 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_spline, interpolating from the k_input grid to the new k_t values. 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)$.
  3. Perform the weighted sum (quadrature) over the μ_o nodes to get the observed multipoles $P_\ell(k_o)$ on the k_output grid.

This function is the standard, performant implementation for applying AP and RSD 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_spline: Interpolates the true multipole moments.
  • _Pk_recon: Reconstructs the true power spectrum on a grid.
  • gausslobatto: Function used to get quadrature nodes and weights.
source
Effort.window_convolutionMethod
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

References

  • The methodology for this type of window measurement is discussed in: arXiv:1810.05051
source
Effort.window_convolutionMethod
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

References

  • The methodology for this type of window measurement is discussed in: arXiv:1810.05051
source