1. Reference

The unequalpy package contains the following modules:

1.1. Matter power spectrum

This module computes the equal- and unequal-time matter power spectrum for different formalisms: Standard Perturbation Theory and Effective Field Theory.


unequalpy.matter.matter_power_spectrum_1loop(wavenumber, growth, powerk, counterterm=0, model='spt')[source]

One-loop matter power spectrum. This function computes the one-loop matter power spectrum in real space in standard perturbation theory at a single redshift, as described in [1]_ and [2]_.

Parameters
  • wavenumber ((nk,) array_like) – Array of wavenumbers in units of \([{\rm Mpc}^{-1}]\) at which to evaluate the matter power spectrum.

  • growth ((nz,) array_like) – Array of linear growth function at one single redshift.

  • powerk (tuple, function) – Tuple of functions for the linear, 22-type and 13-type power spectra at redshift zero.

  • counterterm ((nz,) array_like) – Array of counterterms dealing with deviations from perfect fluid, as described in equation 2.55 in [2], in units of \([{\rm Mpc}^{2}]\). Default is 0.

  • model (string, optional) – You can choose from two perturbation frameworks: {‘spt’: standard perturbation theory} described in [1] or {‘eft’: effective field theory} described in [2]. Default is ‘spt’.

Returns

power_spectrum – The matter power spectrum evaluated at the input redshifts and wavenumbers for the given cosmology. Units of \([{\rm Mpc}^{3}]\).

Return type

(nz,nk) array_like

Examples

>>> import numpy as np
>>> from astropy.cosmology import FlatLambdaCDM
>>> from skypy.power_spectrum import growth_function
>>> cosmo = FlatLambdaCDM(H0=67.11, Ob0=0.049, Om0= 0.2685)

We use precomputed values from the FAST-PT code:

>>> d = np.loadtxt('../Pfastpt.txt',unpack=True)
>>> ks = d[:, 0]
>>> pk, p22, p13 = d[:, 1], d[:, 2], d[:, 3]
>>> p11_int = interp1d( ks, pk, fill_value="extrapolate")
>>> p22_int = interp1d( ks, p22, fill_value="extrapolate")
>>> p13_int = interp1d( ks, p13, fill_value="extrapolate")
>>> powerk = (p11_int, p22_int, p13_int)

The normalised growth function from SkyPy:

>>> g0 = growth_function(0, cosmo)
>>> D0 = growth_function(0, cosmo) / g0

The best-fit counterterm using the Quijote simulations:

>>> ct2 = -0.4

And finally, the SPT and EFT equal-time matter power spectra:

>>> pspt = matter_power_spectrum_1loop(k, D0, powerk)
>>> peft = matter_power_spectrum_1loop(k, D0, ct2, model='eft')

References

..[1] Jeong, D. and Komatsu, E., 2006, astro-ph/0604075. ..[2] de la Bella, L. et al., 2017, doi:10.1088/1475-7516/2017/11/039.

unequalpy.matter.matter_unequal_time_power_spectrum(wavenumber, growth1, growth2, powerk, counterterm1=0, counterterm2=0, model='spt')[source]

Unequal-time one-loop matter power spectrum. This function computes the unequal-time one-loop matter power spectrum in real space in standard perturbation theory, as described in [1]_.

Parameters
  • wavenumber ((nk,) array_like) – Array of wavenumbers in units of \([{\rm Mpc}^{-1}]\) at which to evaluate the matter power spectrum.

  • growth1 ((nz1,), (nz2,) array_like) – Arrays of linear growth function at two redshifts.

  • growth2 ((nz1,), (nz2,) array_like) – Arrays of linear growth function at two redshifts.

  • powerk (tuple, function) – Tuple of functions for the linear, 22-type and 13-type power spectra at redshift zero.

  • counterterm1 ((nz1,), (nz2,) array_like) – Array of counterterms dealing with deviations from perfect fluid, as described in equation 2.55 in [1], in units of \([{\rm Mpc}^{2}]\). Default is 0.

  • counterterm2 ((nz1,), (nz2,) array_like) – Array of counterterms dealing with deviations from perfect fluid, as described in equation 2.55 in [1], in units of \([{\rm Mpc}^{2}]\). Default is 0.

  • model (string, optional) – You can choose from two perturbation frameworks: {‘spt’: standard perturbation theory} described in [2] or {‘eft’: effective field theory} described in [2]. Default is ‘spt’.

Returns

