Halo Modeling
Module for DM Halo calculations
- frb.halos.models.init_hmf()[source]
Initialize the Aemulus Halo Mass Function
WARNING: This uses the original version which codes Tinker+2008 We may refactor to use the more accurate, new version
Returns:
- frb.halos.models.frac_in_halos(zvals, Mlow, Mhigh, rmax=1.0)[source]
- Calculate the fraction of matter in collapsed halos
over a mass range and at a given redshift
Note that the fraction of DM associated with these halos will be scaled down by an additional factor of f_diffuse
Requires Aemulus HMF to be installed
- Parameters:
zvals – ndarray
Mlow – float In h^-1 units already so this will be applied for the halo mass function
Mhigh – float In h^-1 units already
rmax – float Extent of the halo in units of rvir WARNING: This calculation assumes a single concentration for all halos
- Returns:
- ndarray
rho_halo / rho_m
- Return type:
ratios
- frb.halos.models.halo_incidence(Mlow, zFRB, radius=None, hmfe=None, Mhigh=1e+16, nsample=20, cumul=False)[source]
Calculate the (approximate) average number of intersections to halos of a given minimum mass to a given zFRB.
Requires Aemulus HMF to be installed
- Parameters:
Mlow – float Mass of minimum halo in Solar masses The code deals with h^-1 factors so that you do not The minimum value is 2e10
zFRB – float Redshift of the FRB
radius –
Quantity, optional The calculation will specify this radius as rvir derived from
Mlow unless this is specified. And this rvir will vary with redshift
hmfe (hmf.hmf_emulator, optional) – Halo mass function emulator from Aeumulus
Mhigh – float, optional Mass of maximum halo in Solar masses
nsammple – int, optional Number of samplings in redshift 20 should be enough
cumul – bool, optional Return the cumulative quantities instead
- Returns:
If cumul is False Navg: float
Number of average intersections
elif cumul is True zeval: ndarray Ncumul: ndarray
- frb.halos.models.build_grid(z_FRB=1.0, ntrial=10, seed=12345, Mlow=10000000000.0, r_max=2.0, outfile=None, dz_box=0.1, dz_grid=0.01, f_hot=0.75, verbose=True)[source]
Generate a universe of dark matter halos with DM measurements Mainly an internal function for generating useful output grids.
Requires the Aemulus Halo Mass function
- Parameters:
z_FRB – float, optional
ntrial – int, optional
seed – int, optional
Mlow – float, optional h^-1 mass
r_max – float, optional Extent of the halo in units of rvir
outfile – str, optional Write
dz_box – float, optional Size of the slice of the universe for each sub-calculation
dz_grid – float, optional redshift spacing in the DM grid
f_hot – float Fraction of the cosmic fraction of matter in diffuse gas (for DM)
- Returns:
ndarray (ntrial, nz) halo_tbl: Table
Table of all the halos intersected
- Return type:
DM_grid
- frb.halos.models.stellarmass_from_halomass(log_Mhalo, z=0, params=None)[source]
Stellar mass from Halo Mass from Moster+2013 https://doi.org/10.1093/mnras/sts261
- frb.halos.models.halomass_from_stellarmass(log_mstar, z=0, randomize=False)[source]
Halo mass from Stellar mass (Moster+2013). Inverts the function stellarmass_from_halomass numerically.
- Parameters:
log_mstar (float or numpy.ndarray) – log_10 stellar mass in solar mass units.
z (float, optional) – galaxy redshift
- Returns:
- log_10 halo mass
in solar mass units.
- Return type:
log_Mhalo (float)
- frb.halos.models.stellarmass_from_halomass_kravtsov(log_mhalo)[source]
Stellar mass from Halo Mass from Kravtsov+2018. https://ui.adsabs.harvard.edu/abs/2018AstL…44….8K/abstract
Caution: This relation is valid for low z (z~0). Higher z values may require a scaled relation. :param log_mhalo: log_10 halo mass :type log_mhalo: float
- Returns:
log_10 galaxy
- Return type:
log_mstar (float)
- frb.halos.models.halomass_from_stellarmass_kravtsov(log_mstar)[source]
Inverts the function frb.halos.models.stellarmass_from_halomass_kravtsov. :param log_mstar: log_10 stellar mass :type log_mstar: float or numpy.ndarray
- Returns:
log_10 halo mass
- Return type:
log_mhalo (float)
- class frb.halos.models.ModifiedNFW(log_Mhalo=12.2, c=7.67, f_hot: float = None, log_MCGM: float = None, alpha=0.0, y0=1.0, z=0.0, cosmo=FlatLambdaCDM(name='Planck18', H0=<Quantity 67.66 km / (Mpc s)>, Om0=0.30966, Tcmb0=<Quantity 2.7255 K>, Neff=3.046, m_nu=<Quantity [0., 0., 0.06] eV>, Ob0=0.04897), fb=0.15814118710844152, **kwargs)[source]
Bases:
objectGenerate a modified NFW model, e.g. Mathews & Prochaska 2017 for the hot, virialized gas.
- Parameters:
log_Mhalo – float, optional log10 of the Halo mass [dark matter + baryons] (solar masses)
log_MCGM – float, optional log10 of the CGM mass (solar masses) If provided, this sets f_hot [this is its only use]
c – float, optional concentration of the halo
f_hot – float, optional Fraction of the total “expected” baryons in this hot (aka CGM) phase Will likely use this for all diffuse gas
alpha – float, optional Parameter to modify NFW profile power-law
y0 – float, optional Parameter to modify NFW profile position.
z – float, optional Redshift of the halo
cosmo – astropy cosmology, optional Cosmology of the universe.
- H0
Quantity; Hubble constant
- fb
float; Cosmic fraction of baryons (stars+dust+gas) in the entire halo Default to 0.16
- r200
Quantity Virial radius
- rho0
Quantity Density normalization
- M_b
Quantity Mass in baryons of the
- __init__(log_Mhalo=12.2, c=7.67, f_hot: float = None, log_MCGM: float = None, alpha=0.0, y0=1.0, z=0.0, cosmo=FlatLambdaCDM(name='Planck18', H0=<Quantity 67.66 km / (Mpc s)>, Om0=0.30966, Tcmb0=<Quantity 2.7255 K>, Neff=3.046, m_nu=<Quantity [0., 0., 0.06] eV>, Ob0=0.04897), fb=0.15814118710844152, **kwargs)[source]
- setup_param(cosmo)[source]
Setup key parameters of the model
- Parameters:
cosmo – astropy cosmology Cosmology of the universe
- fy_b(y)[source]
Enclosed mass function for the baryons
- Parameters
y: float or ndarray
- Returns:
f_y – Enclosed mass
- Return type:
float or ndarray
- ne(xyz)[source]
Calculate n_e from n_H with a correction for Helium Assume 25% mass is Helium and both electrons have been stripped
- Parameters:
xyz (ndarray (3, npoints)) – Coordinate(s) in kpc
- Returns:
n_e – electron density in cm**-3
- Return type:
float or ndarray
- nH(xyz)[source]
Calculate the Hydrogen number density Includes a correction for Helium
- Parameters:
xyz (ndarray) – Coordinate(s) in kpc
- Returns:
nH – Density in cm**-3
- Return type:
float or ndarray
- rho_b(xyz)[source]
Mass density in baryons in the halo; modified
- Parameters:
xyz (ndarray) – Position (assumes kpc)
- Returns:
rho – Density in g / cm**-3
- Return type:
Quantity
- Ne_Rperp(Rperp, step_size=<Quantity 0.1 kpc>, rmax=1.0, add_units=True, cumul=False)[source]
Calculate N_e at an input impact parameter Rperp Just a simple sum in steps of step_size
This integrates through the entire halo. Use half if this is for the host
- Parameters:
- Returns:
if cumul –
- zval: ndarray (kpc)
z-values where z=0 is the midplane
- Ne_cumul: ndarray
Cumulative Ne values (pc cm**-3)
else –
- Ne: Quantity
Column density of total electrons
- RM_Rperp(Rperp, Bparallel, step_size=<Quantity 0.1 kpc>, rmax=1.0, add_units=True, cumul=False, zmax=None)[source]
Calculate RM at an input impact parameter Rperp Just a simple sum in steps of step_size Assumes a constant Magnetic field
- Parameters:
Rperp (Quantity) – Impact parameter, typically in kpc
(Quantity) (Bparallel) – Magnetic field
step_size (Quantity, optional) – Step size used for numerical integration (sum)
rmax (float, optional) – Maximum radius for integration in units of r200
add_units (bool, optional) – Speed up calculations by avoiding units
cumul (bool, optional)
zmax (float, optional) – Maximum distance along the sightline to integrate. Default is rmax*rvir
- Returns:
if cumul –
- zval: ndarray (kpc)
z-values where z=0 is the midplane
- Ne_cumul: ndarray
Cumulative Ne values (pc cm**-3)
else –
- RM: Quantity
Column density of total electrons
- mass_r(r, step_size=<Quantity 0.1 kpc>)[source]
Calculate baryonic halo mass (not total) to a given radius Just a simple sum in steps of step_size
- Parameters:
r (Quantity) – Radius, typically in kpc
step_size (Quantity, optional) – Step size used for numerical integration (sum)
- Returns:
Mr – Enclosed baryonic mass within r Msun units
- Return type:
Quantity
- class frb.halos.models.MB04(Rc=<Quantity 167. kpc>, log_Mhalo=12.2, c=7.67, f_hot=0.75, **kwargs)[source]
Bases:
ModifiedNFWHalo based on the Maller & Bullock (2004) model of virialized halo gas.
- Parameters:
Rc – Quantity cooling radius
- class frb.halos.models.YF17(log_Mhalo=12.18, c=7.67, f_hot=0.75, **kwargs)[source]
Bases:
ModifiedNFWFaerman et al (2017) model of the Milky Way
For the un-normalized density profile, we adopt the average of the warm and hot components in
- class frb.halos.models.MB15(log_Mhalo=12.18, c=7.67, f_hot=0.75, **kwargs)[source]
Bases:
ModifiedNFWEncodes the Galactic halo profile from Miller & Bregman 2015, ApJ, 800, 14 https://ui.adsabs.harvard.edu/abs/2015ApJ…800…14M/abstract
The default normalization and beta values are taken from their Table 2, last row. The models presented there do not appear to vary too much.
- class frb.halos.models.MilkyWay(log_Mhalo=12.18, c=7.67, f_hot=0.75, alpha=2, y0=2, **kwargs)[source]
Bases:
ModifiedNFWFiducial model for the Galaxy
Halo mass follows latest constraints
Density profile is similar to Maller & Bullock 2004
- class frb.halos.models.M31(log_Mhalo=12.18, c=7.67, alpha=2, y0=2, f_hot=0.75, **kwargs)[source]
Bases:
ModifiedNFWPreferred model for M31
Taking mass from van der Marel 2012
- class frb.halos.models.LMC(log_Mhalo=np.float64(11.0), c=12.1, f_hot=0.75, alpha=2, y0=2, **kwargs)[source]
Bases:
ModifiedNFWPreferred model for LMC
Taking data from D’Onghia & Fox ARAA 2016
- Mass updated according to
https://ui.adsabs.harvard.edu/abs/2019MNRAS.487.2685E/abstract
- class frb.halos.models.SMC(log_Mhalo=np.float64(9.380211241711606), c=15.0, f_hot=0.75, alpha=2, y0=2, **kwargs)[source]
Bases:
ModifiedNFWPreferred model for SMC
Taking data from D’Onghia & Fox ARAA 2016
- class frb.halos.models.M33(log_Mhalo=np.float64(11.698970004336019), c=8.36, f_hot=0.75, alpha=2, y0=2, **kwargs)[source]
Bases:
ModifiedNFWPreferred model for SMC
Taking data from Corbelli 2006
- class frb.halos.models.ICM(log_Mhalo=np.float64(14.698970004336019), c=5, f_hot=0.7, **kwargs)[source]
Bases:
ModifiedNFWIntracluster medium (ICM) model following the analysis of Vikhilnin et al. 2006
We scale the model to the profile fitted to A907
- setup_param(cosmo=None)[source]
Setup key parameters of the model
- Parameters:
cosmo – astropy cosmology Cosmology of the universe
- class frb.halos.models.Virgo(log_Mhalo=np.float64(14.880136251483002), **kwargs)[source]
Bases:
ICMParameterization of Virgo following the Planck Collaboration paper: A&A 596 A101 (2016)
- frb.halos.models.rad3d2(xyz)[source]
Calculate radius squared to x,y,z inputted Assumes the origin is 0,0,0
- Parameters:
xyz (Tuple or ndarray)
- Returns:
rad3d
- Return type:
float or ndarray
Module for DM halo calculations and galaxy halo models.
Constants
- frb.halos.models.m_p = np.float64(1.67262192369e-24)
float64(value=0, /) –
Double-precision floating-point number type, compatible with Python
floatand Cdouble.- Character code:
'd'- Canonical name:
numpy.double
- Alias on this platform (Linux x86_64):
numpy.float64: 64-bit precision floating-point number type: sign bit, 11 bits exponent, 52 bits mantissa.
Proton mass in CGS units for optimized calculations.
Functions
init_hmf
- frb.halos.models.init_hmf()[source]
Initialize the Aemulus Halo Mass Function
WARNING: This uses the original version which codes Tinker+2008 We may refactor to use the more accurate, new version
Returns:
Initialize the Aemulus Halo Mass Function.
Warning
Uses the original version which codes Tinker+2008. May be refactored to use the more accurate new version.
Returns:
hmf_emulator- Initialized HMF emulator object
frac_in_halos
- frb.halos.models.frac_in_halos(zvals, Mlow, Mhigh, rmax=1.0)[source]
- Calculate the fraction of matter in collapsed halos
over a mass range and at a given redshift
Note that the fraction of DM associated with these halos will be scaled down by an additional factor of f_diffuse
Requires Aemulus HMF to be installed
- Parameters:
zvals – ndarray
Mlow – float In h^-1 units already so this will be applied for the halo mass function
Mhigh – float In h^-1 units already
rmax – float Extent of the halo in units of rvir WARNING: This calculation assumes a single concentration for all halos
- Returns:
- ndarray
rho_halo / rho_m
- Return type:
ratios
Calculate the fraction of matter in collapsed halos over a mass range and at a given redshift.
Note
The fraction of DM associated with these halos will be scaled down by an additional factor of f_diffuse.
Parameters:
zvals: ndarray - Redshift valuesMlow: float - Minimum halo mass in h^-1 unitsMhigh: float - Maximum halo mass in h^-1 unitsrmax: float, optional - Extent of halo in units of rvir (default: 1.0)
Returns:
ratios: ndarray - rho_halo / rho_m
Warning
This calculation assumes a single concentration for all halos.
stellarmass_from_halomass
- frb.halos.models.stellarmass_from_halomass(log_Mhalo, z=0, params=None)[source]
Stellar mass from Halo Mass from Moster+2013 https://doi.org/10.1093/mnras/sts261
- Parameters:
- Returns:
- log_10 galaxy stellar mass
in solar mass units.
- Return type:
log_mstar (float)
Stellar mass from Halo Mass using Moster+2013 relation.
Parameters:
log_Mhalo: float - log_10 halo mass in solar mass unitsz: float, optional - Halo redshift (default: 0)params: list, optional - Custom model parameters
Returns:
log_mstar: float - log_10 galaxy stellar mass in solar mass units
Reference: https://doi.org/10.1093/mnras/sts261
halomass_from_stellarmass
- frb.halos.models.halomass_from_stellarmass(log_mstar, z=0, randomize=False)[source]
Halo mass from Stellar mass (Moster+2013). Inverts the function stellarmass_from_halomass numerically.
- Parameters:
log_mstar (float or numpy.ndarray) – log_10 stellar mass in solar mass units.
z (float, optional) – galaxy redshift
- Returns:
- log_10 halo mass
in solar mass units.
- Return type:
log_Mhalo (float)
Halo mass from Stellar mass (Moster+2013). Numerically inverts stellarmass_from_halomass.
Parameters:
log_mstar: float or ndarray - log_10 stellar mass in solar mass unitsz: float, optional - Galaxy redshift (default: 0)randomize: bool, optional - Add scatter to the relation
Returns:
log_Mhalo: float - log_10 halo mass in solar mass units
stellarmass_from_halomass_kravtsov
- frb.halos.models.stellarmass_from_halomass_kravtsov(log_mhalo)[source]
Stellar mass from Halo Mass from Kravtsov+2018. https://ui.adsabs.harvard.edu/abs/2018AstL…44….8K/abstract
Caution: This relation is valid for low z (z~0). Higher z values may require a scaled relation. :param log_mhalo: log_10 halo mass :type log_mhalo: float
- Returns:
log_10 galaxy
- Return type:
log_mstar (float)
Stellar mass from Halo Mass using Kravtsov+2018 relation.
Caution
This relation is valid for low z (z~0). Higher z values may require a scaled relation.
Parameters:
log_mhalo: float - log_10 halo mass
Returns:
log_mstar: float - log_10 galaxy stellar mass
Reference: https://ui.adsabs.harvard.edu/abs/2018AstL…44….8K/abstract
halomass_from_stellarmass_kravtsov
- frb.halos.models.halomass_from_stellarmass_kravtsov(log_mstar)[source]
Inverts the function frb.halos.models.stellarmass_from_halomass_kravtsov. :param log_mstar: log_10 stellar mass :type log_mstar: float or numpy.ndarray
- Returns:
log_10 halo mass
- Return type:
log_mhalo (float)
Inverts stellarmass_from_halomass_kravtsov function.
Parameters:
log_mstar: float or ndarray - log_10 stellar mass
Returns:
log_mhalo: float - log_10 halo mass
rad3d2
- frb.halos.models.rad3d2(xyz)[source]
Calculate radius squared to x,y,z inputted Assumes the origin is 0,0,0
- Parameters:
xyz (Tuple or ndarray)
- Returns:
rad3d
- Return type:
float or ndarray
Calculate radius to x,y,z coordinates. Assumes origin is (0,0,0).
Parameters:
xyz: tuple or ndarray - 3D coordinates
Returns:
rad3d: float or ndarray - 3D radius squared
Classes
ModifiedNFW
- class frb.halos.models.ModifiedNFW(log_Mhalo=12.2, c=7.67, f_hot: float = None, log_MCGM: float = None, alpha=0.0, y0=1.0, z=0.0, cosmo=FlatLambdaCDM(name='Planck18', H0=<Quantity 67.66 km / (Mpc s)>, Om0=0.30966, Tcmb0=<Quantity 2.7255 K>, Neff=3.046, m_nu=<Quantity [0., 0., 0.06] eV>, Ob0=0.04897), fb=0.15814118710844152, **kwargs)[source]
Bases:
objectGenerate a modified NFW model, e.g. Mathews & Prochaska 2017 for the hot, virialized gas.
- Parameters:
log_Mhalo – float, optional log10 of the Halo mass [dark matter + baryons] (solar masses)
log_MCGM – float, optional log10 of the CGM mass (solar masses) If provided, this sets f_hot [this is its only use]
c – float, optional concentration of the halo
f_hot – float, optional Fraction of the total “expected” baryons in this hot (aka CGM) phase Will likely use this for all diffuse gas
alpha – float, optional Parameter to modify NFW profile power-law
y0 – float, optional Parameter to modify NFW profile position.
z – float, optional Redshift of the halo
cosmo – astropy cosmology, optional Cosmology of the universe.
- H0
Quantity; Hubble constant
- fb
float; Cosmic fraction of baryons (stars+dust+gas) in the entire halo Default to 0.16
- r200
Quantity Virial radius
- rho0
Quantity Density normalization
- M_b
Quantity Mass in baryons of the
Generate a modified NFW model for hot, virialized gas (e.g. Mathews & Prochaska 2017).
Parameters:
log_Mhalo: float, optional - log10 of halo mass in solar masses (default: 12.2)c: float, optional - Concentration of the halo (default: 7.67)f_hot: float, optional - Fraction of baryons in hot phase (default: 0.75)alpha: float, optional - Parameter to modify NFW profile power-law (default: 0)y0: float, optional - Parameter to modify NFW profile position (default: 1)z: float, optional - Redshift of the halo (default: 0)cosmo: astropy cosmology, optional - Cosmology of the universe
Key Attributes:
H0: Quantity - Hubble constantfb: float - Cosmic fraction of baryons (default: 0.16)r200: Quantity - Virial radiusrho0: Quantity - Density normalizationM_b: Quantity - Mass in baryons
Methods:
- setup_param(cosmo)[source]
Setup key parameters of the model
- Parameters:
cosmo – astropy cosmology Cosmology of the universe
Setup key parameters of the model.
- Ne_Rperp(Rperp, step_size=<Quantity 0.1 kpc>, rmax=1.0, add_units=True, cumul=False)[source]
Calculate N_e at an input impact parameter Rperp Just a simple sum in steps of step_size
This integrates through the entire halo. Use half if this is for the host
- Parameters:
- Returns:
if cumul –
- zval: ndarray (kpc)
z-values where z=0 is the midplane
- Ne_cumul: ndarray
Cumulative Ne values (pc cm**-3)
else –
- Ne: Quantity
Column density of total electrons
Calculate column density along a perpendicular path through the halo.
- __init__(log_Mhalo=12.2, c=7.67, f_hot: float = None, log_MCGM: float = None, alpha=0.0, y0=1.0, z=0.0, cosmo=FlatLambdaCDM(name='Planck18', H0=<Quantity 67.66 km / (Mpc s)>, Om0=0.30966, Tcmb0=<Quantity 2.7255 K>, Neff=3.046, m_nu=<Quantity [0., 0., 0.06] eV>, Ob0=0.04897), fb=0.15814118710844152, **kwargs)[source]
- setup_param(cosmo)[source]
Setup key parameters of the model
- Parameters:
cosmo – astropy cosmology Cosmology of the universe
- fy_b(y)[source]
Enclosed mass function for the baryons
- Parameters
y: float or ndarray
- Returns:
f_y – Enclosed mass
- Return type:
float or ndarray
- ne(xyz)[source]
Calculate n_e from n_H with a correction for Helium Assume 25% mass is Helium and both electrons have been stripped
- Parameters:
xyz (ndarray (3, npoints)) – Coordinate(s) in kpc
- Returns:
n_e – electron density in cm**-3
- Return type:
float or ndarray
- nH(xyz)[source]
Calculate the Hydrogen number density Includes a correction for Helium
- Parameters:
xyz (ndarray) – Coordinate(s) in kpc
- Returns:
nH – Density in cm**-3
- Return type:
float or ndarray
- rho_b(xyz)[source]
Mass density in baryons in the halo; modified
- Parameters:
xyz (ndarray) – Position (assumes kpc)
- Returns:
rho – Density in g / cm**-3
- Return type:
Quantity
- Ne_Rperp(Rperp, step_size=<Quantity 0.1 kpc>, rmax=1.0, add_units=True, cumul=False)[source]
Calculate N_e at an input impact parameter Rperp Just a simple sum in steps of step_size
This integrates through the entire halo. Use half if this is for the host
- Parameters:
- Returns:
if cumul –
- zval: ndarray (kpc)
z-values where z=0 is the midplane
- Ne_cumul: ndarray
Cumulative Ne values (pc cm**-3)
else –
- Ne: Quantity
Column density of total electrons
- RM_Rperp(Rperp, Bparallel, step_size=<Quantity 0.1 kpc>, rmax=1.0, add_units=True, cumul=False, zmax=None)[source]
Calculate RM at an input impact parameter Rperp Just a simple sum in steps of step_size Assumes a constant Magnetic field
- Parameters:
Rperp (Quantity) – Impact parameter, typically in kpc
(Quantity) (Bparallel) – Magnetic field
step_size (Quantity, optional) – Step size used for numerical integration (sum)
rmax (float, optional) – Maximum radius for integration in units of r200
add_units (bool, optional) – Speed up calculations by avoiding units
cumul (bool, optional)
zmax (float, optional) – Maximum distance along the sightline to integrate. Default is rmax*rvir
- Returns:
if cumul –
- zval: ndarray (kpc)
z-values where z=0 is the midplane
- Ne_cumul: ndarray
Cumulative Ne values (pc cm**-3)
else –
- RM: Quantity
Column density of total electrons
- mass_r(r, step_size=<Quantity 0.1 kpc>)[source]
Calculate baryonic halo mass (not total) to a given radius Just a simple sum in steps of step_size
- Parameters:
r (Quantity) – Radius, typically in kpc
step_size (Quantity, optional) – Step size used for numerical integration (sum)
- Returns:
Mr – Enclosed baryonic mass within r Msun units
- Return type:
Quantity
ICM
- class frb.halos.models.ICM(log_Mhalo=np.float64(14.698970004336019), c=5, f_hot=0.7, **kwargs)[source]
Bases:
ModifiedNFWIntracluster medium (ICM) model following the analysis of Vikhilnin et al. 2006
We scale the model to the profile fitted to A907
Intracluster Medium model, child of ModifiedNFW.
Implements the Miller & Bregman 2015 ICM model for galaxy clusters.
Methods:
- nH(xyz)[source]
Scale by He
- Parameters:
xyz
Returns:
Calculate the number density of Hydrogen.
Parameters:
xyz: ndarray - Coordinates in kpc
Returns:
ndarray- Number density in units of 1/cm³
- setup_param(cosmo=None)[source]
Setup key parameters of the model
- Parameters:
cosmo – astropy cosmology Cosmology of the universe
MilkyWay
- class frb.halos.models.MilkyWay(log_Mhalo=12.18, c=7.67, f_hot=0.75, alpha=2, y0=2, **kwargs)[source]
Bases:
ModifiedNFWFiducial model for the Galaxy
Halo mass follows latest constraints
Density profile is similar to Maller & Bullock 2004
Fiducial model for the Galaxy. Halo mass follows latest constraints.
Density profile is similar to Maller & Bullock 2004.
M31
- class frb.halos.models.M31(log_Mhalo=12.18, c=7.67, alpha=2, y0=2, f_hot=0.75, **kwargs)[source]
Bases:
ModifiedNFWPreferred model for M31
Taking mass from van der Marel 2012
Preferred model for M31. Mass from van der Marel 2012.
Attributes:
distance: Quantity - Distance from Sun (752 kpc)coord: SkyCoord - Coordinates of M31
Methods:
- DM_from_Galactic(scoord, **kwargs)[source]
Calculate DM through M31’s halo from the Sun given a direction
- Parameters:
scoord – SkyCoord Coordinates of the sightline
**kwargs – Passed to Ne_Rperp
- Returns:
- Quantity
Dispersion measure through M31’s halo
- Return type:
Calculate DM through M31’s halo from the Sun given a direction.
Parameters:
scoord: SkyCoord - Coordinates of the sightline
Returns:
DM: Quantity - Dispersion measure through M31’s halo
Examples
Basic halo model usage:
from frb.halos.models import ModifiedNFW, halomass_from_stellarmass
from astropy import units as u
# Create a halo model
halo = ModifiedNFW(log_Mhalo=12.5, z=0.3, f_hot=0.6)
# Calculate DM at 100 kpc offset
offset = 100 * u.kpc
dm_contribution = halo.Ne_Rperp(offset)
print(f"DM contribution: {dm_contribution}")
Stellar-halo mass relations:
from frb.halos.models import stellarmass_from_halomass, halomass_from_stellarmass
# Convert halo mass to stellar mass
log_mhalo = 12.0 # log solar masses
log_mstar = stellarmass_from_halomass(log_mhalo, z=0.5)
# Invert the relation
log_mhalo_recovered = halomass_from_stellarmass(log_mstar, z=0.5)
print(f"Halo mass: {log_mhalo}")
print(f"Stellar mass: {log_mstar:.2f}")
print(f"Recovered halo mass: {log_mhalo_recovered:.2f}")
Galaxy-specific models:
from frb.halos.models import MilkyWay, M31
from astropy.coordinates import SkyCoord
from astropy import units as u
# Milky Way model
mw = MilkyWay(log_Mhalo=12.1, f_hot=0.7)
# M31 model with DM calculation
m31 = M31(log_Mhalo=12.2)
target_coord = SkyCoord('01h33m51s', '+30d39m37s')
dm_m31 = m31.DM_from_Galactic(target_coord)
print(f"DM through M31: {dm_m31}")
Halo mass function calculations:
from frb.halos.models import frac_in_halos
import numpy as np
# Calculate matter fraction in halos
redshifts = np.array([0.1, 0.5, 1.0])
mass_min = 1e11 # Solar masses
mass_max = 1e15
fractions = frac_in_halos(redshifts, mass_min, mass_max)
for z, frac in zip(redshifts, fractions):
print(f"z={z:.1f}: {frac:.3f} of matter in halos")
Advanced usage with scatter:
# Include scatter in stellar-halo mass relation
log_mstar = 10.5
z = 0.3
# Multiple realizations with scatter
log_mhalos = []
for i in range(100):
log_mhalo = halomass_from_stellarmass(log_mstar, z=z, randomize=True)
log_mhalos.append(log_mhalo)
mean_mhalo = np.mean(log_mhalos)
std_mhalo = np.std(log_mhalos)
print(f"Mean log(Mhalo): {mean_mhalo:.2f} ± {std_mhalo:.2f}")