CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
cbl::cosmology::HaloProfile Class Reference

The class HaloProfile. More...

#include "Headers/HaloProfile.h"

Public Member Functions

Constructors/destructors
 HaloProfile ()=default
 default constructor
 
 HaloProfile (const cbl::cosmology::Cosmology cosmology, const double redshift, const double conc, const double Mass, const double Delta, const std::string profile_author, const std::string halo_def, const double trunc_fact=3., const bool miscentering=false, const bool single_profile=false, const double sigma_off=0.1, const double f_off=1.)
 constructor that sets the halo profile parameters. If the miscentering is considered, \(\Sigma(R)\) and \(\Delta\Sigma(R)\) are expressed as (see e.g. Bellagamba et al. 2019): More...
 
 HaloProfile (const cbl::cosmology::Cosmology cosmology, const double redshift, const std::string cM_author, const double Mass, const double Delta, const std::string profile_author, const std::string halo_def, const double trunc_fact=3., const bool miscentering=false, const bool single_profile=false, const double sigma_off=0.1, const double f_off=1.)
 constructor that sets the halo profile parameters, assuming a concentration-mass relation. If the miscentering is considered, \(\Sigma(R)\) and \(\Delta\Sigma(R)\) are expressed as (see e.g. Bellagamba et al. 2019): More...
 
 ~HaloProfile ()=default
 default destructor
 
Member functions used to get the private members or

to compute the cluster halo profiles

double concentration ()
 get the concentration More...
 
double concentration2 (const double Vmax, const double Rmax) const
 compute the halo concentration More...
 
double Delta ()
 get the overdensity More...
 
std::vector< double > density_profile_3D (const std::vector< double > rad)
 the 3D halo density profile. More...
 
double rho_s ()
 the characteristic density of the 3D halo density profile. More...
 
std::vector< double > Sigma (const std::vector< double > rad)
 the total surface density profile of the halo. More...
 
std::vector< double > Sigma_mis (const std::vector< double > rad)
 the miscentering contribution to the total surface density profile of the halo. More...
 
std::vector< double > Sigma_cen (const std::vector< double > rad)
 the centered contribution to the total surface density profile of the halo. More...
 
std::vector< double > Sigma_2h (const std::vector< double > rad, const std::string bias_author, const std::string method_Pk="EisensteinHu", const std::string interp_type="Linear", const bool NL=false)
 the 2-halo contribution to the total surface density profile of the halo, computed by assuming a halo bias model. More...
 
std::vector< double > Sigma_2h (const std::vector< double > rad, const double bias, const std::string method_Pk="EisensteinHu", const std::string interp_type="Linear", const bool NL=false)
 the 2-halo contribution to the total surface density profile of the halo. More...
 
std::vector< double > DeltaSigma (const std::vector< double > rad)
 the total excess surface density profile of the halo. More...
 
std::vector< double > DeltaSigma_mis (const std::vector< double > rad)
 the miscentering contribution to the excess surface density profile of the halo. More...
 
std::vector< double > DeltaSigma_cen (const std::vector< double > rad)
 the centered contribution to the excess surface density profile of the halo. More...
 
std::vector< double > DeltaSigma_2h (const std::vector< double > rad, const std::string bias_author, const std::string method_Pk="EisensteinHu", const std::string interp_type="Linear", const bool NL=false)
 the 2-halo contribution to the excess surface density profile of the halo, computed by assuming a halo bias model. More...
 
std::vector< double > DeltaSigma_2h (const std::vector< double > rad, const double bias, const std::string method_Pk="EisensteinHu", const std::string interp_type="Linear", const bool NL=false)
 the 2-halo contribution to the excess surface density profile of the halo. More...
 
double density_profile_FourierSpace (const double kk)
 the Fourier transform of the normalised halo density More...
 
double Mass_Delta (const double Mass, const double Delta_in, const double Delta_out, const double conc, const bool is_input_conc, const double rRmin_guess=1.e-3, const double rRmax_guess=10.) const
 The halo mass converted to a different value of \(\Delta\), assuming the Navarro-Frenk-White density profile. More...
 
Member functions used to set the private members
void set_cosmology (const cbl::cosmology::Cosmology cosmology)
 set the cosmological model More...
 
void set_mass (const double mass=par::defaultDouble)
 set the private member m_mass More...
 
void set_concentration (const double conc=par::defaultDouble)
 set the private member m_concentration More...
 
void set_f_off (const double f_off=par::defaultDouble)
 set the private member m_f_off More...
 
void set_sigma_off (const double sigma_off=par::defaultDouble)
 set the private member m_sigma_off More...
 
void set_trunc_fact (const double trunc_fact=par::defaultDouble)
 set the private member m_trunc_fact More...
 

Private Member Functions

Protected member functions
void m_set_profile (const cbl::cosmology::Cosmology cosmology, const double redshift, const double conc, const double Mass, const double Delta, const std::string profile_author, const std::string halo_def, const double trunc_fact=3., const bool miscentering=false, const bool single_profile=false, const double sigma_off=0.1, const double f_off=1.)
 private function that sets the halo profile parameters. If the miscentering is considered, \(\Sigma(R)\) and \(\Delta\Sigma(R)\) are expressed as (see e.g. Bellagamba et al. 2019): More...
 