power_spectrum – The unequal-time matter power spectrum evaluated at the input redshifts and wavenumbers for the given cosmology, in units of \([{\rm Mpc}^{3}]\).

Return type

(nz1, nz2, nk) array_like

Examples

>>> import numpy as np
>>> from astropy.cosmology import FlatLambdaCDM
>>> from skypy.power_spectrum import growth_function
>>> cosmo = FlatLambdaCDM(H0=67.11, Ob0=0.049, Om0= 0.2685)

We use precomputed values from the FAST-PT code:

>>> d = np.loadtxt('../Pfastpt.txt',unpack=True)
>>> ks = d[:, 0]
>>> pk, p22, p13 = d[:, 1], d[:, 2], d[:, 3]
>>> p11_int = interp1d( ks, pk, fill_value="extrapolate")
>>> p22_int = interp1d( ks, p22, fill_value="extrapolate")
>>> p13_int = interp1d( ks, p13, fill_value="extrapolate")
>>> powerk = (p11_int, p22_int, p13_int)

The normalised growth function from SkyPy:

>>> g0 = growth_function(0, cosmo)
>>> D0 = growth_function(0, cosmo)]) / g0
>>> D1 = growth_function(1, cosmo)]) / g0

The best-fit counterterm using the Quijote simulations:

>>> ct0, ct1 = -0.4, -0.2

And finally, the SPT and EFT unequal-time matter power spectra:

>>> pspt = matter_unequal_time_power_spectrum(0.1, D0, D1, powerk)
>>> peft = matter_unequal_time_power_spectrum(0.1, D0, D1, powerk,
...     ct0, ct1, 'eft')

References

..[1] de la Bella, L. and Tessore, N. and Bridle, S., 2020,

arXiv 2011.06185.

1.2. Approximation to the power spectrum

This module prepares the power spectra for the lensing analysis, using different approximations to the unequal-time power spectrum.


unequalpy.approximation.geometric_approx(power_spectrum)[source]

Geometric approximation for the power spectrum. This function computes the unequal-time geometric mean approximation for any power spectrum, as described in equation 1.3 in [1]_.

Parameters

power_spectrum ((nz, k) array_like) – Array of power spectra at different redshifts.

Returns

power_power – The geometric approximation of the unequal-time power spectrum evaluated at the input redshifts and wavenumbers for the given cosmology. Units of \({\rm Mpc}^{3}\).

Return type

(nz, nz, nk) array_like

Examples

>>> import numpy as np
>>> from astropy.cosmology import FlatLambdaCDM
>>> from skypy.power_spectrum import growth_function
>>> from unequalpy.approximation import geometric_approx as Pgeom
>>> cosmo = FlatLambdaCDM(H0=67.11, Ob0=0.049, Om0= 0.2685)

We use precomputed values from the FAST-PT code:

>>> d = np.loadtxt('../Pfastpt.txt',unpack=True)
>>> ks = d[:, 0]
>>> pk, p22, p13 = d[:, 1], d[:, 2], d[:, 3]
>>> p11_int = interp1d( ks, pk, fill_value="extrapolate")
>>> p22_int = interp1d( ks, p22, fill_value="extrapolate")
>>> p13_int = interp1d( ks, p13, fill_value="extrapolate")
>>> powerk = (p11_int, p22_int, p13_int)

The normalised growth function from SkyPy:

>>> g0 = growth_function(0, cosmo)
>>> z = np.array([0,1])
>>> Dz = growth_function(z, cosmo)]) / g0

The equal-time matter power spectrum

>>> pet = P1loop(ks, Dz, powerk)

And finally, the geometric approximation to the unequal-time matter power spectrum:

>>> pg_spt = Pgeom(pet)

References

..[1] de la Bella, L. and Tessore, N. and Bridle, S., 2020,

arXiv 2011.06185.

unequalpy.approximation.midpoint_approx(wavenumber, growth_mean, powerk)[source]

Midpoint approximation for the power spectrum. This function computes the unequal-time midpoint approximation for any power spectrum, as described in equation 2.14 in [1]_.

Parameters
  • wavenumber ((nk,) array_like) – Array of wavenumbers in units of \({\rm Mpc}^{-1}\) at which to evaluate the matter power spectrum.

  • growth_mean ((nz1, nz2), array_like) – Growth function array evaluated at the midpoint redshift.

  • powerk (tuple, function) – Tuple of functions for the linear, 22-type and 13-type power spectra at redshift zero.

