API reference

This section documents the functions intended for internal usage by the package.

Index

Background

Effort._FFunction
_F(y)

Arguments

  • y: The value of the parameter y for which the integral is calculated.

Returns

The value of the definite integral for the given y.

Details

The integrand is defined as: $f(x, y) = x^2 \cdot \sqrt{x^2 + y^2} / (1 + e^x)$

The integration is performed over the domain (0, Inf) for the variable x. A relative tolerance of 1e-12 is used for the integration solver.

source
Effort._dFdyFunction
_dFdy(y)

Calculates the definite integral of the function $f(x, y) = x^2 / ((1 + e^x) \cdot \sqrt{x^2 + y^2})$ with respect to x from 0 to Inf, and then multiplies the result by y.

This function is the derivative of the integral function _F(y) with respect to y.

Arguments

  • y: The value of the parameter y used in the integrand and as a multiplicative factor.

Returns

The value of the definite integral multiplied by y for the given y.

source
Effort._get_yFunction
_get_y(mν, a; kB=8.617342e-5, Tν=0.71611 * 2.7255)

Calculates the dimensionless parameter y used in the integral function _F(y).

The parameter y is calculated based on the neutrino mass, scale factor, Boltzmann constant, and neutrino temperature according to the formula:

y = mν * a / (kB * Tν)

Arguments

  • : Neutrino mass (in units where kB and are defined).
  • a: Scale factor.

Keyword Arguments

  • kB: Boltzmann constant (default: 8.617342e-5 eV/K).
  • : Neutrino temperature (default: 0.71611 * 2.7255 K).

Returns

The calculated dimensionless parameter y.

source
Effort._ΩνE2Method
_ΩνE2(a, Ωγ0, mν; kB=8.617342e-5, Tν=0.71611 * 2.7255, Neff=3.044)

Calculates the energy density of relic neutrinos, scaled by the critical density, at a given scale factor a, for a single neutrino mass.

This function accounts for the contribution of a single neutrino mass to the total energy density. It uses _F(y) to incorporate the effect of neutrino mass and temperature.

Arguments

  • a: The scale factor.
  • Ωγ0: The photon density parameter today.
  • : The neutrino mass (a single value).

Keyword Arguments

  • kB: Boltzmann constant (default: 8.617342e-5 eV/K).
  • : Neutrino temperature (default: 0.71611 * 2.7255 K).
  • Neff: Effective number of neutrino species (default: 3.044).

Returns

The calculated neutrino energy density parameter ΩνE2 at scale factor a for the given mass.

Details

The calculation involves a factor Γν derived from Neff and the ratio of neutrino to photon temperatures. The main term is proportional to Ωγ0 / a^4 multiplied by F_interpolant(_get_y(mν, a)).

The parameter y passed to F_interpolant is calculated using _get_y(mν, a).

Formula

The formula used is: ΩνE2 = (15 / π^4) * Γν^4 * (Ωγ0 / a^4) * F(y) where Γν = (4/11)^(1/3) * (Neff/3)^(1/4) and y = mν * a / (kB * Tν).

See Also

source
Effort._ΩνE2Method
_ΩνE2(a, Ωγ0, mν::AbstractVector; kB=8.617342e-5, Tν=0.71611 * 2.7255, Neff=3.044)

Calculates the energy density of relic neutrinos, scaled by the critical density, at a given scale factor a, for a vector of neutrino masses.

This function accounts for the combined contribution of multiple neutrino masses to the total energy density by summing the individual contributions. It uses the F_interpolant function (which is equivalent to _F(y)) for each mass.

Arguments

  • a: The scale factor.
  • Ωγ0: The photon density parameter today.
  • : A vector of neutrino masses (AbstractVector).

Keyword Arguments

  • kB: Boltzmann constant (default: 8.617342e-5 eV/K).
  • : Neutrino temperature (default: 0.71611 * 2.7255 K).
  • Neff: Effective number of neutrino species (default: 3.044).

Returns

The calculated total neutrino energy density parameter ΩνE2 at scale factor a for the sum of contributions from all masses in the vector.

Details

The calculation involves a factor Γν derived from Neff and the ratio of neutrino to photon temperatures. The main term is proportional to Ωγ0 / a^4 multiplied by the sum of F_interpolant(_get_y(mν_i, a)) for each mass mν_i in the input vector .

The parameter y passed to F_interpolant for each mass is calculated using _get_y(mν_i, a).

Formula

The formula used is: ΩνE2 = (15 / π^4) * Γν^4 * (Ωγ0 / a^4) * Σ F(y_i) where Γν = (4/11)^(1/3) * (Neff/3)^(1/4) and y_i = mν_i * a / (kB * Tν).

See Also

source
Effort._dΩνE2daMethod
_dΩνE2da(a, Ωγ0, mν; kB=8.617342e-5, Tν=0.71611 * 2.7255, Neff=3.044)

Calculates the derivative of the neutrino energy density parameter _ΩνE2 with respect to the scale factor a, for a single neutrino mass.

This function computes the derivative of the expression for _ΩνE2 by applying the chain rule, involving both _F(y) and _dFdy(y) functions.

Arguments

  • a: The scale factor.
  • Ωγ0: The photon density parameter today.
  • : The neutrino mass (a single value).

Keyword Arguments

  • kB: Boltzmann constant (default: 8.617342e-5 eV/K).
  • : Neutrino temperature (default: 0.71611 * 2.7255 K).
  • Neff: Effective number of neutrino species (default: 3.044).

Returns

The calculated derivative d(ΩνE2)/da at scale factor a for the given mass.

Details

The calculation is based on the derivative of the _ΩνE2 formula with respect to a.

See Also

source
Effort._dΩνE2daMethod
_dΩνE2da(a, Ωγ0, mν::AbstractVector; kB=8.617342e-5, Tν=0.71611 * 2.7255, Neff=3.044)

Calculates the derivative of the neutrino energy density parameter _ΩνE2 with respect to the scale factor a, for a vector of neutrino masses.

This function computes the derivative of the expression for _ΩνE2 by summing the derivatives of the contributions from each individual neutrino mass. It uses the _F(y) and _dFdy(y) functions for each mass.

Arguments

  • a: The scale factor.
  • Ωγ0: The photon density parameter today.
  • : A vector of neutrino masses (AbstractVector).

Keyword Arguments

  • kB: Boltzmann constant (default: 8.617342e-5 eV/K).
  • : Neutrino temperature (default: 0.71611 * 2.7255 K).
  • Neff: Effective number of neutrino species (default: 3.044).

Returns

The calculated total derivative d(ΩνE2)/da at scale factor a for the sum of contributions from all masses in the vector.

Details

The calculation sums the derivatives of the individual neutrino mass contributions to _ΩνE2 with respect to a.

See Also

source
Effort._a_zFunction
_a_z(z)

Calculates the cosmological scale factor a from the redshift z.

The relationship between scale factor and redshift is given by $a = 1 / (1 + z)$.

Arguments

  • z: The redshift (scalar or array).

Returns

The corresponding scale factor a (scalar or array).

Formula

The formula used is: $a = 1 / (1 + z)$

source
Effort._ρDE_aFunction
_ρDE_a(a, w0, wa)

Calculates the evolution of the dark energy density parameter relative to its value today, as a function of the scale factor a.

This function implements the standard parametrization for the dark energy equation of state w(a) = w0 + wa*(1-a).