double m_return_set_concentration ()
 Return the concentration set by the user. More...
 
double m_concentration_Duffy ()
 The concentration-mass relation by Duffy et al. (2008): More...
 
double m_rho_s_NFW ()
 Characteristic density of the NFW profile. More...
 
double m_rho_s_NFW_trunc ()
 Characteristic density of the truncated NFW profile. More...
 
std::vector< double > m_rho_NFW (const std::vector< double > rad)
 The NFW 3D density profile. More...
 
std::vector< double > m_rho_NFW_trunc (const std::vector< double > rad)
 The truncated NFW 3D density profile. More...
 
std::vector< double > m_Sigma_NFW (const std::vector< double > rad)
 Surface density of the NFW profile. More...
 
std::vector< double > m_Sigma_mean_NFW (const std::vector< double > rad)
 Mean surface density of the NFW profile. More...
 
std::vector< double > m_Sigma_NFW_trunc (const std::vector< double > rad)
 Surface density of the truncated NFW profile. More...
 
std::vector< double > m_Sigma_mean_NFW_trunc (const std::vector< double > rad)
 Mean surface density of the truncated NFW profile. More...
 
std::vector< double > m_Sigma_mis (const std::vector< double > rad)
 the mis-centered surface density profile of the halo (see e.g. Bellagamba et al. 2019). More...
 
std::vector< double > m_Sigma_including_miscentering (const std::vector< double > rad)
 Surface density including the miscentering contribution. More...
 
std::vector< double > m_DeltaSigma_cen (const std::vector< double > rad)
 Excess surface density without the miscentering contribution. More...
 
std::vector< double > m_DeltaSigma_mis (const std::vector< double > rad)
 the mis-centered excess surface density profile of the halo (see e.g. Bellagamba et al. 2019). More...
 
std::vector< double > m_DeltaSigma_including_miscentering (const std::vector< double > rad)
 Excess surface density including the miscentering contribution. More...
 
void m_set_cM_relation (const std::string cM_author="Duffy")
 set the concentration-mass relation. A function is set, computing the concentration of a dark matter halo of a given a mass, at a given redshift; the models implemented are the following: More...
 

Private Attributes

std::shared_ptr< cosmology::Cosmologym_cosmology = NULL
 pointer to the input cosmology
 
std::string m_profile_author = par::defaultString
 cluster density profile author
 
std::string m_halo_def
 halo (overdensity) definition
 
double m_Delta
 overdensity factor
 
std::function< double(const double, cbl::cosmology::Cosmology, const double)> m_Delta_func
 function returning the overdensity factor
 
double m_trunc_fact
 truncation factor $F_t$ defining the truncation radius, that is \(r_t = F_t*r_{\Delta}\)
 
bool m_miscentering
 if true, account for the miscentering
 
double m_sigma_off
 rms of the off-centered cluster distribution
 
double m_f_off
 fraction of off-centered clusters
 
double m_sigma_off_threshold
 lower threshold for sigma_off
 
double m_f_off_threshold
 lower threshold for f_off
 
double m_min_Smis
 minimum value where the off-centered Sigma(R) is computed
 
double m_max_Smis
 maxmimum value where the off-centered Sigma(R) is computed
 
bool m_single_profile
 if true, the miscetering is related to a single cluster. Otherwise, a population of clusters is considered
 
std::vector< double >(HaloProfile::* m_rho_ptr )(const std::vector< double >)
 cluster 3D density profile function
 
double(HaloProfile::* m_rho_s_ptr )()
 function returning the characteristic density of the cluster profile
 
std::vector< double >(HaloProfile::* m_Sigma_mis_ptr )(const std::vector< double >)
 function returning the miscentering contribution to the total cluster surface density for a given radius
 
std::vector< double >(HaloProfile::* m_Sigma_cen_ptr )(const std::vector< double >)
 function returning the centered cluster surface density for a given radius
 
std::vector< double >(HaloProfile::* m_Sigma_mean_cen_ptr )(const std::vector< double >)
 function returning the centered mean cluster surface density for a given radius
 
std::vector< double >(HaloProfile::* m_Sigma_ptr )(const std::vector< double >)
 function returning the cluster surface density for a given radius
 
std::vector< double >(HaloProfile::* m_DeltaSigma_mis_ptr )(const std::vector< double >)
 function returning the miscentering contribution to the total cluster excess surface density for a given radius
 
std::vector< double >(HaloProfile::* m_DeltaSigma_cen_ptr )(const std::vector< double >)
 function returning the centered cluster excess surface density for a given radius
 
std::vector< double >(HaloProfile::* m_DeltaSigma_ptr )(const std::vector< double >)
 function returning the excess cluster surface density for a given radius
 
bool m_isSet_cM_relation = false
 if true, the concentration-mass relation is set
 
double(HaloProfile::* m_return_concentration )()
 pointer returning the concentration, the one set by the user or from a concentration-mass relation function
 
double m_A = par::defaultDouble
 parameter A in the c-M relation
 
double m_B = par::defaultDouble
 parameter B in the c-M relation
 
double m_C = par::defaultDouble
 parameter C in the c-M relation
 
double m_mass = par::defaultDouble
 halo mass
 
double m_concentration = par::defaultDouble
 halo concentration
 
