powspec_calc - All things power spectrum

powspechi.powspec_calc.alm_dNdphi(l, m, etacut=0.9, vns=array([1., 1., 1., 1.]), psis=array([0., 0., 0., 0.]), gsim=<function fconst>, *args, **kwargs)

Calculate the \(a_{\ell m}\) coefficients of a function of type \(f(\theta, \phi) = g(\theta) \cdot h(\phi)\), where:

\[h(\phi) = \frac{1}{2\pi} \left[ 1 + 2\sum_{n = 1}^{\infty} v_n \cos[n(\phi - \Psi_n)] \right].\]
Parameters
  • l (int, scalar) – The multipole moment \(\ell\) associated with the polar angle \(\theta\)

  • m (int, scalar) – The mode associated with the azimuthal angle \(\phi\)

  • etacut (float, scalar, optional) – The limit imposed on pseudorapidity, i.e., \(|\eta|\) < etacut. If there is no limit, just set it to None. Default: 0.9.

  • vns (float, optional) – The array representing \(v_n\), with \(n > 0\). Default: array([1., 1., 1., 1.]).

  • psis (float, optional) – The array representing \(\Psi_n\), with \(n > 0\). Default: array([0., 0., 0., 0.])

  • gsim (function, optional) – The polar function \(g(\theta)\). Default: monte_carlos.fconst.

  • *args – Arguments to be passed to gsim

  • **kwargs – Keyword-only arguments to be passed to gsim

Returns

a_lm – The coefficient with indices l and m from the decomposition of \(f(\theta, \phi)\) in spherical harmonics.

Return type

complex, scalar

Notes

It should be remarked that if the default values of vns and psis are used, one should get in return the values

\[b_{\ell m} \sim \int_{q_i}^{q_f} \sin{\theta} g(\theta) P_{\ell m}(\cos{\theta})d\theta,\]

where \((q_i, q_f)\) is the interval in \(\theta\) corresponding to the imposed \(\eta\) limit and \(P_{\ell m}\) are the associated Legendre polynomials.

powspechi.powspec_calc.avcls_zvtx(avcls, nevts)

Calculate the weighted average of the averaged spectra from distinct event ensembles.

Parameters
  • avcls (dict) – A dictionary following a certain hierarchy: ‘vtx’ \(\to\) ‘full’/’mdz’ \(\to\) list[mean_array, err_array]. The first key ‘vtx’ stands for the vertex interval to which the averaged spectrum belongs. The sub-dictionary associated with the key ‘vtx’ is the standard format for averaged spectrum found throughout this package.

  • nevts (dict) – A dictionary whose ‘vtx’ keys are the same as avcls, while its values correspond to the total number of events in each ensemble.

Returns

mean_zvtx – A dictionary with the same format as the standard averaged power spectrum, where keys ‘full’ and ‘mdz’ correspond to the full spectrum and \(C^{m\neq0}_{\ell}\), respectively.

Return type

dict

powspechi.powspec_calc.cls_calc(lsize, alms, *args, **kwargs)

Calculate the angular power spectrum analytically from a function of \(\ell, m\).

Parameters
  • lsize (int, scalar) – The maximum value for the multipole moment

  • alms (function) – The function which calculates the \(a_{\ell m}\) coefficients.

  • *args – Arguments to be passed to alms.

  • **kwargs – Keyword-only arguments to be passed to alms.

Returns

cls_true – A dictionary in the typical power spectrum format of the powspechi package.

Return type

dict

Notes

It is recommended to use cls_calc with alm_dNdphi as the alms parameter. One may create their own alms function to analytically calculate the angular power spectrum of functions not belonging to the type \(f(\theta, \phi) = g(\theta)\cdot h(\phi)\). However, keep in mind the execution time.

powspechi.powspec_calc.isobackground(clsres_file, skip=True)

From a special type of file create a dictionary containing \(\langle N_{\ell} \rangle\), i.e., an average power spectrum used to correct for the ensemble multiplicity.

Parameters
  • clsres_file (string) – A file containing the average power spectrum \(\langle N_{\ell} \rangle\). It has four columns which follow the order: ‘full’ ‘err_full’ ‘mdz’ ‘err_mdz’. Refer to maps2cld to see the meaning of ‘full’ and ‘mdz’. As for the prefix ‘err’, it indicates the error on the mean of its corresponding spectrum.

  • skip (bool, optional) – If True it skips the first line of the file, whereas if set to False no line will be skipped. Default: True.

Returns

clsres – A dictionary with keys ‘full’ and ‘mdz’, whose values are lists with the full spectrum and the same when \(m\neq0\). For each of these lists, the index 0 contains the mean, while index 1 contains the error on the mean. Both quantities are ndarrays.

Return type

dict

Notes