Returns

power – The midpoint power spectrum evaluated at the input redshifts and wavenumbers for the given cosmology. Units of \({\rm Mpc}^{3}\).

Return type

(nz1, nz2, nk), array_like

Examples

>>> import numpy as np
>>> from astropy.cosmology import FlatLambdaCDM
>>> from skypy.power_spectrum import growth_function
>>> from unequalpy.approximation import midpoint_approx
>>> from unequalpy.approximation import growth_midpoint
>>> cosmo = FlatLambdaCDM(H0=67.11, Ob0=0.049, Om0= 0.2685)

We use precomputed values from the FAST-PT code:

>>> d = np.loadtxt('../Pfastpt.txt',unpack=True)
>>> ks = d[:, 0]
>>> pk, p22, p13 = d[:, 1], d[:, 2], d[:, 3]
>>> p11_int = interp1d( ks, pk, fill_value="extrapolate")
>>> p22_int = interp1d( ks, p22, fill_value="extrapolate")
>>> p13_int = interp1d( ks, p13, fill_value="extrapolate")
>>> powerk = (p11_int, p22_int, p13_int)

The normalised growth function from SkyPy at some midoint redshift between z1 and z2:

>>> z1, z2 = 0, 2
>>> D12  = growth_midpoint(z1, z2, growth_function, cosmo)

Finally, the midpoint approximation to the unequal-time matter power spectrum:

>>> pspt = midpoint_approx(ks, D12, powerk)

References

..[1] de la Bella, L. and Tessore, N. and Bridle, S., 2020,

arXiv 2011.06185.

unequalpy.approximation.growth_midpoint(redshift1, redshift2, growth_function, *args)[source]

Growth function at the midpoint. This function computes the linear growth function evaluated at the mean of two given redshift values.

Parameters
  • redshift1 ((nz1,), (nz2,), array_like) – Array of wavenumbers in units of \({\rm Mpc}^{-1}\) at which to evaluate the matter power spectrum.

  • redshift2 ((nz1,), (nz2,), array_like) – Array of wavenumbers in units of \({\rm Mpc}^{-1}\) at which to evaluate the matter power spectrum.

  • growth_function (function) – Method to compute the growth function.

  • *args (tuple) – Arguments for the growth function method.

Returns

growth – The growth function evaluated at the midpoint redshift.

Return type

(nz1, nz2), array_like

Examples

>>> import numpy as np
>>> from astropy.cosmology import FlatLambdaCDM
>>> from skypy.power_spectrum import growth_function
>>> from unequalpy.approximation import growth_midpoint
>>> cosmo = FlatLambdaCDM(H0=67.11, Ob0=0.049, Om0= 0.2685)

The normalised growth function from SkyPy (you can choose any other method) at some midoint redshift between z1 and z2:

>>> z1, z2 = 0, 2
>>> D12  = growth_midpoint(z1, z2, growth_function, cosmo)

1.3. Lens filters

This module computes different filters for the computation of the angular power spectrum for the lensing potential


unequalpy.lens_filter.lensing_efficiency(x, zx, nz)[source]

Lensing efficiency. This function computes the general form of the lensing efficiency function given in [1]_.

Parameters
  • x ((nx,) array_like) – Array of comoving distances at which evaluate the lensing efficiency function.

  • zx (array_like) – Array of redshift as a function of comoving distance.

  • nz ((nx, ) array_like) – Array of redshift distribution of galaxies.

Returns

efficiency – Array of lensing efficiency values.

Return type

(nx,) array_like

Examples

>>> import numpy as np
>>> from astropy.cosmology import FlatLambdaCDM
>>> cosmo = FlatLambdaCDM(H0=67.11, Ob0=0.049, Om0= 0.2685)
>>> z_list = np.linspace(0, 1100,num=1000000)
>>> xl = cosmo.comoving_distance(z_list).value
>>> qf = lensing_efficiency(xl, 1)

References

..[1] Lemos, P. and Challinor, A. and Efstathiou, G., 2017,

arXiv: 1704.01054.

unequalpy.lens_filter.lensing_efficiency_cmb(x, xs)[source]

Parametric lensing efficiency cmb. This function computes the cmb lensing efficiency function given in [1]_.

Parameters
  • x ((nx,) array_like) – Array of comoving distances at which evaluate the lensing efficiency function.

  • xs (float) – Value of the comoving distance at the last scatering surface.

Returns

efficiency – Array of lensing efficiency values.

