pyccl.halos.profiles.ia module

class pyccl.halos.profiles.ia.SatelliteShearHOD(self, *, mass_def, concentration, a1h=0.001, b=-2, lmax=6, log10Mmin_0=12.0, log10Mmin_p=0.0, siglnM_0=0.4, siglnM_p=0.0, log10M0_0=7.0, log10M0_p=0.0, log10M1_0=13.3, log10M1_p=0.0, alpha_0=1.0, alpha_p=0.0, bg_0=1.0, bg_p=0.0, bmax_0=1.0, bmax_p=0.0, a_pivot=1.0, ns_independent=False, integration_method='FFTLog', rmin=0.001, N_r=512, N_jn=10000)[source]

Bases: HaloProfileHOD

Halo HOD class that calculates the satellite galaxy intrinsic shear field in real and fourier space, according to Fortuna et al. 2021.. Can be used to compute halo model-based intrinsic alignment (angular) power spectra.

The satellite intrinsic shear profile in real space is assumed to be

\[\gamma^I(r)=a_{1\mathrm{h}}\left(\frac{r}{r_\mathrm{vir}} \right)^b \sin^b\theta,\]

where \(a_{1\mathrm{h}}\) is the amplitude of intrinsic alignments on the 1-halo scale, \(b\) the index defining the radial dependence and \(\theta\) the angle defining the projection of the semi-major axis of the galaxy along the line of sight.

Parameters:
  • mass_def (MassDef or str) – a mass definition object, or a name string.

  • concentration (Concentration) – concentration-mass relation to use with this profile.

  • a1h (float) – Amplitude of the satellite intrinsic shear profile.

  • b (float) – Power-law index of the satellite intrinsic shear profile. If zero, the profile is assumed to be constant inside the halo.

  • lmax (int) – Maximum multipole to be summed in the plane-wave expansion (Eq. (C1) in Fortuna et al. 2021, default=6).

  • log10Mmin_0 (float) – offset parameter for \(\log_{10}M_{\rm min}\).

  • log10Mmin_p (float) – tilt parameter for \(\log_{10}M_{\rm min}\).

  • siglnM_0 (float) – offset parameter for \(\sigma_{{\rm ln}M}\).

  • siglnM_p (float) – tilt parameter for \(\sigma_{{\rm ln}M}\).

  • log10M0_0 (float) – offset parameter for \(\log_{10}M_0\).

  • log10M0_p (float) – tilt parameter for \(\log_{10}M_0\).

  • log10M1_0 (float) – offset parameter for \(\log_{10}M_1\).

  • log10M1_p (float) – tilt parameter for \(\log_{10}M_1\).

  • alpha_0 (float) – offset parameter for \(\alpha\).

  • alpha_p (float) – tilt parameter for \(\alpha\).

  • bg_0 (float) – offset parameter for \(\beta_g\).

  • bg_p (float) – tilt parameter for \(\beta_g\).

  • bmax_0 (float) – offset parameter for \(\beta_{\rm max}\).

  • bmax_p (float) – tilt parameter for \(\beta_{\rm max}\).

  • a_pivot (float) – pivot scale factor \(a_*\).

  • ns_independent (bool) – drop requirement to only form satellites when centrals are present.

  • integration_method (str) – Method used to obtain the fourier transform of the profile. Can be 'FFTLog', 'simpson' or 'spline'.

  • rmin (float) – For 'simpson' or 'spline' integration, minimum value of physical radius used to carry out the radial integral (in Mpc).

  • N_r (int) – For 'simpson' or 'spline' integration, number of points to be used when sampling the radial integral (in log space).

  • N_jn (int) – For 'simpson' or 'spline' integration, number of points to be used when sampling the spherical Bessel functions, that are later used to interpolate. Interpolating the Bessel functions increases the speed of the computations compared to explicitly evaluating them, without significant loss of accuracy.

update_parameters(*, a1h=None, b=None, lmax=None, log10Mmin_0=None, log10Mmin_p=None, siglnM_0=None, siglnM_p=None, log10M0_0=None, log10M0_p=None, log10M1_0=None, log10M1_p=None, alpha_0=None, alpha_p=None, bg_0=None, bg_p=None, bmax_0=None, bmax_p=None, a_pivot=None, ns_independent=None, rmin=None, N_r=None, N_jn=None)[source]

Update any of the parameters associated with this profile. Any parameter set to None won’t be updated.

Parameters:
  • a1h (float) – Amplitude of the satellite intrinsic shear profile.

  • b (float) – Power-law index of the satellite intrinsic shear profile. If zero, the profile is assumed to be constant inside the halo.

  • lmax (int) – Maximum multipole to be summed in the plane-wave expansion.

  • log10Mmin_0 (float) – offset parameter for \(\log_{10}M_{\rm min}\).

  • log10Mmin_p (float) – tilt parameter for \(\log_{10}M_{\rm min}\).

  • siglnM_0 (float) – offset parameter for \(\sigma_{{\rm ln}M}\).

  • siglnM_p (float) – tilt parameter for \(\sigma_{{\rm ln}M}\).

  • log10M0_0 (float) – offset parameter for \(\log_{10}M_0\).

  • log10M0_p (float) – tilt parameter for \(\log_{10}M_0\).

  • log10M1_0 (float) – offset parameter for \(\log_{10}M_1\).

  • log10M1_p (float) – tilt parameter for \(\log_{10}M_1\).

  • alpha_0 (float) – offset parameter for \(\alpha\).

  • alpha_p (float) – tilt parameter for \(\alpha\).

  • bg_0 (float) – offset parameter for \(\beta_g\).

  • bg_p (float) – tilt parameter for \(\beta_g\).

  • bmax_0 (float) – offset parameter for \(\beta_{\rm max}\).

  • bmax_p (float) – tilt parameter for \(\beta_{\rm max}\).

  • a_pivot (float) – pivot scale factor \(a_*\).

  • ns_independent (bool) – drop requirement to only form satellites when centrals are present

  • rmin (float) – For simpson or spline integration, minimum value of physical radius used to carry out the radial integral (in Mpc).

  • N_r (int) – For simpson or spline integration, number of points to be used when sampling the radial integral (in logarithmic space).

  • N_jn (int) – For simpson or spline integration, number of points to be used when sampling the spherical Bessel functions, that are later used to interpolate. Interpolating the Bessel functions increases the speed of the computations compared to explicitly evaluating them, without significant loss of accuracy.

get_normalization(cosmo, a, *, hmc)[source]

Returns the normalization of this profile, which is the mean galaxy number density.

Parameters:
Returns:

normalization factor of this profile.

Return type:

float

gamma_I(r, r_vir)[source]

Returns the intrinsic satellite shear,

\[\gamma^I(r)=a_{1\mathrm{h}} \left(\frac{r}{r_\mathrm{vir}}\right)^b.\]

If \(b\) is 0, then only the value of the amplitude \(a_\mathrm{1h}\) is returned. In addition, according to Fortuna et al. 2021, we use a constant value of 0.06 Mpc for \(r<0.06\) Mpc and set a maximum of 0.3 for \(\gamma^I(r)\).