abcmb.spectrum module

class abcmb.spectrum.SpectrumSolver(ellmin=2, ellmax=2500, lensing=False, k_axis_transfer=Array([1.00000000e-04, ..., 4.00000000e-01], dtype=float64), k_axis_Pk_output=Array([0.0001, ..., 0.1], dtype=float64), k_pivot=0.05, scale_sw=1, scale_isw=1, scale_dop=1, scale_pol=1)[source]

Bases: Module

CMB angular power spectrum computation.

Computes temperature and polarization angular power spectra by integrating transfer functions over wavenumber and time.

Attributes:

  • ells jnp.array

    Multipole values for output power spectra

  • ells_indices jnp.array

    Indices into bessel_l_tab corresponding to ells

  • lensing_ells jnp.array

    Extended multipole range for lensing calculations

  • lensing_ells_indices jnp.array

    Indices into bessel_l_tab for lensing multipoles

  • lensing_mus jnp.array

    Used for lensing, the Gauss-Legendre quadrature roots for the correlation function -> Cl integral.

  • lensing_ws jnp.array

    Used for lensing, the Gauss-Legendre quadrature weights for the correlation function -> Cl integral.

  • lensing bool

    Whether to include gravitational lensing effects

  • k_axis_transfer jnp.array

    Wavenumber grid for transfer function integration (units: Mpc^{-1})

  • k_axis_Pk_output jnp.array

    Wavenumber grid for matter power spectrum output (units: Mpc^{-1})

  • k_pivot float

    Pivot scale for primordial power spectrum normalization (units: Mpc^{-1}, default: 0.05)

  • scale_sw float

    Multiplicative factor for Sachs-Wolfe term (default: 1.0)

  • scale_isw float

    Multiplicative factor for integrated Sachs-Wolfe term (default: 1.0)

  • scale_dop float

    Multiplicative factor for Doppler term (default: 1.0)

  • scale_pol float

    Multiplicative factor for polarization term (default: 1.0)

Methods:

  • primordial_spectrum :

    Compute primordial power spectrum

  • Pk_lin :

    Compute linear matter power spectrum

  • get_Cl :

    Compute angular power spectra for multiple ℓ

  • Cl_one_ell :

    Compute angular power spectrum for single ℓ

  • integrand_T0 :

    Compute SW+ISW temperature source integrand

  • integrand_T1 :

    Compute ISW temperature source integrand

  • integrand_T2 :

    Compute polarization temperature source integrand

  • integrand_E :

    Compute E-mode polarization source integrand

Cl_one_ell(idx, PT, BG, params)[source]

Computes angular power spectrum for single multipole.

Integrates transfer functions over wavenumber.

Parameters:

  • idx int

    Index into bessel_l_tab for multipole ℓ

  • PT perturbations.PerturbationTable

    Perturbation evolution table

  • BG background.Background

    Background cosmology module

  • params dict

    Dictionary of input and derived parameters

Returns:

tuple

(C_ℓ^TT, C_ℓ^TE, C_ℓ^EE) angular power spectra

Pk_cb(k, z, PT, params)[source]

Compute linear Baryon+DarkMatter power spectrum at wavenumbers k and redshift z. Does not include any other massive species present.

Parameters:

  • k float or array

    Wavenumber (Mpc^{-1})

  • z float

    Redshift to evaluate.

  • PT perturbations.PerturbationTable

    Perturbation evolution table

  • params dict

    Dictionary of input and derived parameters

Returns:

float or array

Linear Baryon+DarkMatter power spectrum P_cb(k, z), units Mpc^3

Pk_lin(k, z, PT, params)[source]

Compute linear matter power spectrum at wavenumbers k and redshift z.

Parameters:

  • k float or array

    Wavenumber (Mpc^{-1})

  • z float

    Redshift to evaluate.

  • PT perturbations.PerturbationTable

    Perturbation evolution table

  • params dict

    Dictionary of input and derived parameters

Returns:

float or array

Linear matter power spectrum P(k, z), units Mpc^3

get_Cl(PT, BG, params)[source]

Compute angular power spectra for multiple multipoles.

Parameters:

  • PT perturbations.PerturbationTable

    Perturbation evolution table

  • BG background.Background

    Background cosmology module

  • params dict

    Dictionary of input and derived parameters

Returns:

tuple

(ClTT, ClTE, ClEE) angular power spectra

lensed_Cls(ells, ClTT_unlensed, ClTE_unlensed, ClEE_unlensed, PT, BG, params)[source]

Compute lensed CMB power spectra.

Applies gravitational lensing corrections to unlensed temperature and polarization power spectra using Wigner rotation matrices.

Parameters:

  • ells array

    Multipole values

  • ClTT_unlensed array

    Unlensed temperature power spectrum

  • ClTE_unlensed array

    Unlensed temperature-E-mode cross spectrum

  • ClEE_unlensed array

    Unlensed E-mode polarization power spectrum

  • PT perturbations.PerturbationTable

    Perturbation evolution table

  • BG background.Background

    Background cosmology module

  • params dict

    Dictionary of input and derived parameters

Returns:

tuple

(ClTT, ClTE, ClEE) lensed power spectra

lensing_Cl(ells, PT, BG, params)[source]

Angular lensing power spectrum at multipole ell.

IMPORTANT: Assumes Limber approximation throughout, even at ell=2.

Eq.(3.14) in astro-ph/0601594, except shifts ell -> ell+1/2 to match CLASS.

Parameters:

  • ell float or array

    Multipole

  • PT perturbations.PerturbationTable

    Perturbation evolution table

  • BG background.Background

    Background cosmology module

  • params dict

    Dictionary of input and derived parameters

Returns:

float or array

Angular lensing matter power spectrum Cl^phiphi, dimensionless.

lensing_power_spectrum(k, lna, PT, BG, params)[source]

Computes the lensing power spectrum at wavenumbers k and redshift z. Eq.(3.15) in astro-ph/0601594

Parameters:

  • k float or array

    Wavenumber (Mpc^{-1})

  • lna float

    Scale factor

  • PT perturbations.PerturbationTable

    Perturbation evolution table

  • BG background.Background

    Background cosmology module

  • params dict

    Dictionary of input and derived parameters

Returns:

float or array

Lensing matter power spectrum P(k, z), dimensionless.

primordial_spectrum(k, params)[source]

Compute primordial curvature power spectrum.

Parameters:

  • k float or array

    Wavenumber (units: Mpc^{-1})

  • params dict

    Dictionary of input and derived parameters

Returns:

float or array

Primordial power spectrum P_R(k), units Mpc^3

abcmb.spectrum.phi0(i, x)[source]

New method for computing phi0 (or just jl) We tabulated the bessel function between its smallest value (~1.e-10) out to the fifth local maximum. This is a different interval for each l, but we kept identical shape so it can be a large 2D array. If the incoming argument is within this interval, we use fast_interp. Otherwise we use the large x expansion above.

abcmb.spectrum.phi1(i, x)[source]

New method for computing phi1, or jl’. We tabulated the bessel function between its smallest value (~1.e-10) out to the fifth local maximum. This is a different interval for each l, but we kept identical shape so it can be a large 2D array. If the incoming argument is within this interval, we use fast_interp. Otherwise we use the large x expansion above.

abcmb.spectrum.phi2(i, x)[source]

New method for computing phi2 = (3 jl’’ + jl)/2 We tabulated the bessel function between its smallest value (~1.e-10) out to the fifth local maximum. This is a different interval for each l, but we kept identical shape so it can be a large 2D array. If the incoming argument is within this interval, we use fast_interp. Otherwise we use the large x expansion above.