pyccl.correlation module¶
Correlation functon computations.
- Choices of algorithms used to compute correlation functions:
- ‘Bessel’ is a direct integration using Bessel functions. ‘FFTLog’ is fast using a fast Fourier transform. ‘Legendre’ uses a sum over Legendre polynomials.
-
pyccl.correlation.
correlation
(cosmo, ell, C_ell, theta, corr_type='gg', method='fftlog')[source]¶ Compute the angular correlation function.
Parameters: - cosmo (
Cosmology
) – A Cosmology object. - ell (array_like) – Multipoles corresponding to the input angular power spectrum.
- C_ell (array_like) – Input angular power spectrum.
- theta (float or array_like) – Angular separation(s) at which to calculate the angular correlation function (in degrees).
- corr_type (string) – Type of correlation function. Choices: ‘GG’ (galaxy-galaxy), ‘GL’ (galaxy-shear), ‘L+’ (shear-shear, xi+), ‘L-‘ (shear-shear, xi-).
- method (string, optional) – Method to compute the correlation function. Choices: ‘Bessel’ (direct integration over Bessel function), ‘FFTLog’ (fast integration with FFTLog), ‘Legendre’ ( brute-force sum over Legendre polynomials).
Returns: - Value(s) of the correlation function at the input
angular separations.
Return type: float or array_like
- cosmo (
-
pyccl.correlation.
correlation_3d
(cosmo, a, r)[source]¶ Compute the 3D correlation function.
Parameters: Returns: Value(s) of the correlation function at the input distance(s).
-
pyccl.correlation.
correlation_3dRsd
(cosmo, a, s, mu, beta, use_spline=True)[source]¶ Compute the 3DRsd correlation function using linear approximation with multipoles.
Parameters: - cosmo (
Cosmology
) – A Cosmology object. - a (float) – scale factor.
- s (float or array_like) – distance(s) at which to calculate the 3DRsd correlation function (in Mpc).
- mu (float) – cosine of the angle at which to calculate the 3DRsd correlation function (in Radian).
- beta (float) – growth rate divided by galaxy bias.
- use_spline – switch that determines whether the RSD correlation function is calculated using global splines of multipoles.
Returns: Value(s) of the correlation function at the input distance(s) & angle.
- cosmo (
-
pyccl.correlation.
correlation_3dRsd_avgmu
(cosmo, a, s, beta)[source]¶ Compute the 3DRsd correlation function averaged over mu at constant s.
Parameters: Returns: Value(s) of the correlation function at the input distance(s) & angle.
-
pyccl.correlation.
correlation_multipole
(cosmo, a, beta, l, s)[source]¶ Compute the correlation multipoles.
Parameters: Returns: Value(s) of the correlation function at the input distance(s).