Return type

(nx,) array_like

Examples

>>> import numpy as np
>>> from astropy.cosmology import FlatLambdaCDM
>>> cosmo = FlatLambdaCDM(H0=67.11, Ob0=0.049, Om0= 0.2685)
>>> z_list = np.linspace(0, 1100,num=1000000)
>>> xs = 14000
>>> xl = cosmo.comoving_distance(z_list).value
>>> qf = lensing_efficiency_cmb(xl, xs)

References

..[1] Hanson, D. and Challinor, A. and Lewis, A., 2010,

doi: 10.1007/s10714-010-1036-y.

unequalpy.lens_filter.filter_shear(x, zx, lens_efficiency, cosmology)[source]

Filter for shear. This function computes the filter for the shear power spectra, described in [1]_.

Parameters
  • x ((nx,) array_like) – Array of comoving distances at which evaluate the lensing efficiency function.

  • zx (array_like) – Array of redshift as a function of comoving distance.

  • lens_efficiency ((nx,) array_like) – Array of lensing efficiency values.

  • cosmology (astropy.cosmology.Cosmology) – Cosmology object providing methods for the evolution history of omega_matter and omega_lambda with redshift.

Returns

filter – Array of filter values.

Return type

(nx,) array_like

Examples

>>> import numpy as np
>>> from astropy.cosmology import FlatLambdaCDM
>>> cosmo = FlatLambdaCDM(H0=67.11, Ob0=0.049, Om0= 0.2685)
>>> z_list = np.linspace(0, 1100,num=1000000)
>>> xl = cosmo.comoving_distance(z_list).value
>>> q1 = lensing_efficiency(xl, 1)
>>> f1 = filter(xl, z_list, q1, cosmology)

References

..[1] Lemos, P. and Challinor, A. and Efstathiou, G., 2017,

arXiv: 1704.01054.

unequalpy.lens_filter.filter_convergence(x, zx, lens_efficiency, cosmology)[source]

Filter for convergence. This function computes the filter for the convergenge power spectra, described in [1]_.

Parameters
  • x ((nx,) array_like) – Array of comoving distances at which evaluate the lensing efficiency function.

  • zx (array_like) – Array of redshift as a function of comoving distance.

  • lens_efficiency ((nx,) array_like) – Array of lensing efficiency values.

  • cosmology (astropy.cosmology.Cosmology) – Cosmology object providing methods for the evolution history of omega_matter and omega_lambda with redshift.

Returns

filter – Array of filter values.

Return type

(nx,) array_like

Examples

>>> import numpy as np
>>> from astropy.cosmology import FlatLambdaCDM
>>> cosmo = FlatLambdaCDM(H0=67.11, Ob0=0.049, Om0= 0.2685)
>>> z_list = np.linspace(0, 1100,num=1000000)
>>> xl = cosmo.comoving_distance(z_list).value
>>> q1 = lensing_efficiency(xl, 1)
>>> f1 = filter(xl, z_list, q1, cosmology)

References

..[1] Kilbinger, M., 2014, doi: 10.1088/0034-4885/78/8/086901.

unequalpy.lens_filter.filter_galaxy_clustering(x, zx, nz, linear_bias, cosmology)[source]

Filter for galaxy clustering. This function computes filter for the galaxy clustering power spectra, described in [1]_.

Parameters
  • x ((nx,) array_like) – Array of comoving distances at which evaluate the lensing efficiency function.

  • zx (array_like) – Array of redshift as a function of comoving distance.

  • nz ((nx,) array_like) – Array of redshift distribution of galaxies.

  • linear_bias (float) – Float for the value of the linear bias parameter.

  • cosmology (astropy.cosmology.Cosmology) – Cosmology object providing methods for the evolution history of omega_matter and omega_lambda with redshift.

Returns

filter – Array of filter values.

Return type

(nx,) array_like

Examples

>>> import numpy as np
>>> from astropy.cosmology import FlatLambdaCDM
>>> cosmo = FlatLambdaCDM(H0=67.11, Ob0=0.049, Om0= 0.2685)
>>> z_list = np.linspace(0, 1100,num=1000000)
>>> xl = cosmo.comoving_distance(z_list).value
>>> n1 = redshift_distribution_galaxies(xl, 1)
>>> f1 = filter(xl, z_list, n1, cosmo)

References

..[1] Kilbinger, M., 2014, doi: 10.1088/0034-4885/78/8/086901.