Arguments

  • a: The scale factor (scalar or array).
  • w0: The present-day value of the dark energy equation of state parameter.
  • wa: The derivative of the dark energy equation of state parameter with respect to (1-a).

Returns

The dark energy density parameter relative to its value today, ρ_DE(a) / ρ_DE(a=1), at the given scale factor a (scalar or array).

Formula

The formula used is: $\rho_\mathrm{DE}(a) / \rho_\mathrm{DE}(a=1) = a^(-3 * (1 + w0 + wa)) * e^{3 * wa * (a - 1)}$

This function uses broadcasting (@.) to handle scalar or array inputs for a.

See Also

  • _ρDE_z(z, w0, wa): Calculates the dark energy density evolution as a function of redshift z.
source
Effort._ρDE_zFunction
_ρDE_z(z, w0, wa)

Calculates the evolution of the dark energy density parameter relative to its value today, as a function of the redshift z.

This function implements the standard parametrization for the dark energy equation of state $w(a) = w0 + wa*(1-a)$, converted to depend on redshift z.

Arguments

  • z: The redshift (scalar or array).
  • w0: The present-day value of the dark energy equation of state parameter.
  • wa: The derivative of the dark energy equation of state parameter with respect to (1-a).

Returns

The dark energy density parameter relative to its value today, ρ_DE(z) / ρ_DE(z=0), at the given redshift z (scalar or array).

Formula

The formula used is: $\rho_\mathrm{DE}(z) / \rho_\mathrm{DE}(z=0) = (1 + z)^(3 * (1 + w0 + wa)) * e^{-3 * wa * z / (1 + z)}$

This function uses broadcasting (@.) to handle scalar or array inputs for z.

See Also

  • _ρDE_a(a, w0, wa): Calculates the dark energy density evolution as a function of scale factor a.
source
Effort._dρDEdaFunction
_dρDEda(a, w0, wa)

Calculates the derivative of the dark energy density parameter evolution, d(ρ_DE(a)/ρ_DE(a=1))/da, with respect to the scale factor a.

This function computes the derivative of the formula implemented in _ρDE_a(a, w0, wa).

Arguments

  • a: The scale factor (scalar or array).
  • w0: The present-day value of the dark energy equation of state parameter.
  • wa: The derivative of the dark energy equation of state parameter with respect to (1-a).

Returns

The calculated derivative of the dark energy density parameter evolution with respect to a at the given scale factor a (scalar or array).

Formula

The formula used is: $\frac{d}{da} \left( \frac{\rho_{\text{DE}}(a)}{\rho_{\text{DE}}(a=1)} \right) = 3 \left( -\frac{1 + w_0 + w_a}{a} + w_a \right) \frac{\rho_{\text{DE}}(a)}{\rho_{\text{DE}}(a=1)}$

This function uses broadcasting (@.) to handle scalar or array inputs for a.

See Also

