tpcf_multipole

halotools.mock_observables.tpcf_multipole(s_mu_tcpf_result, mu_bins, order=0)[source]

Calculate the multipoles of the two point correlation function after first computing s_mu_tpcf.

Parameters:
s_mu_tcpf_resultnp.ndarray

2-D array with the two point correlation function calculated in bins of \(s\) and \(\mu\). See s_mu_tpcf.

mu_binsarray_like

array of \(\mu = \cos(\theta_{\rm LOS})\) bins for which s_mu_tcpf_result has been calculated. Must be between [0,1].

orderint, optional

order of the multpole returned.

Returns:
xi_lnp.array

multipole of s_mu_tcpf_result of the indicated order.

Examples

For demonstration purposes we create a randomly distributed set of points within a periodic cube of length 250 Mpc/h.

>>> Npts = 100
>>> Lbox = 250.
>>> x = np.random.uniform(0, Lbox, Npts)
>>> y = np.random.uniform(0, Lbox, Npts)
>>> z = np.random.uniform(0, Lbox, Npts)

We transform our x, y, z points into the array shape used by the pair-counter by taking the transpose of the result of numpy.vstack. This boilerplate transformation is used throughout the mock_observables sub-package:

>>> sample1 = np.vstack((x,y,z)).T

First, we calculate the correlation function using s_mu_tpcf.

>>> from halotools.mock_observables import s_mu_tpcf
>>> s_bins  = np.linspace(0.01, 25, 10)
>>> mu_bins = np.linspace(0, 1, 15)
>>> xi_s_mu = s_mu_tpcf(sample1, s_bins, mu_bins, period=Lbox)

Then, we can calculate the quadrapole of the correlation function:

>>> xi_2 = tpcf_multipole(xi_s_mu, mu_bins, order=2)