While the correction for the \(m\neq0\) average spectrum is simply \(\langle N^{m\neq0}_{\ell} \rangle\), in the case of the full spectrum, even and odd \(\ell\) modes are treated differently. Odd modes are corrected with odd \(\langle N_{\ell} \rangle\), whereas even modes are corrected with even \(\langle N^{m\neq0}_{\ell} \rangle\). The reason lies in considering the artificial spectrum features which arise from limited sky maps. If \(\langle C_{\ell} \rangle\) is simply subtracted by \(\langle N_{\ell} \rangle\), then such features will disappear, thus the resulting spectrum will not faithfully reproduce the expected full spectrum under said circumstances.

powspechi.powspec_calc.lns(nside)

Create a multipole (\(\ell\)) array based on the chosen resolution.

Parameters

nside (int, scalar) – A parameter related to the chosen HEALPix map resolution

Returns

ln – A 1-D array of int type that satisfies the chosen nside

Return type

int, ndarray

powspechi.powspec_calc.maps2cld(maps)

Calculate the angular power spectrum of a given map or maps.

Parameters

maps (array_like) – A single map or array/list of maps. It must be a HEALPix map, i.e., the number of indices must correspond to a nside value.

Returns

  • cld (dict) – A dictionary whose keys correspond to the ‘full’ power spectrum and the same without the \(a_{\ell 0}\) modes, denoted ‘mdz’. The values of cld are ndarrays with dimensions dependent on the number of entry maps and their resolution.

  • averd (dict) – If more than one map is given, the averaged power spectrum is calculated. Its keys are also ‘full’ and ‘mdz’. Its values are lists of arrays: index 0 corresponds to the mean cld value, while index 1 is the error on the mean.

Notes

A ‘full’ angular power spectrum has the following expression:

\[C_{\ell} = \frac{1}{2\ell + 1}\sum_{m = -\ell}^{m = \ell} |a_{\ell m}|^2,\]

while ‘mdz’, which stands for \(m\neq0\) has the form

\[C^{m\neq0}_{\ell} = C_{\ell} - \frac{1}{2\ell + 1} |a_{\ell 0}|^2,\]

\(a_{\ell m}\) are the coefficients associated with the spherical harmonics \(Y_{\ell m}\).

powspechi.powspec_calc.subisocorr(averd, isobkg)

Subtract the average spectrum calculated through HEALPix \(\langle C_{\ell} \rangle\) from the spectrum of ensemble multiplicity \(\langle N_{\ell}\rangle\).

Parameters
  • averd (dict) – A dictionary containing the power spectra \(\langle C_{\ell} \rangle\) and \(\langle C^{m\neq0}_{\ell} \rangle\). They should be contained in a list with index 0 for the mean and index 1 for its error. Such lists should be values corresponding to different keys. Their recommended names are ‘full’ and ‘mdz’, respectively.

  • isobkg (dict) – A dictionary following the same format, i.e., same keys and list types, as averd. It should contain the averaged spectrum used to correct for the ensemble’s multiplicity distribution, \(\langle N_{\ell} \rangle\).

Returns

averd_sic – A dictionary following the same format as averd. It contains the corrected averaged spectra \(\langle S_{\ell}\rangle\) and \(\langle S^{m\neq0}_{\ell}\rangle\), as well as their propagated error.

Return type

dict

powspechi.powspec_calc.vns_calc(n, averd, blms, mixed=True)

Calculate the \(v_n\) coefficients of a particle distribution \(f(\theta, \phi) = g(\theta) h(\phi)\) using the angular power spectrum method.

Parameters
  • n (int, scalar) – The index of \(v_n\). Should be larger than 0.

  • averd (dict) – The averaged angular power spectrum of a distribution of type \(f(\theta, \phi) = g(\theta) h(\phi)\). It follows the standard spectrum format.

  • blms (float, array_like) – The coefficients associated with the polar function \(g(\theta)\). See alm_dNdphi for its expression. The array indices should correspond to \(b_{nn}\).

  • mixed (bool, optional) – If True, the values of \(v_1,v_2\) are considered in the calculation of \(v_3,v_4\), respectively. In that case, the values \(b_{31}\) and \(b_{42}\) should be appended to the blms array. Default: True.

Returns

  • vn (float, scalar) – The resulting value of \(v_n\).

  • err (float, scalar) – The error associated with the calculation of vn.

Raises

IndexError – If one desires to calculate vn for \(n = 0\) or \(n > 4\).

Notes

The expressions for \(v_n\) are the following 1:

\[\begin{split}|v_n|^2 &= \frac{2n + 1}{2} \cdot \frac{C^{m\neq0}_n}{|b_{nn}|^2} \cdot \frac{|b_{00}|^2}{C_0}, \\ |v_n|^2 &= \frac{1}{|b_{nn}|^2} \left[ \frac{2n + 1}{2} \cdot C^{m\neq0}_n - \frac{2n - 3}{2} \cdot \frac{|b_{nn-2}|^2}{|b_{n-2n-2}|^2} \cdot C^{m\neq0}_{n-2} \right] \frac{|b_{00}|^2}{C_0},\end{split}\]

valid for \(n = 1, 2\) and \(n = 3, 4\), respectively.

References

1
  1. Machado, “Heavy ion anisotropies: a closer look at the angular power spectrum”, arXiv:1907.00413 [hep-ph] (2019).