double m_redshift = par::defaultDouble
 halo redshift
 

Detailed Description

The class HaloProfile.

Definition at line 52 of file HaloProfile.h.

Constructor & Destructor Documentation

◆ HaloProfile() [1/2]

cbl::cosmology::HaloProfile::HaloProfile ( const cbl::cosmology::Cosmology  cosmology,
const double  redshift,
const double  conc,
const double  Mass,
const double  Delta,
const std::string  profile_author,
const std::string  halo_def,
const double  trunc_fact = 3.,
const bool  miscentering = false,
const bool  single_profile = false,
const double  sigma_off = 0.1,
const double  f_off = 1. 
)

constructor that sets the halo profile parameters. If the miscentering is considered, \(\Sigma(R)\) and \(\Delta\Sigma(R)\) are expressed as (see e.g. Bellagamba et al. 2019):

\(\Sigma_{\rm 1h}(R)=(1-f_{\rm off})\Sigma_{\rm cen}(R)+ f_{\rm off}\Sigma_{\rm off}(R),\,\,\,\,(1)\)

\(\Delta\Sigma_{\rm 1h}(R)=(1-f_{\rm off})\Delta\Sigma_{\rm cen}(R)+ f_{\rm off}\Delta\Sigma_{\rm off}(R).\,\,\,\,(2)\)

Parameters
cosmologythe cosmology
redshiftthe redshift
conchalo concentration
Masshalo mass
Deltaoverdensity factor which needs to be multiplied to the critical density in order to define an overdensity
profile_authorthe density profile author(s); available options are: "NFW" \(\rightarrow\) Navarro-Frenk-White profile; "NFW_trunc" \(\rightarrow\) truncated Navarro-Frenk-White profile
halo_defthe halo definition; available options are: "vir" \(\rightarrow\) all matter within the radius \(r_{vir}\) for which the mean internal density is cbl::cosmology::Cosmology::Delta_c times the critical density \(\rho_{crit}=3H^2/8\pi G\); "critical" \(\rightarrow\) all matter within the radius \(r_{200}\) for which the mean internal density is 200 times the critical density; "mean" \(\rightarrow\) all matter within the radius \(r_{\Delta}\) for which the mean internal density is \(\Delta\Omega_m\) times the critical mean background density
trunc_factfactor \(F_t\) defining the truncation radius, that is \(r_t = F_tr_{\Delta}\)
miscenteringif true, account for the miscentering
single_profileif true, the miscetering is related to a single cluster. Otherwise, a population of clusters is considered (the latter is the case of stacked weak lensing analyses)
sigma_offif single_profile=false, this is the standard deviation of the miscentered cluster population, \(\sigma_{\rm off}\), in Mpc \(/h\). Otherwise, it is the miscentering of a single cluster profile
f_offfraction of miscentered clusters, \(f_{\rm off}\). Used only if single_profile=false
Warning
If the miscentering is considered, the methods HaloProfile::Sigma() and HaloProfile::DeltaSigma() include the miscentering contribution (i.e. Eq. 1 and Eq. 2).
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities (e.g. \(\Sigma\), \(\Delta\Sigma\)) are in units of \(h\)

Definition at line 45 of file HaloProfile.cpp.

◆ HaloProfile() [2/2]

cbl::cosmology::HaloProfile::HaloProfile ( const cbl::cosmology::Cosmology  cosmology,
const double  redshift,
const std::string  cM_author,
const double  Mass,
const double  Delta,
const std::string  profile_author,
const std::string  halo_def,
const double  trunc_fact = 3.,
const bool  miscentering = false,
const bool  single_profile = false,
const double  sigma_off = 0.1,
const double  f_off = 1. 
)

constructor that sets the halo profile parameters, assuming a concentration-mass relation. If the miscentering is considered, \(\Sigma(R)\) and \(\Delta\Sigma(R)\) are expressed as (see e.g. Bellagamba et al. 2019):

\(\Sigma_{\rm 1h}(R)=(1-f_{\rm off})\Sigma_{\rm cen}(R)+ f_{\rm off}\Sigma_{\rm off}(R),\,\,\,\,(1)\)

\(\Delta\Sigma_{\rm 1h}(R)=(1-f_{\rm off})\Delta\Sigma_{\rm cen}(R)+ f_{\rm off}\Delta\Sigma_{\rm off}(R).\,\,\,\,(2)\)

