Public API
Models
Models in Slomo.jl are specified as new types whose attributes track the parameters of the model.
For instance[1], we can make a new model representing an NFW halo as
halo = Slomo.Halos.NFWModel(25.3, 3.7e6)The NFWModel has two parameters, rs and rhos, representing the scale radius and scale density of the halo in units of kpc and solar masses per cubic kpc, respectively. We can then query the enclosed mass of this model at a particular radius with the mass function, e.g.,
mass_at_5kpc = mass(halo, 5)Models are immutable composite types (i.e. structs in the usual Julia syntax). For those familiar with Python, these work like NamedTuples.
While it isn't possible to change the value of a parameter in a Model, you can easily make a copy of the model with different parameters using the Slomo.Models.update function:
denser_halo = update(halo, rhos=4e6)The Models module defines the basic functionality for working with DensityModels (a subtype of Models), including functions for evaluating the volume density, surface density, and enclosed mass profiles. All subtypes of DensityModel should at least implement the density() method. Most models implemented in Slomo already implement custom mass methods for computational convenience.
Slomo.Models — ModuleBasic model types and methods
Slomo.Models.density — Methoddensity(model::DensityModel, r)Local volume density at r.
Slomo.Models.density2d — Methoddensity2d(model::DensityModel, R)Local surface density at R. If not defined for a subtype of DensityModel, then calculate numerically as the Abel transform from the volume density profile.
Slomo.Models.mass — Methodmass(model::DensityModel, r)Enclosed mass within r. If not defined for a subtype of DensityModel, then calculate numerically as the integral of the volume density.
Slomo.Models.update — Methodupdate(model)Create a new model from the specified parameters.
Anisotropy
Slomo.Anisotropy — ModuleCollection of anisotropy models.
IsotropicModel: isotropic modelConstantBetaModel: constant beta modelRSBetaModel: flexible beta model from Read & Stegerbeta: orbital anisotropy profile
Slomo.Anisotropy.ConstantBetaModel — TypeConstantBetaModel(beta)Model that assumes beta is constant with galactocentric radius.
Slomo.Anisotropy.IsotropicModel — TypeIsotropicModel()Parameter-less model representing an isotropic system (i.e. beta = 0).
Slomo.Anisotropy.RSBetaModel — TypeRSBetaModel(beta0, betaInf, rbeta, nbeta)Flexible velocity anisotropy model from Read & Steger 2017 (eq. 9)
beta0: inner asymptotic anisotropybetaInf: outer asymptotic anisotropyrbeta: transition radiusnbeta: transition sharpness, higher is a faster transition
beta0, betaInf, nbeta = 0, 1, 2 corresponds to the Osipkov-Merritt profile
beta0, betaInf, nbeta = 0, 0.5, 1 corresponds to the Mamon-Lokas profile
Slomo.Anisotropy.beta — Methodbeta(model::AnisotropyModel, r)Orbital anisotropy as a function of radius. The orbital anisotropy is defined as
Tracers
Slomo.Tracers — ModuleCollection of tracer density models.
SersicModel : Representation of a Sersic surface density profile
Slomo.Tracers.SersicModel — TypeSersicModel(Re, n, Mtot)Representation of a Sersic surface density profile.
Re: the effective radius in kpc.n: the Sersic indexMtot: the total mass associated to the profile.
For the purposes of using the SersicModel to track the tracer population's density distribution, Mtot is somewhat arbitrary. For using it as a component of the total mass profile, this should take units of solar masses.
Halos
Slomo.Halos — ModuleDark matter halo models and relations.
Significant credit to Benedikt Diemer for developing Colossus.
While much of this is inspired by Colossus; any implementation faults are my own.
The default halo mass definition is "200c", i.e. the virial radius of the halo is defined to be the radius that encloses a density that is 200 times the critical density of the universe.
The default cosmology (from Slomo.CosmologyTools) is currently set to the values from Planck 2018 (h = 0.6766, Neff = 3.046, OmegaM = 0.3111).
Common HaloModel methods
scale_radius: scale radius of the halovirial_radius: virial radius of the halovirial_mass: virial mass of the haloconcentration: halo concentration (defined as rvir / rs)
HaloModel subtypes
NFWModel: Navarro, Frenk, and White (1997) profileCoreNFWModel: Core-NFW (Read et al. 2016) profileGNFWModel: generalized NFW model (i.e., αβγ where α = 1, β = 3)ABGModel: α-β-γ profile (double power law)EinastoModel: Einasto (1965) profileSolNFWModel: Soliton + NFW model (Marsh & Pop 2015)SolABGModel: Soliton + αβγ model (Wasserman et al. 2019)
Halo relations
hmcr: halo mass-concentration relation from Dutton & Maccio (2014)shmr: stellar-halo mass relation from Rodriguez-Puebla et al. (2017)abg_from_logshm: αβγ scaling with log(Mstar / Mhalo) from DiCintio et al. (2014)
Slomo.Halos.ABGModel — TypeABGModel(rs, rhos, alpha, beta, gamma)αβγ double-power law density model.
rs: scale radius in kpcrhos: scale density in Msun / kpc3alpha: transition sharpnessbeta: negative outer log slopegamma: negative inner log slope
Slomo.Halos.CoreNFWModel — TypeCoreNFWModel(rs, rhos, rc, nc)coreNFW halo density model from Read et. al (2016)
Read, Agertz, and Collins 2016
rs: scale radius in kpcrhos: scale density in Msun / kpc3rc: core radius in kpcnc: core index (0 -> no core, 1 -> full core)
Slomo.Halos.EinastoModel — TypeEinastoModel(rs, rhos, alpha)Einasto halo model
rs: scale radius in kpcrhos: scale density in Msun / kpc3alpha: Einasto shape index, lower is more cored
Slomo.Halos.GNFWModel — TypeGNFWModel(rs, rhos, gamma)Generalized NFW (gNFW) halo density model. This is a specific instance of a generalized Hernquist (i.e. αβγ) model with the transition parameter α = 1, the outer slope β = 3, and the inner slope, γ, left as a free parameter.
rs: scale radius in kpcrhos: scale density in Msun / kpc3gamma: negative of inner density log slope
Slomo.Halos.NFWModel — TypeNFWModel(rs, rhos)NFW halo density model. Default constructor will make a halo with Mvir = 1e12 Msun.
rs: scale radius in kpcrhos: scale density in Msun / kpc3
Slomo.Halos.SolABGModel — TypeSolABGModel(rs, rhos, alpha, beta, gamma, rsol, rhosol, repsilon)Soliton + αβγ model. For consistency, the αβγ density profile must match the soliton density profile at repsilon. If that is not the case, the constructor will update repsilon.
rs: NFW scale radius in kpcrhos: NFW scale density in Msun / kpc3alpha: transition sharpnessbeta: negative outer log slopegamma: negative inner log slopersol: soliton scale radius in kpcrhosol: soliton scale density in kpcrepsilon: transition radius
Slomo.Halos.SolNFWModel — TypeSolNFWModel(rs, rhos, rsol, rhosol, repsilon)Soliton NFW model. For consistency, the NFW density must match the soliton density at repsilon. If that is not the case, the constructor will update repsilon.
rs: NFW scale radius in kpcrhos: NFW scale density in Msun / kpc3rsol: soliton scale radius in kpcrhosol: soliton scale density in kpcrepsilon: transition radius
Slomo.Halos.ABG_from_virial — MethodABG_from_virial(Mvir, cvir, alpha, beta, gamma; <keyword arguments>)Construct a ABGModel halo from the virial mass, concentration, and shape parameters.
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensity
Slomo.Halos.ABG_from_virial — MethodABG_from_virial(Mvir, cvir, Mstar; <keyword arguments>)Construct a ABGModel halo from the virial mass, concentration, and stellar mass, using the scaling relations from DiCintio et al. 2014
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensity
Slomo.Halos.CoreNFW_from_virial — MethodCoreNFW_from_virial(Mvir, cvir, Re, t_sf; <keyword arguments>)Construct a CoreNFWModel from the scaling relations in Read et al. (2016).
Mvir: virial mass in Msuncvir: halo concentrationRe: effective radius of stellar distribution, in kpct_sf: time since start of star formation, in Gyr
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensity
Slomo.Halos.Einasto_from_virial — MethodEinasto_from_virial(Mvir, cvir, alpha; <keyword arguments>)Construct an Einasto halo from the halo mass, concentration, and shape parameter.
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensity
Slomo.Halos.GNFW_from_virial — MethodGNFW_from_virial(Mvir, cvir, gamma; <keyword arguments>)Construct a gNFW halo from the virial mass and concentration.
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensity
Slomo.Halos.NFW_from_virial — MethodNFW_from_virial(Mvir, cvir; <keyword arguments>)Construct an NFW halo from the virial mass and concentration.
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensity
Slomo.Halos.SolABG_from_virial — MethodSolABG_from_virial(Mvir, cvir, alpha, beta, gamma, m22; <keyword arguments>)Construct a Soliton + αβγ halo from the halo mass, concentration, shape parameters, and axion mass.
Arguments
rsol = nothing: for default, use the soliton size-halo mass scaling relationshipmdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensityx_start::Real = 2.0: factor of radius scale to start the search for the rootmaxevals:Int = 100: maximum number of function evalutions for the root finding
Slomo.Halos.SolNFW_from_virial — MethodSolNFW_from_virial(Mvir, cvir, m22; <keyword arguments>)Construct a Solition-NFW halo from the halo mass, concentration, and axion mass.
Arguments
rsol = nothing: for default, use the soliton size-halo mass scaling relationshipmdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensityx_start::Real = 2.0: factor of radius scale to start the search for the rootmaxevals:Int = 100: maximum number of function evalutions for the root finding
Slomo.Halos.abg_from_logshm — Methodabg_from_logshm(logshm)αβγ scalings with log10(Mstar / Mhalo) from DiCintio+2014a.
Slomo.Halos.concentration — Methodconcentration(halo::HaloModel; <keyword arguments>)Compute the halo concentration.
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensityxstart::Real = 10.0: factor of the scale radius to start the search for the rootrtol::Real = 1e-2: relative error tolerance for the root findingmaxevals:Int = 100: maximum number of function evalutions for the root finding
Slomo.Halos.hmcr — Methodhmcr(Mvir; <keyword arguments>)Compute the expected halo concentration from halo mass using the relation of Dutton & Maccio 2014 (equations 7, 10-13).
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensity
Slomo.Halos.hmcr_prior — Methodhmcr_prior(Mvir, cvir; <keyword arguments>)Compute the probability of drawing a halo with virial parameters (Mvir, cvir), using the halo mass-concentration relation of Dutton & Maccio 2014.
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensitysigma_logc::Real = 0.16: scatter in log_concentration at fixed Mvir (in dex)
Slomo.Halos.scale_radius — Methodscale_radius(halo::HaloModel)Return the scale radius of the halo. The convention adopted here is that the scale radius occurs where the logarithmic slope of the halo density profile is equal to -2. For the standard NFW model this is equal to rs, but other halo models may compute some other value to stay consistent with this convention (e.g., the ABGModel).
Slomo.Halos.shmr — Methodshmr(Mvir; <keyword arguments>)Compute the expected stellar mass from halo mass using the relation of Rodriguez-Puebla et al. 2017.
See equations 25-33 for the functional form and 49 - 55 for the parameters.
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensity
Slomo.Halos.shmr_prior — Methodshmr_prior(Mvir, Mstar; <keyword arguments>)Compute the probability of drawing a halo mass and stellar mass pair (Mvir, Mstar).
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensitysigma_logMstar::Real = 0.15scatter in log_Mstar (in dex)
Slomo.Halos.virial_mass — Methodvirial_mass(Rvir::Real; <keyword arguments>)Compute the virial mass from the virial radius.
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensity
Slomo.Halos.virial_mass — Methodvirial_mass(halo::HaloModel; <keyword arguments>)Compute the virial mass via root finding.
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensityxstart::Real = 10.0: factor of the scale radius to start the search for the rootrtol::Real = 1e-2: relative error tolerance for the root findingmaxevals:Int = 100: maximum number of function evalutions for the root finding
Slomo.Halos.virial_radius — Methodvirial_radius(halo::HaloModel, <keyword arguments>)Compute the virial radius of the halo via root finding.
Arguments
mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensityz::Real = 0.0: redshift at which to evaluate the overdensityxstart::Real = 10.0: factor of the scale radius to start the search for the rootrtol::Real = 1e-2: relative error tolerance for the root findingmaxevals:Int = 100: maximum number of function evalutions for the root finding
Jeans modeling
Slomo.Jeans — ModuleJeans model with line-of-sight velocity dispersion predictions.
JeansModel: collection of mass model, tracer density model, and anisotropy modelsigma_los: calculate the line-of-sight velocity dispersionkappa_los: calculate the line-of-sight velocity kurtosis
Slomo.Jeans.JeansModel — TypeJeansModel(mass_model, tracer_model, anisotropy_model)Represents the key ingredients for the Jeans modeling. The first argument can be an array of DensityModel instances, in which case we will consider the mass model to be the sum of all its components. The second argument refers to the density profile of the kinematic tracer population (e.g., the stars or star clusters).
Slomo.Jeans.kappa_los — Methodkappa_los(model::JeansModel, R; <keyword arguments>, parameters...)Calculate the line-of-sight velocity kurtosis by numerically integrating the spherical Jeans equations.
Examples
model = JeansModel(NFWModel(), SersicModel(), ConstantBetaModel())
R = collect(1:0.5:10)
kappa = kappa_los(model, R)Arguments
n_interp::Int = 10: number of points on the interpolation gridfudge::Real = 1e-6: fudge factor for tweaking integration boundsinterp::Bool = true: whether or not to interpolate from a grid of Rreturn_sigma::Bool = false: whether or not to additionally return the velocity dispersionrmax::Real = 1e2: maximum radius (in kpc) for integrating the Jeans equationparameters: keywords used to update JeansModel parameter values
This is only valid for the ConstantBetaModel anisotropy model. Using some other anisotropy model will give nonsense results.
See also: sigma_los
Slomo.Jeans.sigma_los — Methodsigma_los(model::JeansModel, R; <keyword arguments>, parameters...)Calculate the line-of-sight velocity dispersion by numerically integrating the spherical Jeans equation.
Examples
model = JeansModel(NFWModel(), SersicModel(), ConstantBetaModel())
R = collect(1:0.5:10)
sigma = sigma_los(model, R)Arguments
n_interp::Int = 10: number of points on the interpolation gridfudge::Real = 1e-6: fudge factor for tweaking integration boundsinterp::Bool = true: whether or not to interpolate from a grid of Rrmax::Real = 1e2: maximum radius (in kpc) for integrating the Jeans equationparameters: keywords used to update JeansModel parameter values
See also: sigma_los_parallel, kappa_los
Slomo.Jeans.sigma_los_parallel — Methodsigma_los_parallel(model::JeansModel, R, parameter_sets; <keyword arguments>)Calculate the line-of-sight velocity dispersion by numerically integrating the spherical Jeans equation in parallel for different sets of parameters.
parameter_sets should be a list of dictionaries where each dictionary contains an update to the model parameters.
Examples
model = JeansModel(NFWModel(), SersicModel(), ConstantBetaModel())
R = collect(1:0.5:10)
betas = collect(0.1:0.1:0.9)
parameter_sets = [Dict(:beta => b) for b in betas]
sigmas = sigma_los(model, R, parameter_sets)Arguments
n_interp::Int = 10: number of points on the interpolation gridfudge::Real = 1e-6: fudge factor for tweaking integration boundsinterp::Bool = true: whether or not to interpolate from a grid of Rrmax::Real = 1e2: maximum radius (in kpc) for integrating the Jeans equation
See also: sigma_los
- 1Pun intended