source
Effort._E_aMethod
_E_a(a, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the normalized Hubble parameter, E(a), at a given scale factor a.

E(a) describes the expansion rate of the universe relative to the Hubble constant today, incorporating contributions from different energy density components: radiation (photons and massless neutrinos), cold dark matter and baryons, dark energy, and massive neutrinos.

Arguments

  • a: The scale factor (scalar or array).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : The total neutrino mass (or a vector of masses), used in the calculation of the massive neutrino energy density. Defaults to 0.0 (massless neutrinos).
  • w0: The present-day value of the dark energy equation of state parameter w(a) = w0 + wa*(1-a). Defaults to -1.0 (ΛCDM).
  • wa: The derivative of the dark energy equation of state parameter with respect to (1-a). Defaults to 0.0 (ΛCDM).

Returns

The calculated normalized Hubble parameter E(a) (scalar or array).

Details

The calculation includes:

  • Photon density Ωγ0 = 2.469e-5 / h^2.
  • Massless neutrino density Ων0 (calculated from _ΩνE2 at a=1).
  • Dark energy density ΩΛ0 (calculated to ensure a flat universe: 1 - Ωγ0 - Ωcb0 - Ων0).
  • Massive neutrino density _ΩνE2(a, Ωγ0, mν).
  • Dark energy evolution _ρDE_a(a, w0, wa) (density relative to today's dark energy density).

The formula used is: E(a) = sqrt(Ωγ0 * a^-4 + Ωcb0 * a^-3 + ΩΛ0 * ρDE(a) + ΩνE2(a)) where ρDE(a) is the dark energy density relative to its value today, and ΩνE2(a) is the massive neutrino energy density parameter at scale factor a.

This function uses broadcasting (@.) to handle scalar or array inputs for a.

See Also

source
Effort._E_aMethod
_E_a(a, w0wacosmo::w0waCDMCosmology)

Calculates the normalized Hubble parameter, E(a), at a given scale factor a, using parameters extracted from a w0waCDMCosmology struct.

This method is a convenience wrapper around the main _E_a function. It extracts the cold dark matter and baryon density (Ωcb0), Hubble parameter (h), neutrino mass (), and dark energy parameters (w0, wa) from the provided cosmology struct and passes them to the primary _E_a method.

Arguments

  • a: The scale factor (scalar or array).
  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

The calculated normalized Hubble parameter E(a) (scalar or array).

Details

The parameters Ωcb0, h, , w0, and wa are extracted from the w0wacosmo struct. Ωcb0 is calculated as (w0wacosmo.ωb + w0wacosmo.ωc) / w0wacosmo.h^2.

This method calls the primary _E_a(a, Ωcb0, h; mν, w0, wa) method internally.

See Also

  • _E_a(a, Ωcb0, h; mν, w0, wa): The primary method that performs the calculation.
  • w0waCDMCosmology: The struct type containing the cosmological parameters.
source
Effort._E_zMethod
_E_z(z, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the normalized Hubble parameter, E(z), as a function of redshift z.

This function is the redshift-dependent counterpart to _E_a(a, Ωcb0, h; mν, w0, wa). It first converts z to the scale factor a using _a_z(z) and then calls the _E_a function.

Arguments

  • z: The redshift (scalar or array).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

The calculated normalized Hubble parameter E(z) (scalar or array).

See Also

source
Effort._E_zMethod
_E_z(z, w0wacosmo::w0waCDMCosmology)

Calculates the normalized Hubble parameter, E(z), as a function of redshift z, using parameters extracted from a w0waCDMCosmology struct.

This function is the redshift-dependent counterpart to _E_a(a, w0wacosmo::w0waCDMCosmology). It's a convenience method that extracts parameters from the struct and calls the primary _E_z(z, Ωcb0, h; mν, w0, wa) method.

Arguments

  • z: The redshift (scalar or array).
  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

The calculated normalized Hubble parameter E(z) (scalar or array).

See Also

source
Effort._dlogEdlogaFunction
_dlogEdloga(a, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the logarithmic derivative of the normalized Hubble parameter, $\frac{d(\log E)}{d(\log a)}$, with respect to the logarithm of the scale factor a.

This quantity is useful in cosmological calculations, particularly when analyzing the growth of structure. It is derived from the derivative of E(a) with respect to a.

Arguments

  • a: The scale factor (scalar or array).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

The calculated value of $\frac{d(\log E)}{d(\log a)}$ at the given scale factor a (scalar or array).

Details

The calculation involves the derivative of the _E_a function with respect to a. The formula is derived from $\frac{d(\log E)}{d(\log a)} = \frac{a}{E} \frac{dE}{da}$. The derivative dE/da involves terms related to the derivatives of the density components with respect to a, including _dρDEda(a, w0, wa) and _dΩνE2da(a, Ωγ0, mν).

This function uses broadcasting (@.) to handle scalar or array inputs for a.

See Also

source
Effort._ΩmaMethod
_Ωma(a, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the total matter density parameter, Ω_m(a), at a given scale factor a.

This represents the combined density of cold dark matter and baryons relative to the critical density at scale factor a.

Arguments

  • a: The scale factor (scalar or array).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es) (used in the calculation of _E_a).
  • w0: Dark energy equation of state parameter (used in the calculation of _E_a).
  • wa: Dark energy equation of state parameter derivative (used in the calculation of _E_a).

Returns

The calculated total matter density parameter Ω_m(a) at the given scale factor a (scalar or array).

Formula

The formula used is: $\Omega_m(a) = \frac{\Omega_{\text{cb}0} a^{-3}}{E(a)^2}$ where $E(a)$ is the normalized Hubble parameter calculated using _E_a(a, Ωcb0, h; mν, w0, wa).

This function uses broadcasting (@.) to handle scalar or array inputs for a.

See Also

source
Effort._ΩmaMethod
_Ωma(a, w0wacosmo::w0waCDMCosmology)

Calculates the total matter density parameter, Ω_m(a), at a given scale factor a, using parameters extracted from a w0waCDMCosmology struct.

This method is a convenience wrapper around the primary _Ωma(a, Ωcb0, h; mν, w0, wa) function. It extracts the cold dark matter and baryon density (Ωcb0), Hubble parameter (h), neutrino mass (), and dark energy parameters (w0, wa) from the provided cosmology struct and passes them to the primary _Ωma method.

Arguments

  • a: The scale factor (scalar or array).
  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

The calculated total matter density parameter Ω_m(a) at the given scale factor a (scalar or array).

Details

The parameters Ωcb0, h, , w0, and wa are extracted from the w0wacosmo struct. Ωcb0 is calculated as (w0wacosmo.ωb + w0wacosmo.ωc) / w0wacosmo.h^2.

This method calls the primary _Ωma(a, Ωcb0, h; mν, w0, wa) method internally.

See Also

  • _Ωma(a, Ωcb0, h; mν, w0, wa): The primary method that performs the calculation.
  • w0waCDMCosmology: The struct type containing the cosmological parameters.
source
Effort._r̃_z_checkFunction
_r̃_z_check(z, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the conformal distance r̃(z) to a given redshift z using numerical integration.

This is a "check" version, typically slower but potentially more accurate, used for verifying results from faster methods. The conformal distance is the integral of 1/E(z) with respect to z.

Arguments

  • z: The redshift (scalar).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

The calculated conformal distance r̃(z) (scalar).

Details

The function calculates the integral $\int_0^z \frac{dz'}{E(z')}$ where $E(z')$ is the normalized Hubble parameter at redshift $z'$, calculated using _E_a after converting $z'$ to scale factor using _a_z. The integration is performed using IntegralProblem and the QuadGKJL() solver.

Formula

The conformal distance is defined as: $\tilde{r}(z) = \int_0^z \frac{dz'}{E(z')}$

See Also

  • _r̃_z: The standard, faster method for calculating conformal distance.
  • _E_a: Calculates the normalized Hubble parameter as a function of scale factor.
  • _a_z: Converts redshift to scale factor.
  • _r_z_check: Calculates the comoving distance using this check version.
source
Effort._r̃_zMethod
_r̃_z(z, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the conformal distance r̃(z) to a given redshift z using Gauss-Legendre quadrature.

This is the standard, faster method for calculating the conformal distance, which is the integral of 1/E(z) with respect to z.

Arguments

  • z: The redshift (scalar or array).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

The calculated conformal distance r̃(z) (scalar or array).

Details

The function approximates the integral $\int_0^z \frac{dz'}{E(z')}$ using Gauss-Legendre quadrature with a specified number of points (here, 9). It uses _transformed_weights to get the quadrature points and weights over the interval [0, z]. The integrand $1/E(z')$ is evaluated at these points using _E_a (after converting z' to a with _a_z), and the result is a weighted sum.

Formula

The conformal distance is defined as: $\tilde{r}(z) = \int_0^z \frac{dz'}{E(z')}$ This function computes this integral numerically.

See Also

source
Effort._r̃_zMethod
_r̃_z(z, w0wacosmo::w0waCDMCosmology)

Calculates the conformal distance r̃(z) to a given redshift z, using parameters extracted from a w0waCDMCosmology struct.

This method is a convenience wrapper around the primary _r̃_z(z, Ωcb0, h; mν, w0, wa) function. It extracts the necessary cosmological parameters from the provided struct.

Arguments

  • z: The redshift (scalar or array).
  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

The calculated conformal distance r̃(z) (scalar or array).

Details

The parameters Ωcb0, h, , w0, and wa are extracted from the w0wacosmo struct. Ωcb0 is calculated as (w0wacosmo.ωb + w0wacosmo.ωc) / w0wacosmo.h^2.

This method calls the primary _r̃_z(z, Ωcb0, h; mν, w0, wa) method internally.

See Also

  • _r̃_z(z, Ωcb0, h; mν, w0, wa): The primary method for calculating conformal distance.
  • w0waCDMCosmology: The struct type containing the cosmological parameters.
  • _r_z: Calculates the comoving distance.
source
Effort._r_z_checkMethod
_r_z_check(z, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the comoving distance r(z) to a given redshift z using the "check" version of the conformal distance calculation.

The comoving distance is related to the conformal distance by a factor involving the speed of light and the Hubble parameter today. This version uses the slower, potentially more accurate _r̃_z_check for the conformal distance.

Arguments

  • z: The redshift (scalar).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

The calculated comoving distance r(z) (scalar).

Details

The comoving distance is calculated by scaling the conformal distance obtained from _r̃_z_check(z, Ωcb0, h; mν, w0, wa) by the factor $c_0 / (100 h)$, where $c_0$ is the speed of light (in units consistent with h).

Formula

The comoving distance is defined as: $r(z) = \frac{c_0}{100 h} \tilde{r}(z)$ This function uses $\tilde{r}(z) = \text{_r̃_z_check}(z, \dots)$.

See Also

  • _r̃_z_check: The slower, check version of the conformal distance calculation.
  • _r_z: The standard, faster method for calculating comoving distance.
source
Effort._r_zMethod
_r_z(z, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the comoving distance r(z) to a given redshift z using the standard conformal distance calculation.

The comoving distance is related to the conformal distance by a factor involving the speed of light and the Hubble parameter today. This version uses the standard, faster _r̃_z for the conformal distance.

Arguments

  • z: The redshift (scalar or array).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

The calculated comoving distance r(z) (scalar or array).

Details

The comoving distance is calculated by scaling the conformal distance obtained from _r̃_z(z, Ωcb0, h; mν, w0, wa) by the factor $c_0 / (100 h)$, where $c_0$ is the speed of light (in units consistent with h).

Formula

The comoving distance is defined as: $r(z) = \frac{c_0}{100 h} \tilde{r}(z)$ This function uses $\tilde{r}(z) = \text{_r̃_z}(z, \dots)$.

See Also

source
Effort._r_zMethod
_r_z(z, w0wacosmo::w0waCDMCosmology)

Calculates the comoving distance r(z) to a given redshift z, using parameters extracted from a w0waCDMCosmology struct.

This method is a convenience wrapper around the primary _r_z(z, Ωcb0, h; mν, w0, wa) function. It extracts the necessary cosmological parameters from the provided struct.

Arguments

  • z: The redshift (scalar or array).
  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

The calculated comoving distance r(z) (scalar or array).

Details

The parameters Ωcb0, h, , w0, and wa are extracted from the w0wacosmo struct. Ωcb0 is calculated as (w0wacosmo.ωb + w0wacosmo.ωc) / w0wacosmo.h^2.

This method calls the primary _r_z(z, Ωcb0, h; mν, w0, wa) method internally.

See Also

  • _r_z(z, Ωcb0, h; mν, w0, wa): The primary method for calculating comoving distance.
  • w0waCDMCosmology: The struct type containing the cosmological parameters.
  • _r̃_z: Calculates the conformal distance.
source
Effort._d̃A_zMethod
_d̃A_z(z, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the conformal angular diameter distance d̃_A(z) to a given redshift z.

The conformal angular diameter distance is defined as the conformal comoving distance divided by (1 + z).

Arguments

  • z: The redshift (scalar or array).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

The calculated conformal angular diameter distance d̃_A(z) (scalar or array).

Details

The function calculates the conformal comoving distance using _r̃_z(z, Ωcb0, h; mν, w0, wa) and then divides by (1 + z).

Formula

The formula used is: $\tilde{d}_A(z) = \frac{\tilde{r}(z)}{1 + z}$ where $\tilde{r}(z)$ is the conformal comoving distance.

See Also

source
Effort._d̃A_zMethod
_d̃A_z(z, w0wacosmo::w0waCDMCosmology)

Calculates the conformal angular diameter distance d̃_A(z) to a given redshift z, using parameters extracted from a w0waCDMCosmology struct.

This method is a convenience wrapper around the primary _d̃A_z(z, Ωcb0, h; mν, w0, wa) function. It extracts the necessary cosmological parameters from the provided struct.

Arguments

  • z: The redshift (scalar or array).
  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

The calculated conformal angular diameter distance d̃_A(z) (scalar or array).

Details

The parameters Ωcb0, h, , w0, and wa are extracted from the w0wacosmo struct. Ωcb0 is calculated as (w0wacosmo.ωb + w0wacosmo.ωc) / w0wacosmo.h^2.

This method calls the primary _d̃A_z(z, Ωcb0, h; mν, w0, wa) method internally.

See Also

  • _d̃A_z(z, Ωcb0, h; mν, w0, wa): The primary method for calculating conformal angular diameter distance.
  • w0waCDMCosmology: The struct type containing the cosmological parameters.
  • _dA_z: Calculates the standard angular diameter distance.
source
Effort._dA_zMethod
_dA_z(z, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the angular diameter distance d_A(z) to a given redshift z.

The angular diameter distance is defined as the comoving distance divided by (1 + z).

Arguments

  • z: The redshift (scalar or array).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

The calculated angular diameter distance d_A(z) (scalar or array).

Details

The function calculates the comoving distance using _r_z(z, Ωcb0, h; mν, w0, wa) and then divides by (1 + z).

Formula

The formula used is: $d_A(z) = \frac{r(z)}{1 + z}$ where $r(z)$ is the comoving distance.

See Also

source
Effort._dA_zMethod
_dA_z(z, w0wacosmo::w0waCDMCosmology)

Calculates the angular diameter distance d_A(z) to a given redshift z, using parameters extracted from a w0waCDMCosmology struct.

This method is a convenience wrapper around the primary _dA_z(z, Ωcb0, h; mν, w0, wa) function. It extracts the necessary cosmological parameters from the provided struct.

Arguments

  • z: The redshift (scalar or array).
  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

The calculated angular diameter distance d_A(z) (scalar or array).

Details

The parameters Ωcb0, h, , w0, and wa are extracted from the w0wacosmo struct. Ωcb0 is calculated as (w0wacosmo.ωb + w0wacosmo.ωc) / w0wacosmo.h^2.

This method calls the primary _dA_z(z, Ωcb0, h; mν, w0, wa) method internally.

See Also

  • _dA_z(z, Ωcb0, h; mν, w0, wa): The primary method for calculating angular diameter distance.
  • w0waCDMCosmology: The struct type containing the cosmological parameters.
  • _d̃A_z: Calculates the conformal angular diameter distance.
source
Effort._growth!Function
_growth!(du, u, p, loga)

Defines the in-place right-hand side of the second-order ordinary differential equation for the linear growth factor D(a), with log(a) as the independent variable.

The state vector u is [D(log a), dD/d(log a)]. This function calculates the derivatives du = [dD/d(log a), d^2D/d(log a)^2] based on the growth equation.

Arguments

  • du: The output vector where the calculated derivatives are stored (modified in-place).
  • u: The current state vector [D(log a), dD/d(log a)].
  • p: A vector of parameters [Ωcb0, mν, h, w0, wa].
  • loga: The natural logarithm of the scale factor, log(a).

Returns

Modifies the du vector in-place.

Details

The function solves the second-order differential equation for the linear growth factor, often written as:

\[\frac{d^2 D}{d(\ln a)^2} + \left(2 + \frac{d \ln E}{d \ln a}\right) \frac{d D}{d \ln a} - \frac{3}{2} \Omega_m(a) D = 0\]

where $E(a)$ is the normalized Hubble parameter and $\Omega_m(a)$ is the matter density parameter.

The terms $\frac{d \ln E}{d \ln a}$ and $\Omega_m(a)$ are calculated using _dlogEdloga and _Ωma respectively, with parameters extracted from p.

The system of first-order ODEs implemented is: $\frac{d u[1]}{d(\ln a)} = u[2]$ $\frac{d u[2]}{d(\ln a)} = -\left(2 + \frac{d \ln E}{d \ln a}\right) u[2] + \frac{3}{2} \Omega_m(a) u[1]$

See Also

  • _growth_solver: Functions that solve this ODE.
  • _dlogEdloga: Calculates the logarithmic derivative of E(a).
  • _Ωma: Calculates the matter density parameter at scale factor a.
  • _E_a: Related normalized Hubble parameter.
source
Effort._growth_solverMethod
_growth_solver(Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Solves the ODE for the linear growth factor D(a) and its derivative dD/d(log a) over a fixed range of log(a), typically from an early time to slightly past a=1.

This function sets up and solves the _growth! ODE using a standard solver.

Arguments

  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

A DifferentialEquations.jl solution object containing the values of D(log a) and dD/d(log a) over the solved log(a) range.

Details

The ODE is solved from log(amin) to log(1.01), where amin = 1/139. Initial conditions u₀ = [amin, amin] are used, corresponding to D(a) ≈ a at early times. The problem is solved using the Tsit5() solver with a relative tolerance of 1e-5.

See Also

source
Effort._growth_solverMethod
_growth_solver(w0wacosmo::w0waCDMCosmology)

Solves the ODE for the linear growth factor D(a) and its derivative dD/d(log a) using parameters extracted from a w0waCDMCosmology struct.

This method is a convenience wrapper around the primary _growth_solver(Ωcb0, h; mν, w0, wa) function. It extracts the necessary cosmological parameters from the provided struct.

Arguments

  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

A DifferentialEquations.jl solution object containing the values of D(log a) and dD/d(log a) over the solved log(a) range.

Details

The parameters Ωcb0, h, , w0, and wa are extracted from the w0wacosmo struct. Ωcb0 is calculated as (w0wacosmo.ωb + w0wacosmo.ωc) / w0wacosmo.h^2.

This method calls the primary _growth_solver(Ωcb0, h; mν, w0, wa) method internally.

See Also

source
Effort._growth_solverMethod
_growth_solver(z, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Solves the ODE for the linear growth factor D(a) and its derivative dD/d(log a) and returns the solution evaluated specifically at the given redshift(s) z.

This function solves the _growth! ODE over a range of log(a) and then extracts the solution values corresponding to the provided redshift(s).

Arguments

  • z: The redshift or an array of redshifts at which to save the solution.
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

A 2xN array (where N is the number of redshifts in z) containing the solution. The first row contains the growth factor D(z), and the second row contains the derivative dD/d(log a) evaluated at redshift z.

Details

The ODE is solved from log(amin) to log(1.01), where amin = 1/139. Initial conditions u₀ = [amin, amin] are used, corresponding to D(a) ≈ a at early times. The problem is solved using the Tsit5() solver with a relative tolerance of 1e-5. The solution is saved specifically at the log(a) values corresponding to the input redshifts z, obtained using _a_z.

See Also

source
Effort._growth_solverMethod
_D_z(z, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the linear growth factor D(z) for a single redshift z.

The linear growth factor describes how density perturbations grow in the linear regime of structure formation. It is obtained by solving a second-order ODE.

Arguments

  • z: The redshift (scalar).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

The calculated linear growth factor D(z) (scalar).

Details

This function solves the growth ODE using the _growth_solver(Ωcb0, h; mν, w0, wa) method, which solves over a fixed range of log(a). It then evaluates the solution at the log(a) value corresponding to the input redshift z (obtained via _a_z) to get the value of D(z).

See Also

source
Effort._D_zMethod
_growth_solver(z, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Solves the ODE for the linear growth factor D(a) and its derivative dD/d(log a) and returns the solution evaluated specifically at the given redshift(s) z.

This function solves the _growth! ODE over a range of log(a) and then extracts the solution values corresponding to the provided redshift(s).

Arguments

  • z: The redshift or an array of redshifts at which to save the solution.
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

A 2xN array (where N is the number of redshifts in z) containing the solution. The first row contains the growth factor D(z), and the second row contains the derivative dD/d(log a) evaluated at redshift z.

Details

The ODE is solved from log(amin) to log(1.01), where amin = 1/139. Initial conditions u₀ = [amin, amin] are used, corresponding to D(a) ≈ a at early times. The problem is solved using the Tsit5() solver with a relative tolerance of 1e-5. The solution is saved specifically at the log(a) values corresponding to the input redshifts z, obtained using _a_z.

See Also

source
Effort._D_zMethod
_D_z(z::AbstractVector, Ωcb0, h; mν=0.0, w0=-1.0, wa=0.0)

Calculates the linear growth factor D(z) for a vector of redshifts z.

The linear growth factor describes how density perturbations grow in the linear regime of structure formation. It is obtained by solving a second-order ODE.

Arguments

  • z: A vector of redshifts (AbstractVector).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

A vector containing the calculated linear growth factor D(z) for each redshift in the input vector z.

Details

This function solves the growth ODE using the _growth_solver(z, Ωcb0, h; mν, w0, wa) method, which solves the ODE and saves the solution specifically at the log(a) values corresponding to the input redshifts z. It then extracts the first row of the solution (which contains the D(z) values) and reverses it.

See Also

source
Effort._D_zMethod
_D_z(z, w0wacosmo::w0waCDMCosmology)

Calculates the linear growth factor D(z) for a given redshift or vector of redshifts z, using parameters extracted from a w0waCDMCosmology struct.

This method is a convenience wrapper around the primary _D_z(z, Ωcb0, h; mν, w0, wa) or _D_z(z::AbstractVector, Ωcb0, h; mν, w0, wa) functions. It extracts the necessary cosmological parameters from the provided struct and calls the appropriate method based on whether z is a scalar or a vector.

Arguments

  • z: The redshift or an array of redshifts.
  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

The calculated linear growth factor D(z) (scalar or array).

Details

The parameters Ωcb0, h, , w0, and wa are extracted from the w0wacosmo struct. Ωcb0 is calculated as (w0wacosmo.ωb + w0wacosmo.ωc) / w0wacosmo.h^2.

This method calls either _D_z(z, Ωcb0, h; mν, w0, wa) or _D_z(z::AbstractVector, Ωcb0, h; mν, w0, wa) internally, depending on the type of z.

See Also

source
Effort._f_zMethod
_f_z(z::AbstractVector, Ωcb0, h; mν=0, w0=-1.0, wa=0.0)

Calculates the linear growth rate f(z) for a vector of redshifts z.

The linear growth rate is defined as $f(z) = \frac{d \ln D}{d \ln a} = \frac{dD/d(\ln a)}{D}$, where D(z) is the linear growth factor and a is the scale factor.

Arguments

  • z: A vector of redshifts (AbstractVector).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

A vector containing the calculated linear growth rate f(z) for each redshift in the input vector z.

Details

This function uses the _growth_solver(z, Ωcb0, h; mν, w0, wa) method to solve the growth ODE and obtain the growth factor D(z) and its derivative with respect to log(a), dD/d(log a), at the specified redshifts. It then calculates f(z) as the ratio of dD/d(log a) to D(z) at each redshift. The result is reversed before returning.

Formula

The formula used is: $f(z) = \frac{d \ln D}{d \ln a} = \frac{dD/d(\ln a)}{D}$

See Also

source
Effort._f_zMethod
_f_z(z, Ωcb0, h; mν=0, w0=-1.0, wa=0.0)

Calculates the linear growth rate f(z) for a single redshift z.

The linear growth rate is defined as $f(z) = \frac{d \ln D}{d \ln a} = \frac{dD/d(\ln a)}{D}$, where D(z) is the linear growth factor and a is the scale factor.

Arguments

  • z: The redshift (scalar).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

The calculated linear growth rate f(z) (scalar).

Details

This function uses the _growth_solver(z, Ωcb0, h; mν, w0, wa) method to solve the growth ODE and obtain the growth factor D(z) and its derivative with respect to log(a), dD/d(log a), at the specified redshift. It then calculates f(z) as the ratio of dD/d(log a) to D(z).

Formula

The formula used is: $f(z) = \frac{d \ln D}{d \ln a} = \frac{dD/d(\ln a)}{D}$

See Also

source
Effort._f_zMethod
_f_z(z, w0wacosmo::w0waCDMCosmology)

Calculates the linear growth rate f(z) for a given redshift or vector of redshifts z, using parameters extracted from a w0waCDMCosmology struct.

This method is a convenience wrapper around the primary _f_z(z, Ωcb0, h; mν, w0, wa) or _f_z(z::AbstractVector, Ωcb0, h; mν, w0, wa) functions. It extracts the necessary cosmological parameters from the provided struct and calls the appropriate method based on whether z is a scalar or a vector.

Arguments

  • z: The redshift or an array of redshifts.
  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

The calculated linear growth rate f(z) (scalar or array).

Details

The parameters Ωcb0, h, , w0, and wa are extracted from the w0wacosmo struct. Ωcb0 is calculated as (w0wacosmo.ωb + w0wacosmo.ωc) / w0wacosmo.h^2.

This method calls either _f_z(z, Ωcb0, h; mν, w0, wa) or _f_z(z::AbstractVector, Ωcb0, h; mν, w0, wa) internally, depending on the type of z.

See Also

source
Effort._D_f_zMethod
_D_f_z(z, Ωcb0, h; mν=0, w0=-1.0, wa=0.0)

Calculates both the linear growth factor D(z) and the linear growth rate f(z) for a vector of redshifts z.

This function is a convenience to get both quantities from a single ODE solution. The growth rate is defined as $f(z) = \frac{d \ln D}{d \ln a} = \frac{dD/d(\ln a)}{D}$.

Arguments

  • z: A vector of redshifts (AbstractVector).
  • Ωcb0: The density parameter for cold dark matter and baryons today.
  • h: The Hubble parameter today, divided by 100 km/s/Mpc.

Keyword Arguments

  • : Total neutrino mass(es).
  • w0: Dark energy equation of state parameter.
  • wa: Dark energy equation of state parameter derivative.

Returns

A tuple (D_values, f_values), where D_values is a vector of the linear growth factor D(z) and f_values is a vector of the linear growth rate f(z) for each redshift in the input vector z. Both vectors are reversed before returning.

Details

This function uses the _growth_solver(z, Ωcb0, h; mν, w0, wa) method to solve the growth ODE and obtain the growth factor D(z) and its derivative with respect to log(a), dD/d(log a), at the specified redshifts. It then calculates f(z) as the ratio of dD/d(log a) to D(z) at each redshift. Both the D(z) and calculated f(z) vectors are returned.

Formula

The formula used for f(z) is: $f(z) = \frac{d \ln D}{d \ln a} = \frac{dD/d(\ln a)}{D}$

See Also

source
Effort._D_f_zMethod
_D_f_z(z, w0wacosmo::w0waCDMCosmology)

Calculates both the linear growth factor D(z) and the linear growth rate f(z) for a vector of redshifts z, using parameters extracted from a w0waCDMCosmology struct.

This method is a convenience wrapper around the primary _D_f_z(z, Ωcb0, h; mν, w0, wa) function. It extracts the necessary cosmological parameters from the provided struct.

Arguments

  • z: A vector of redshifts.
  • w0wacosmo: A struct of type w0waCDMCosmology containing the cosmological parameters.

Returns

A tuple (D_values, f_values), where D_values is a vector of the linear growth factor D(z) and f_values is a vector of the linear growth rate f(z) for each redshift in the input vector z.

Details

The parameters Ωcb0, h, , w0, and wa are extracted from the w0wacosmo struct. Ωcb0 is calculated as (w0wacosmo.ωb + w0wacosmo.ωc) / w0wacosmo.h^2.

This method calls the primary _D_f_z(z, Ωcb0, h; mν, w0, wa) method internally.

See Also

  • _D_f_z(z, Ωcb0, h; mν, w0, wa): The primary method for calculating D(z) and f(z).
  • w0waCDMCosmology: The struct type containing the cosmological parameters.
  • _D_z: Calculates the linear growth factor separately.
  • _f_z: Calculates the linear growth rate separately.
source

Projection

Effort._PkμFunction
_Pkμ(k, μ, Int_Mono, Int_Quad, Int_Hexa)

Reconstructs the anisotropic power spectrum $P(k, \mu)$ at a given wavenumber k and cosine of the angle to the line-of-sight μ, using its Legendre multipole moments.

Arguments

  • k: The wavenumber.
  • μ: The cosine of the angle to the line-of-sight.
  • Int_Mono: A function or interpolant that provides the monopole moment $I_0(k)$ at wavenumber k.
  • Int_Quad: A function or interpolant that provides the quadrupole moment $I_2(k)$ at wavenumber k.
  • Int_Hexa: A function or interpolant that provides the hexadecapole moment $I_4(k)$ at wavenumber k.

Returns

The value of the anisotropic power spectrum $P(k, \mu)$ at the given k and μ.

Details

The anisotropic power spectrum is reconstructed as a sum of its multipole moments multiplied by the corresponding Legendre polynomials evaluated at μ. The function uses the 0th, 2nd, and 4th order Legendre polynomials.

Formula

The formula used is:

\[P(k, \mu) = I_0(k) \mathcal{L}_0(\mu) + I_2(k) \mathcal{L}_2(\mu) + I_4(k) \mathcal{L}_4(\mu)\]

where $I_l(k)$ are the multipole moments and $\mathcal{L}_l(\mu)$ are the Legendre polynomials of order $l$.

See Also

  • _Legendre_0: Calculates the 0th order Legendre polynomial.
  • _Legendre_2: Calculates the 2nd order Legendre polynomial.
  • _Legendre_4: Calculates the 4th order Legendre polynomial.
source
Effort._k_trueMethod
_k_true(k_o, μ_o, q_perp, F)

Calculates the true (physical) wavenumber k from the observed wavenumber k_o and observed cosine of the angle to the line-of-sight μ_o.

This transformation accounts for anisotropic effects, likely redshift-space distortions (RSD) or anisotropic cosmological scaling, parameterized by q_perp and F.

Arguments

  • k_o: The observed wavenumber (scalar).
  • μ_o: The observed cosine of the angle to the line-of-sight (scalar).
  • q_perp: A parameter related to perpendicular anisotropic scaling.
  • F: A parameter related to parallel anisotropic scaling (often the growth rate f divided by the anisotropic scaling parameter q_parallel).

Returns

The calculated true wavenumber k (scalar).

Formula

The formula used is:

\[k = \frac{k_o}{q_\perp} \sqrt{1 + \mu_o^2 \left(\frac{1}{F^2} - 1\right)}\]

See Also

source
Effort._k_trueMethod
_k_true(k_o::Array, μ_o::Array, q_perp, F)

Calculates the true (physical) wavenumber k for arrays of observed wavenumbers k_o and observed cosines of the angle to the line-of-sight μ_o.

This method applies the transformation from observed to true wavenumber element-wise or for combinations of input arrays, accounting for anisotropic effects parameterized by q_perp and F.

Arguments

  • k_o: An array of observed wavenumbers.
  • μ_o: An array of observed cosines of the angle to the line-of-sight.
  • q_perp: A parameter related to perpendicular anisotropic scaling.
  • F: A parameter related to parallel anisotropic scaling.

Returns

A vector containing the calculated true wavenumbers k for the given input arrays.

Details

The function calculates k for pairs or combinations of values from the input arrays k_o and μ_o using a formula derived from anisotropic scaling. The calculation involves broadcasting and array operations to handle the array inputs efficiently. The result is flattened into a vector.

Formula

The underlying transformation for each pair of k_o and μ_o is:

\[k = \frac{k_o}{q_\perp} \sqrt{1 + \mu_o^2 \left(\frac{1}{F^2} - 1\right)}\]

See Also

source
Effort._μ_trueMethod
_μ_true(μ_o, F)

Calculates the true (physical) cosine of the angle to the line-of-sight μ from the observed cosine of the angle to the line-of-sight μ_o.

This transformation accounts for anisotropic effects, likely redshift-space distortions (RSD) or anisotropic cosmological scaling, parameterized by F.

Arguments

  • μ_o: The observed cosine of the angle to the line-of-sight (scalar).
  • F: A parameter related to parallel anisotropic scaling (often the growth rate f divided by the anisotropic scaling parameter q_parallel).

Returns

The calculated true cosine of the angle to the line-of-sight μ (scalar).

Formula

The formula used is:

\[\mu = \frac{\mu_o}{F \sqrt{1 + \mu_o^2 \left(\frac{1}{F^2} - 1\right)}}\]

See Also

source
Effort._μ_trueMethod
_μ_true(μ_o::Array, F)

Calculates the true (physical) cosine of the angle to the line-of-sight μ for an array of observed cosines of the angle to the line-of-sight μ_o.

This method applies the transformation from observed to true angle cosine element-wise, accounting for anisotropic effects parameterized by F.

Arguments

  • μ_o: An array of observed cosines of the angle to the line-of-sight.
  • F: A parameter related to parallel anisotropic scaling.

Returns

An array containing the calculated true cosines of the angle to the line-of-sight μ.

Details

The function calculates μ for each value in the input array μ_o using a formula derived from anisotropic scaling. Broadcasting (@.) is used to apply the calculation element-wise.

Formula

The underlying transformation for each μ_o is:

\[\mu = \frac{\mu_o}{F \sqrt{1 + \mu_o^2 \left(\frac{1}{F^2} - 1\right)}}\]

See Also

source
Effort._P_obsFunction
_P_obs(k_o, μ_o, q_par, q_perp, Int_Mono, Int_Quad, Int_Hexa)

Calculates the observed power spectrum $P_{\text{obs}}(k_o, \mu_o)$ at a given observed wavenumber k_o and observed cosine of the angle to the line-of-sight μ_o.

This function transforms the observed coordinates to true (physical) coordinates, calculates the true power spectrum using provided interpolants for the multipole moments, and applies the appropriate scaling factor due to anisotropic effects.

Arguments

  • k_o: The observed wavenumber.
  • μ_o: The observed cosine of the angle to the line-of-sight.
  • q_par: A parameter related to parallel anisotropic scaling.
  • q_perp: A parameter related to perpendicular anisotropic scaling.
  • Int_Mono: An interpolation function for the monopole moment $I_0(k)$ in true k.
  • Int_Quad: An interpolation function for the quadrupole moment $I_2(k)$ in true k.
  • Int_Hexa: An interpolation function for the hexadecapole moment $I_4(k)$ in true k.

Returns

The value of the observed power spectrum $P_{\text{obs}}(k_o, \mu_o)$.

Details

The observed coordinates $(k_o, \mu_o)$ are transformed to true coordinates $(k_t, \mu_t)$ using the _k_true and _μ_true functions, with $F = q_\parallel / q_\perp$. The true power spectrum $P(k_t, \mu_t)$ is then reconstructed using _Pkμ and the provided multipole interpolants. Finally, the result is scaled by $1 / (q_\parallel q_\perp^2)$.

Formula

The formula used is:

\[P_{\text{obs}}(k_o, \mu_o) = \frac{1}{q_\parallel q_\perp^2} P(k_t, \mu_t)\]

where

\[k_t = \text{_k_true}(k_o, \mu_o, q_\perp, F)\]

math \mut = \text{μtrue}(\muo, F)

and

math F = q\parallel / q\perp ```

See Also

  • _k_true: Transforms observed wavenumber to true wavenumber.
  • _μ_true: Transforms observed angle cosine to true angle cosine.
  • _Pkμ: Reconstructs the true power spectrum from multipole moments.
  • interp_Pℓs: Creates the multipole interpolants.
source
Effort.interp_PℓsFunction
interp_Pℓs(Mono_array, Quad_array, Hexa_array, k_grid)

Creates interpolation functions for the monopole, quadrupole, and hexadecapole moments of the power spectrum.

These interpolants can then be used to efficiently evaluate the multipole moments at arbitrary wavenumbers k.

Arguments

  • Mono_array: An array containing the values of the monopole moment $I_0(k)$.
  • Quad_array: An array containing the values of the quadrupole moment $I_2(k)$.
  • Hexa_array: An array containing the values of the hexadecapole moment $I_4(k)$.
  • k_grid: An array containing the corresponding wavenumber k values for the multipole arrays.

Returns

A tuple containing three interpolation functions: (Int_Mono, Int_Quad, Int_Hexa).

Details

The function uses AkimaInterpolation from the Interpolations.jl package to create the interpolants. Extrapolation is set to ExtrapolationType.Extension, which means the interpolant will use the nearest data points to extrapolate outside the provided k_grid range. Note that extrapolation can sometimes introduce errors.

See Also

  • _Pkμ: Uses the interpolation functions to reconstruct the anisotropic power spectrum.
source
Effort.apply_AP_checkMethod
apply_AP_check(k_input::Array, k_output::Array, Mono_array::Array, Quad_array::Array, Hexa_array::Array, q_par, q_perp)

Calculates the observed power spectrum multipole moments (monopole, quadrupole, hexadecapole) on a given observed wavenumber grid k_output, from arrays of true multipole moments provided on an input wavenumber grid k_input, using numerical integration.

This is a check version, intended for verifying results from faster methods. It is significantly slower due to the use of numerical integration over the angle μ.

Arguments

  • k_input: An array of wavenumber values on which the input true multipole moments (Mono_array, Quad_array, Hexa_array) are defined.
  • k_output: An array of observed wavenumber values at which to calculate the output observed multipoles.
  • Mono_array: An array containing the values of the true monopole moment $I_0(k)$ on the k_input grid.
  • Quad_array: An array containing the values of the true quadrupole moment $I_2(k)$ on the k_input grid.
  • Hexa_array: 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.

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 wavenumbers in k_output.

Details

This method first creates interpolation functions for the true multipole moments using interp_Pℓs based on the k_input grid. It then calls the core apply_AP_check(k_grid, int_Mono, int_Quad, int_Hexa, q_par, q_perp) method, passing k_output as the grid at which to calculate the observed multipoles.

This function is a slower check implementation and should not be used in performance-critical code.

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 observed power spectrum $P_{\text{obs}}(k_o, \mu_o)$ is calculated using _P_obs(k_o, μ_o, q_par, q_perp, int_Mono, int_Quad, int_Hexa).

See Also

source
Effort._Pk_reconFunction
_Pk_recon(mono::Matrix, quad::Matrix, hexa::Matrix, l0, l2, l4)

Reconstructs the anisotropic power spectrum $P(k, \mu)$ on a grid of wavenumbers k and cosines of the angle to the line-of-sight μ, using matrices of its Legendre multipole moments and vectors of Legendre polynomial values.

This function is designed to efficiently reconstruct the 2D power spectrum for multiple k and μ values simultaneously, assuming the multipole moments are provided as matrices (e.g., N_k x 1) and Legendre polynomials as vectors (e.g., N_μ).

Arguments

  • mono: A matrix containing the monopole moment $I_0(k)$ values (expected dimensions N_k x 1).
  • quad: A matrix containing the quadrupole moment $I_2(k)$ values (expected dimensions N_k x 1).
  • hexa: A matrix containing the hexadecapole moment $I_4(k)$ values (expected dimensions N_k x 1).
  • l0: A vector containing the 0th order Legendre polynomial $\mathcal{L}_0(\mu)$ values evaluated at the desired μ values (expected dimensions N_μ).
  • l2: A vector containing the 2nd order Legendre polynomial $\mathcal{L}_2(\mu)$ values evaluated at the desired μ values (expected dimensions N_μ).
  • l4: A vector containing the 4th order Legendre polynomial $\mathcal{L}_4(\mu)$ values evaluated at the desired μ values (expected dimensions N_μ).

Returns

A matrix representing the anisotropic power spectrum $P(k, \mu)$ on the N_k x N_μ grid.

Details

The function reconstructs the anisotropic power spectrum using the formula that sums the multipole moments multiplied by the corresponding Legendre polynomials. The matrix and vector operations are broadcast to calculate the result for all combinations of input k (from the rows of the moment matrices) and μ (from the elements of the Legendre polynomial vectors).

Formula

The formula used for each element $(i, j)$ of the output matrix (corresponding to the $i$-th wavenumber and $j$-th angle cosine) is:

\[P(k_i, \mu_j) = I_0(k_i) \mathcal{L}_0(\mu_j) + I_2(k_i) \mathcal{L}_2(\mu_j) + I_4(k_i) \mathcal{L}_4(\mu_j)\]

See Also

source

Utils

Effort._transformed_weightsFunction
_transformed_weights(quadrature_rule, order, a, b)

Transforms the points and weights of a standard quadrature rule from the interval [-1, 1] to a specified interval [a, b].

This is a utility function used to adapt standard quadrature rules (like Gauss-Legendre) for numerical integration over arbitrary intervals [a, b].

Arguments

  • quadrature_rule: A function that takes an order and returns a tuple (points, weights) for the standard interval [-1, 1].
  • order: The order of the quadrature rule (number of points).
  • a: The lower bound of the target interval.
  • b: The upper bound of the target interval.

Returns

A tuple (transformed_points, transformed_weights) for the interval [a, b].

Details

The transformation is applied to the standard points $x_i^{\text{std}}$ and weights $w_i^{\text{std}}$ obtained from the quadrature_rule:

  • Transformed points: $x_i = \frac{b - a}{2} x_i^{\text{std}} + \frac{b + a}{2}$
  • Transformed weights: $w_i = \frac{b - a}{2} w_i^{\text{std}}$

Formula

The transformation formulas are: Points: $x_i = \frac{b - a}{2} x_i^{\text{std}} + \frac{b + a}{2}$ Weights: $w_i = \frac{b - a}{2} w_i^{\text{std}}$

See Also

  • _r̃_z: An example function that uses this utility for numerical integration.
source
Effort._Legendre_0Function
_Legendre_0(x)

Calculates the 0th order Legendre polynomial, $\mathcal{L}_0(x)$.

Arguments

  • x: The input value (typically the cosine of an angle, -1 ≤ x ≤ 1).

Returns

The value of the 0th order Legendre polynomial evaluated at x.

Formula

The formula for the 0th order Legendre polynomial is:

\[\mathcal{L}_0(x) = 1\]

See Also

  • _Legendre_2: Calculates the 2nd order Legendre polynomial.
  • _Legendre_4: Calculates the 4th order Legendre polynomial.
  • _Pkμ: A function that uses Legendre polynomials.
source
Effort._Legendre_2Function
_Legendre_2(x)

Calculates the 2nd order Legendre polynomial, $\mathcal{L}_2(x)$.

Arguments

  • x: The input value (typically the cosine of an angle, -1 ≤ x ≤ 1).

Returns

The value of the 2nd order Legendre polynomial evaluated at x.

Formula

The formula for the 2nd order Legendre polynomial is:

\[\mathcal{L}_2(x) = \frac{1}{2} (3x^2 - 1)\]

See Also

  • _Legendre_0: Calculates the 0th order Legendre polynomial.
  • _Legendre_4: Calculates the 4th order Legendre polynomial.
  • _Pkμ: A function that uses Legendre polynomials.
source
Effort._Legendre_4Function
_Legendre_4(x)

Calculates the 4th order Legendre polynomial, $\mathcal{L}_4(x)$.

Arguments

  • x: The input value (typically the cosine of an angle, -1 ≤ x ≤ 1).

Returns

The value of the 4th order Legendre polynomial evaluated at x.

Formula

The formula for the 4th order Legendre polynomial is:

\[\mathcal{L}_4(x) = \frac{1}{8} (35x^4 - 30x^2 + 3)\]

See Also

  • _Legendre_0: Calculates the 0th order Legendre polynomial.
  • _Legendre_2: Calculates the 2nd order Legendre polynomial.
  • _Pkμ: A function that uses Legendre polynomials.
source
Effort._cubic_splineFunction
_cubic_spline(u, t, new_t::AbstractArray)

A convenience wrapper to create and apply a cubic spline interpolation using DataInterpolations.jl.

This function simplifies the process of creating a CubicSpline interpolant for the data (u, t) and evaluating it at the points new_t.

Arguments

  • u: An array of data values.
  • t: An array of data points corresponding to u.
  • new_t: An array of points at which to interpolate.

Returns

An array of interpolated values corresponding to new_t.

Details

This function is a convenience wrapper around DataInterpolations.CubicSpline(u, t; extrapolation=ExtrapolationType.Extension).(new_t). It creates a cubic spline interpolant with extrapolation enabled using ExtrapolationType.Extension and immediately evaluates it at all points in new_t.

See Also

  • DataInterpolations.CubicSpline: The underlying interpolation function.
  • _quadratic_spline: Wrapper for quadratic spline interpolation.
  • _akima_spline: Wrapper for Akima interpolation.
source
Effort._quadratic_splineFunction
_quadratic_spline(u, t, new_t::AbstractArray)

A convenience wrapper to create and apply a quadratic spline interpolation using DataInterpolations.jl.

This function simplifies the process of creating a QuadraticSpline interpolant for the data (u, t) and evaluating it at the points new_t.

Arguments

  • u: An array of data values.
  • t: An array of data points corresponding to u.
  • new_t: An array of points at which to interpolate.

Returns

An array of interpolated values corresponding to new_t.

Details

This function is a convenience wrapper around DataInterpolations.QuadraticSpline(u, t; extrapolation=ExtrapolationType.Extension).(new_t). It creates a quadratic spline interpolant with extrapolation enabled using ExtrapolationType.Extension and immediately evaluates it at all points in new_t.

See Also

  • DataInterpolations.QuadraticSpline: The underlying interpolation function.
  • _cubic_spline: Wrapper for cubic spline interpolation.
  • _akima_spline: Wrapper for Akima interpolation.
source
Effort._akima_splineFunction
_akima_spline(u, t, new_t::AbstractArray)

A convenience wrapper to create and apply an Akima interpolation using DataInterpolations.jl.

This function simplifies the process of creating an AkimaInterpolation interpolant for the data (u, t) and evaluating it at the points new_t.

Arguments

  • u: An array of data values.
  • t: An array of data points corresponding to u.
  • new_t: An array of points at which to interpolate.

Returns

An array of interpolated values corresponding to new_t.

Details

This function is a convenience wrapper around DataInterpolations.AkimaInterpolation(u, t; extrapolation=ExtrapolationType.Extension).(new_t). It creates an Akima interpolant with extrapolation enabled using ExtrapolationType.Extension and immediately evaluates it at all points in new_t.

See Also

  • DataInterpolations.AkimaInterpolation: The underlying interpolation function.
  • _cubic_spline: Wrapper for cubic spline interpolation.
  • _quadratic_spline: Wrapper for quadratic spline interpolation.
source