Parameters
cosmologythe cosmology
redshiftthe redshift
cM_authorauthor(s) who proposed the concentration-mass relation. Possibilities are: "Duffy" (Duffy et al. 2008)
Masshalo mass
Deltaoverdensity factor which needs to be multiplied to the critical density in order to define an overdensity
profile_authorthe density profile author(s); available options are: "NFW" \(\rightarrow\) Navarro-Frenk-White profile; "NFW_trunc" \(\rightarrow\) truncated Navarro-Frenk-White profile
halo_defthe halo definition; available options are: "vir" \(\rightarrow\) all matter within the radius \(r_{vir}\) for which the mean internal density is cbl::cosmology::Cosmology::Delta_c times the critical density \(\rho_{crit}=3H^2/8\pi G\); "critical" \(\rightarrow\) all matter within the radius \(r_{200}\) for which the mean internal density is 200 times the critical density; "mean" \(\rightarrow\) all matter within the radius \(r_{\Delta}\) for which the mean internal density is \(\Delta\Omega_m\) times the critical mean background density
trunc_factfactor \(F_t\) defining the truncation radius, that is \(r_t = F_tr_{\Delta}\)
miscenteringif true, account for the miscentering
single_profileif true, the miscetering is related to a single cluster. Otherwise, a population of clusters is considered (the latter is the case of stacked weak lensing analyses)
sigma_offif single_profile=false, this is the standard deviation of the miscentered cluster population, \(\sigma_{\rm off}\), in Mpc \(/h\). Otherwise, it is the miscentering of a single cluster profile
f_offfraction of miscentered clusters, \(f_{\rm off}\). Used only if single_profile=false
Warning
the Duffy et al. concentrantion-mass relation refers to the 0<z<2 redshift range, obtained from their full samples (see Table 1 of Duffy et al. 2008).
If the miscentering is considered, the methods HaloProfile::Sigma() and HaloProfile::DeltaSigma() include the miscentering contribution (i.e. Eq. 1 and Eq. 2).
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities (e.g. \(\Sigma\), \(\Delta\Sigma\)) are in units of \(h\)

Definition at line 54 of file HaloProfile.cpp.

Member Function Documentation

◆ concentration()

double cbl::cosmology::HaloProfile::concentration ( )
inline

get the concentration

Returns
the concentration of the halo
Warning
if a concentration-mass relation is set, the output concentration is derived from such relation

Definition at line 567 of file HaloProfile.h.

◆ concentration2()

double cbl::cosmology::HaloProfile::concentration2 ( const double  Vmax,
const double  Rmax 
) const

compute the halo concentration

Author
Carlo Giocoli
cgioc.nosp@m.oli@.nosp@m.gmail.nosp@m..com
Parameters
VmaxVmax
RmaxRmax
Returns
the halo concentration

Definition at line 861 of file HaloProfile.cpp.

◆ Delta()

double cbl::cosmology::HaloProfile::Delta ( )
inline

get the overdensity

Returns
the overdensity which needs to be multiplied to the critical density in order to define an overdensity. Given an overdensity factor \(\Delta\), this function returns \(\Delta\) for a critical overdensity, \(\Delta\Omega_{\rm m}(z)\) for a mean overdensity, and cbl::cosmology::Cosmology::Delta_c for a overdensity defined with respect to the critical density.

Definition at line 592 of file HaloProfile.h.

◆ DeltaSigma()

std::vector< double > cbl::cosmology::HaloProfile::DeltaSigma ( const std::vector< double >  rad)

the total excess surface density profile of the halo.

Parameters
radradius
Returns
the halo excess surface density profile in \(M_{\odot} / \)pc^2
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 757 of file HaloProfile.cpp.

◆ DeltaSigma_2h() [1/2]

std::vector< double > cbl::cosmology::HaloProfile::DeltaSigma_2h ( const std::vector< double >  rad,
const double  bias,
const std::string  method_Pk = "EisensteinHu",
const std::string  interp_type = "Linear",
const bool  NL = false 
)

the 2-halo contribution to the excess surface density profile of the halo.

Parameters
radradius
biasthe halo bias value
method_Pkmethod used for the computation of the power spectrum. Valid choices for method_Pk are: CAMB [http://camb.info/], CLASS [http://class-code.net/], MPTbreeze-v1 [http://arxiv.org/abs/1207.1465], EisensteinHu [http://background.uchicago.edu/~whu/transfer/transferpage.html]
interp_typemethod to interpolate the power spectrum. Possibilities are: "Linear", "Poly", "Spline", "Spline_periodic", "Akima", "Akima_periodic", "Steffen"
NLfalse \(\rightarrow\) linear power spectrum; true \(\rightarrow\) non-linear power spectrum
Returns
the halo excess surface density profile in \(M_{\odot} / \)pc^2
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 790 of file HaloProfile.cpp.

◆ DeltaSigma_2h() [2/2]

std::vector< double > cbl::cosmology::HaloProfile::DeltaSigma_2h ( const std::vector< double >  rad,
const std::string  bias_author,
const std::string  method_Pk = "EisensteinHu",
const std::string  interp_type = "Linear",
const bool  NL = false 
)

the 2-halo contribution to the excess surface density profile of the halo, computed by assuming a halo bias model.

Parameters
radradius
bias_authorauthor(s) who proposed the bias; valid authors are: ST99 (Sheth & Tormen 1999), SMT01 (Sheth, Mo & Tormen 2001), SMT01_WL04 (Sheth, Mo & Tormen 2001 with the correction of Warren 2004), Tinker (Tinker et al. 2010)
method_Pkmethod used for the computation of the power spectrum. Valid choices for method_Pk are: CAMB [http://camb.info/], CLASS [http://class-code.net/], MPTbreeze-v1 [http://arxiv.org/abs/1207.1465], EisensteinHu [http://background.uchicago.edu/~whu/transfer/transferpage.html]
interp_typemethod to interpolate the power spectrum. Possibilities are: "Linear", "Poly", "Spline", "Spline_periodic", "Akima", "Akima_periodic", "Steffen"
NLfalse \(\rightarrow\) linear power spectrum; true \(\rightarrow\) non-linear power spectrum
Returns
the halo excess surface density profile in \(M_{\odot} / \)pc^2
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 835 of file HaloProfile.cpp.

◆ DeltaSigma_cen()

std::vector< double > cbl::cosmology::HaloProfile::DeltaSigma_cen ( const std::vector< double >  rad)

the centered contribution to the excess surface density profile of the halo.

Parameters
radradius
Returns
the halo excess surface density profile in \(M_{\odot} / \)pc^2
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 778 of file HaloProfile.cpp.

◆ DeltaSigma_mis()

std::vector< double > cbl::cosmology::HaloProfile::DeltaSigma_mis ( const std::vector< double >  rad)

the miscentering contribution to the excess surface density profile of the halo.

Parameters
radradius
Returns
the halo excess surface density profile in \(M_{\odot} / \)pc^2
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 766 of file HaloProfile.cpp.

◆ density_profile_3D()

std::vector< double > cbl::cosmology::HaloProfile::density_profile_3D ( const std::vector< double >  rad)

the 3D halo density profile.

Parameters
radradius
Returns
the halo density profile
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 655 of file HaloProfile.cpp.

◆ density_profile_FourierSpace()

double cbl::cosmology::HaloProfile::density_profile_FourierSpace ( const double  kk)

the Fourier transform of the normalised halo density

this function computes the Fourier transform of the normalised density distribution of dark matter haloes; the Navarro-Frenk-White profile is the only one currently implemented (see e.g. eq. 81 of Cooray & Sheth 2002 and 70 of van den Bosch et al. 2012)

\[\tilde{u}_h(k, M_h, z) = \frac{4\pi\rho_sr_s^3}{M_h} \left[\cos\mu\left[{\rm Ci}(\mu+\mu c) - {\rm Ci}(\mu)\right] + \sin\mu\left[{\rm Si}(\mu+\mu c) - {\rm Si}(\mu)\right] - \frac{\sin\mu c}{\mu+\mu c}\right]\]

where

\[\mu\equiv kr_s\,,\]

\[\rho_s = \frac{\rho_{crit}\Delta_c}{3} \frac{c^3}{\ln(1+c)-c/(1+c)}\,,\]

\[{\rm Ci}(x)=-\int_x^\infty\frac{\cos t}{t}\,{\rm d}t\,,\]

\[{\rm Si}(x)=-\int_x^\infty\frac{\sin t}{t}\,{\rm d}t\]

the relation between the halo concentration, \(c=c_{vir}=r_{vir}/r_s\), and halo mass, \(M_h\), is computed by cbl::modelling::twopt::concentration; \(\Delta_c(z)\) is computed by cbl::cosmology::Cosmology::Delta_c and \(\rho_{crit}(z)\) is computed by cbl::cosmology::Cosmology::rho_crit; \(r_{vir}(M_h, z)\) is computed by cbl::cosmology::Cosmology::r_vir

Parameters
kkthe wave vector module at which the model is computed
Returns
the halo density profile

Definition at line 847 of file HaloProfile.cpp.

◆ m_concentration_Duffy()

double cbl::cosmology::HaloProfile::m_concentration_Duffy ( )
private

The concentration-mass relation by Duffy et al. (2008):

\[c(M_h, z) = A(M_h/M_{pivot})^B\,(1+z)^C\]

Returns
the concentration

Definition at line 64 of file HaloProfile.cpp.

◆ m_DeltaSigma_cen()

std::vector< double > cbl::cosmology::HaloProfile::m_DeltaSigma_cen ( const std::vector< double >  rad)
private

Excess surface density without the miscentering contribution.

Parameters
radradius
Returns
the excess surface density

Definition at line 412 of file HaloProfile.cpp.

◆ m_DeltaSigma_including_miscentering()

std::vector< double > cbl::cosmology::HaloProfile::m_DeltaSigma_including_miscentering ( const std::vector< double >  rad)
private

Excess surface density including the miscentering contribution.

Parameters
radradius
Returns
the excess surface density

Definition at line 464 of file HaloProfile.cpp.

◆ m_DeltaSigma_mis()

std::vector< double > cbl::cosmology::HaloProfile::m_DeltaSigma_mis ( const std::vector< double >  rad)
private

the mis-centered excess surface density profile of the halo (see e.g. Bellagamba et al. 2019).

Parameters
radradius
Returns
the halo excess surface density profile in \(M_{\odot} / \)pc^2

Definition at line 429 of file HaloProfile.cpp.

◆ m_return_set_concentration()

double cbl::cosmology::HaloProfile::m_return_set_concentration ( )
inlineprivate

Return the concentration set by the user.

Returns
the concentration

Definition at line 226 of file HaloProfile.h.

◆ m_rho_NFW()

std::vector< double > cbl::cosmology::HaloProfile::m_rho_NFW ( const std::vector< double >  rad)
private

The NFW 3D density profile.

Parameters
radradius
Returns
the NFW density profile

Definition at line 108 of file HaloProfile.cpp.

◆ m_rho_NFW_trunc()

std::vector< double > cbl::cosmology::HaloProfile::m_rho_NFW_trunc ( const std::vector< double >  rad)
private

The truncated NFW 3D density profile.

Parameters
radradius
Returns
the truncated NFW density profile

Definition at line 128 of file HaloProfile.cpp.

◆ m_rho_s_NFW()

double cbl::cosmology::HaloProfile::m_rho_s_NFW ( )
private

Characteristic density of the NFW profile.

Returns
the NFW characteristic density

Definition at line 76 of file HaloProfile.cpp.

◆ m_rho_s_NFW_trunc()

double cbl::cosmology::HaloProfile::m_rho_s_NFW_trunc ( )
private

Characteristic density of the truncated NFW profile.

Returns
the truncated NFW characteristic density

Definition at line 89 of file HaloProfile.cpp.

◆ m_set_cM_relation()

void cbl::cosmology::HaloProfile::m_set_cM_relation ( const std::string  cM_author = "Duffy")
private

set the concentration-mass relation. A function is set, computing the concentration of a dark matter halo of a given a mass, at a given redshift; the models implemented are the following:

  • Duffy et al. 2008:

    \[c(M_h, z) = A(M_h/M_{pivot})^B\,(1+z)^C\]

Parameters
cM_authorauthor(s) who proposed the concentration-mass relation. Possibilities are: "Duffy" (Duffy et al. 2008)
Warning
the Duffy et al. concentrantion-mass relation refers to the 0<z<2 redshift range, obtained from their full samples (see Table 1 of Duffy et al. 2008).

Definition at line 484 of file HaloProfile.cpp.

◆ m_set_profile()

void cbl::cosmology::HaloProfile::m_set_profile ( const cbl::cosmology::Cosmology  cosmology,
const double  redshift,
const double  conc,
const double  Mass,
const double  Delta,
const std::string  profile_author,
const std::string  halo_def,
const double  trunc_fact = 3.,
const bool  miscentering = false,
const bool  single_profile = false,
const double  sigma_off = 0.1,
const double  f_off = 1. 
)
private

private function that sets the halo profile parameters. If the miscentering is considered, \(\Sigma(R)\) and \(\Delta\Sigma(R)\) are expressed as (see e.g. Bellagamba et al. 2019):

\(\Sigma_{\rm 1h}(R)=(1-f_{\rm off})\Sigma_{\rm cen}(R)+ f_{\rm off}\Sigma_{\rm off}(R),\,\,\,\,(1)\)

\(\Delta\Sigma_{\rm 1h}(R)=(1-f_{\rm off})\Delta\Sigma_{\rm cen}(R)+ f_{\rm off}\Delta\Sigma_{\rm off}(R).\,\,\,\,(2)\)

Parameters
cosmologythe cosmology
redshiftthe redshift
conchalo concentration
Masshalo mass
Deltaoverdensity factor which needs to be multiplied to the critical density in order to define an overdensity
profile_authorthe density profile author(s); available options are: "NFW" \(\rightarrow\) Navarro-Frenk-White profile; "NFW_trunc" \(\rightarrow\) truncated Navarro-Frenk-White profile
halo_defthe halo definition; available options are: "vir" \(\rightarrow\) all matter within the radius \(r_{vir}\) for which the mean internal density is cbl::cosmology::Cosmology::Delta_c times the critical density \(\rho_{crit}=3H^2/8\pi G\); "critical" \(\rightarrow\) all matter within the radius \(r_{200}\) for which the mean internal density is 200 times the critical density; "mean" \(\rightarrow\) all matter within the radius \(r_{\Delta}\) for which the mean internal density is \(\Delta\Omega_m\) times the critical mean background density
trunc_factfactor \(F_t\) defining the truncation radius, that is \(r_t = F_tr_{\Delta}\)
miscenteringif true, account for the miscentering
single_profileif true, the miscetering is related to a single cluster. Otherwise, a population of clusters is considered (the latter is the case of stacked weak lensing analyses)
sigma_offif single_profile=false, this is the standard deviation of the miscentered cluster population, \(\sigma_{\rm off}\), in Mpc \(/h\). Otherwise, it is the miscentering of a single cluster profile
f_offfraction of miscentered clusters, \(f_{\rm off}\). Used only if single_profile=false
Warning
If the miscentering is considered, the methods HaloProfile::Sigma() and HaloProfile::DeltaSigma() include the miscentering contribution (i.e. Eq. 1 and Eq. 2).
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities (e.g. \(\Sigma\), \(\Delta\Sigma\)) are in units of \(h\)

Definition at line 541 of file HaloProfile.cpp.

◆ m_Sigma_including_miscentering()

std::vector< double > cbl::cosmology::HaloProfile::m_Sigma_including_miscentering ( const std::vector< double >  rad)
private

Surface density including the miscentering contribution.

Parameters
radradius
Returns
the surface density

Definition at line 392 of file HaloProfile.cpp.

◆ m_Sigma_mean_NFW()

std::vector< double > cbl::cosmology::HaloProfile::m_Sigma_mean_NFW ( const std::vector< double >  rad)
private

Mean surface density of the NFW profile.

Parameters
radradius
Returns
the NFW surface density

Definition at line 178 of file HaloProfile.cpp.

◆ m_Sigma_mean_NFW_trunc()

std::vector< double > cbl::cosmology::HaloProfile::m_Sigma_mean_NFW_trunc ( const std::vector< double >  rad)
private

Mean surface density of the truncated NFW profile.

Parameters
radradius
Returns
the truncated NFW surface density

Definition at line 270 of file HaloProfile.cpp.

◆ m_Sigma_mis()

std::vector< double > cbl::cosmology::HaloProfile::m_Sigma_mis ( const std::vector< double >  rad)
private

the mis-centered surface density profile of the halo (see e.g. Bellagamba et al. 2019).

Parameters
radradius
Returns
the halo surface density profile in \(M_{\odot} / \)pc^2

Definition at line 322 of file HaloProfile.cpp.

◆ m_Sigma_NFW()

std::vector< double > cbl::cosmology::HaloProfile::m_Sigma_NFW ( const std::vector< double >  rad)
private

Surface density of the NFW profile.

Parameters
radradius
Returns
the NFW surface density

Definition at line 149 of file HaloProfile.cpp.

◆ m_Sigma_NFW_trunc()

std::vector< double > cbl::cosmology::HaloProfile::m_Sigma_NFW_trunc ( const std::vector< double >  rad)
private

Surface density of the truncated NFW profile.

Parameters
radradius
Returns
the truncated NFW surface density

Definition at line 209 of file HaloProfile.cpp.

◆ Mass_Delta()

double cbl::cosmology::HaloProfile::Mass_Delta ( const double  Mass,
const double  Delta_in,
const double  Delta_out,
const double  conc,
const bool  is_input_conc,
const double  rRmin_guess = 1.e-3,
const double  rRmax_guess = 10. 
) const

The halo mass converted to a different value of \(\Delta\), assuming the Navarro-Frenk-White density profile.

This function converts a given input mass \(M_\Delta\) to \(M_{\Delta^{new}}\) (e.g. \(M_{500} \rightarrow M_{200}\)).

Specifically, the algorithm currently implemented can be derived as follows. Given the Navarro-Frenk-White profile:

\[ \rho(r) = \frac{\rho_0}{\frac{r}{R_s} \left(1+\frac{r}{R_s} \right)^2} \]

the total halo mass contained within a radius \(R_\Delta\) is:

\[ M_\Delta = \int_0^{R_\Delta} 4\pi r^2\rho(r)dr = 4\pi \rho_0 R_s^3 \left[ \ln(1+c_\Delta)-\frac{c_\Delta}{1+c_\Delta} \right] \]

where the concentration is defined as \(c_\Delta\equiv R_\Delta/R_s\). Thus, we can write:

\[ \frac{\ln(1+c_\Delta) - \frac{c_\Delta}{1+c_\Delta}} {M_\Delta} = \frac{\ln(1+c_{\Delta^{new}}) - \frac{c_{\Delta^{new}}}{1+c_{\Delta^{new}}}} {M_{\Delta^{new}}} \]

\( M_\Delta \) can be written as follows:

\[ M_\Delta = \frac{4}{3}\pi\Delta\rho_{crit}R_\Delta^3 = \frac{\Delta}{\Delta^{new}} \left(\frac{R_\Delta}{R_{\Delta^{new}}} \right)^3 M_{\Delta^{new}} = \frac{\Delta}{\Delta^{new}} x^3 M_{\Delta^{new}} \]

where \( x\equiv R_\Delta/R_{\Delta^{new}} \). Thus we have:

\[ \frac{M_\Delta}{M_{\Delta^{new}}} \left[ \ln(1+c_{\Delta^{new}}) - \frac{c_{\Delta^{new}}}{1+c_{\Delta^{new}}} \right] - \left[ \ln(1+c_\Delta) - \frac{c_\Delta}{1+c_\Delta} \right] = 0 \]

where \( c_{\Delta^{new}} \equiv R_{\Delta^{new}}/R_s = c_\Delta R_{\Delta^{new}}/R_\Delta = c_\Delta/x \). The algorithm solves the above equation as a function of \( x \), providing in output:

\[ M_{\Delta^{new}} = \frac{1}{x^3}\frac{\Delta^{new}}{\Delta} M_\Delta \]

Parameters
Massthe input mass \(M_\Delta\) (e.g. \(M_{500}\))
Delta_inthe input \(\Delta\) (e.g. \(\Delta=500\))
Delta_outthe output \(\Delta^{new}\) (e.g. \(\Delta^{new}=200\))
concthe concentration related to either the input \(\Delta\), or the output \(\Delta^{new}\) (this is specified by is_input_conc)
is_input_conctrue \(\rightarrow\) the given concentration is related to the input \(\Delta\); false \(\rightarrow\) the given concentration is related to the output \(\Delta^{new}\)
rRmin_guessthe minimum guess value of \(r_\Delta/r_{\Delta^{new}}\) used by the gsl minimisation function
rRmax_guessthe maximum guess value of \(r_\Delta/r_{\Delta^{new}}\) used by the gsl minimisation function
Returns
\(M_{\Delta^{new}}\)
Warning
the current implementation assumes the Navarro-Frenk-White profile

Definition at line 878 of file HaloProfile.cpp.

◆ rho_s()

double cbl::cosmology::HaloProfile::rho_s ( )

the characteristic density of the 3D halo density profile.

Returns
the halo density profile

Definition at line 647 of file HaloProfile.cpp.

◆ set_concentration()

void cbl::cosmology::HaloProfile::set_concentration ( const double  conc = par::defaultDouble)
inline

set the private member m_concentration

Parameters
concconcentration

Definition at line 984 of file HaloProfile.h.

◆ set_cosmology()

void cbl::cosmology::HaloProfile::set_cosmology ( const cbl::cosmology::Cosmology  cosmology)

set the cosmological model

Parameters
cosmologythe cosmology

Definition at line 638 of file HaloProfile.cpp.

◆ set_f_off()

void cbl::cosmology::HaloProfile::set_f_off ( const double  f_off = par::defaultDouble)
inline

set the private member m_f_off

Parameters
f_offf_off

Definition at line 991 of file HaloProfile.h.

◆ set_mass()

void cbl::cosmology::HaloProfile::set_mass ( const double  mass = par::defaultDouble)
inline

set the private member m_mass

Parameters
massthe mass of the cluster

Definition at line 977 of file HaloProfile.h.

◆ set_sigma_off()

void cbl::cosmology::HaloProfile::set_sigma_off ( const double  sigma_off = par::defaultDouble)
inline

set the private member m_sigma_off

Parameters
sigma_offsigma_off

Definition at line 998 of file HaloProfile.h.

◆ set_trunc_fact()

void cbl::cosmology::HaloProfile::set_trunc_fact ( const double  trunc_fact = par::defaultDouble)
inline

set the private member m_trunc_fact

Parameters
trunc_facttrunc_fact

Definition at line 1005 of file HaloProfile.h.

◆ Sigma()

std::vector< double > cbl::cosmology::HaloProfile::Sigma ( const std::vector< double >  rad)

the total surface density profile of the halo.

Parameters
radradius
Returns
the halo surface density profile in \(M_{\odot} / \)pc^2
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 664 of file HaloProfile.cpp.

◆ Sigma_2h() [1/2]

std::vector< double > cbl::cosmology::HaloProfile::Sigma_2h ( const std::vector< double >  rad,
const double  bias,
const std::string  method_Pk = "EisensteinHu",
const std::string  interp_type = "Linear",
const bool  NL = false 
)

the 2-halo contribution to the total surface density profile of the halo.

Parameters
radradius
biasthe halo bias value
method_Pkmethod used for the computation of the power spectrum. Valid choices for method_Pk are: CAMB [http://camb.info/], CLASS [http://class-code.net/], MPTbreeze-v1 [http://arxiv.org/abs/1207.1465], EisensteinHu [http://background.uchicago.edu/~whu/transfer/transferpage.html]
interp_typemethod to interpolate the power spectrum. Possibilities are: "Linear", "Poly", "Spline", "Spline_periodic", "Akima", "Akima_periodic", "Steffen"
NLfalse \(\rightarrow\) linear power spectrum; true \(\rightarrow\) non-linear power spectrum
Returns
the halo surface density profile in \(M_{\odot} / \)pc^2
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 700 of file HaloProfile.cpp.

◆ Sigma_2h() [2/2]

std::vector< double > cbl::cosmology::HaloProfile::Sigma_2h ( const std::vector< double >  rad,
const std::string  bias_author,
const std::string  method_Pk = "EisensteinHu",
const std::string  interp_type = "Linear",
const bool  NL = false 
)

the 2-halo contribution to the total surface density profile of the halo, computed by assuming a halo bias model.

Parameters
radradius
bias_authorauthor(s) who proposed the bias; valid authors are: ST99 (Sheth & Tormen 1999), SMT01 (Sheth, Mo & Tormen 2001), SMT01_WL04 (Sheth, Mo & Tormen 2001 with the correction of Warren 2004), Tinker (Tinker et al. 2010)
method_Pkmethod used for the computation of the power spectrum. Valid choices for method_Pk are: CAMB [http://camb.info/], CLASS [http://class-code.net/], MPTbreeze-v1 [http://arxiv.org/abs/1207.1465], EisensteinHu [http://background.uchicago.edu/~whu/transfer/transferpage.html]
interp_typemethod to interpolate the power spectrum. Possibilities are: "Linear", "Poly", "Spline", "Spline_periodic", "Akima", "Akima_periodic", "Steffen"
NLfalse \(\rightarrow\) linear power spectrum; true \(\rightarrow\) non-linear power spectrum
Returns
the halo surface density profile in \(M_{\odot} / \)pc^2
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 745 of file HaloProfile.cpp.

◆ Sigma_cen()

std::vector< double > cbl::cosmology::HaloProfile::Sigma_cen ( const std::vector< double >  rad)

the centered contribution to the total surface density profile of the halo.

Parameters
radradius
Returns
the halo surface density profile in \(M_{\odot} / \)pc^2
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 688 of file HaloProfile.cpp.

◆ Sigma_mis()

std::vector< double > cbl::cosmology::HaloProfile::Sigma_mis ( const std::vector< double >  rad)

the miscentering contribution to the total surface density profile of the halo.

Parameters
radradius
Returns
the halo surface density profile in \(M_{\odot} / \)pc^2
Warning
If in the cbl::cosmology::Cosmology object the \(h\) units are set, both the input quantities (e.g. mass, radius) and the output quantities are in units of \(h\)

Definition at line 673 of file HaloProfile.cpp.


The documentation for this class was generated from the following files: