Private API

Models

Slomo.Models.has_analytic_profileMethod
has_analytic_profile(f, model::Model)

Stupid hack to see if there is an analytic profile defined for a subtype or if it needs to be numerically integrated.

source
Slomo.Models.potentialMethod
potential(model::DensityModel, r)

Gravitational potential at r. Equal to the square of the circular velocity. If not defined for a subtype of DensityModel, then calculate from the enclosed mass.

source

Anisotropy

Tracers

Slomo.Tracers.M_sersicMethod
M_sersic(r, Re, n, Mtot)

Compute the enclosed mass within deprojected radius for the Sersic model.

  • r: deprojected (3d) radii in kpc
  • Re: effective radius in kpc
  • n: index
  • Mtot: total mass
source
Slomo.Tracers.b_cbMethod
b_cb(n)

Compute the "b" parameter in the Sersic function, from the Ciotti & Bertin (1999) approximation. n is the Sersic index.

source
Slomo.Tracers.p_lnMethod
p_ln(n)

Compute the "p" parameter in fitting deprojected Sersic function, from Lima Neto et al. (1999). n is the Sersic index

source
Slomo.Tracers.rho_sersicMethod
rho_sersic(r, Re, n, Mtot)

Compute the volume density profile for Sersic model.

  • r: deprojected (3d) radii in kpc
  • Re: effective radius in kpc
  • n: index
  • Mtot: total mass
source
Slomo.Tracers.s_sersicMethod
s_sersic(R, Re, n, Mtot)

Compute the surface density profile for Sersic model.

  • R: projected (2d) radii in kpc
  • Re: effective radius in kpc
  • n: index
  • Mtot: total mass
source

Halos

Slomo.Halos.M_ABGMethod
M_ABG(r, rs, rhos, alpha, beta, gamma)

Compute the enclosed mass profile for an αβγ model.

  • r: radii in kpc
  • rs: scale radius in kpc
  • rhos: scale density in Msun / kpc3
  • alpha: transition sharpness
  • beta: negative outer log slope
  • gamma: negative inner log slope
source
Slomo.Halos.M_EinastoMethod
M_Einasto(r, rs, rhos, alpha)

Compute the enclosed mass for Einasto model.

  • r: radii in kpc
  • rs: scale radius in kpc
  • rhos: scale density in Msun / kpc3
  • alpha: Einasto shape index, lower is more cored
source
Slomo.Halos.M_GNFWMethod
M_GNFW(r, rs, rhos, gamma)

Compute the enclosed mass for gNFW model.

  • r: radii in kpc
  • rs: scale radius in kpc
  • rhos: scale density in Msun / kpc3
  • gamma: negative of inner density log slope
source
Slomo.Halos.M_solMethod
M_sol(r, rsol, rhosol)

Compute the enclosed mass for the soliton model. It's analytic!

  • rsol: soliton scale radius in kpc
  • rhosol: soliton scale density in kpc
source
Slomo.Halos.Mvir_from_RvirMethod
Mvir_from_Rvir(Rvir; <keyword arguments>)

Compute the virial mass from the virial radius, using the definition

\[M_\mathrm{vir} = 4\pi \Delta\rho / (3 R_\mathrm{vir}^3)\]

Arguments

  • mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")
  • cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensity
  • z::Real = 0.0: redshift at which to evaluate the overdensity
source
Slomo.Halos.Rvir_from_MvirMethod
Rvir_from_Mvir(Mvir; <keyword arguments>)

Compute the virial radius from the virial mass, using the definition

\[M_\mathrm{vir} = 4\pi \Delta\rho / (3 R_\mathrm{vir}^3)\]

Arguments

  • mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")
  • cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensity
  • z::Real = 0.0: redshift at which to evaluate the overdensity
source
Slomo.Halos.alpha_peakMethod
alpha_peak(Mvir; <keyword arguments>)

Compute the Einasto shape parameter from peak height using the relation of Gao+2008.

Arguments

  • mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")
  • cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensity
  • z::Real = 0.0: redshift at which to evaluate the overdensity
source
Slomo.Halos.alpha_peak_priorMethod
alpha_peak_prior(Mvir, alpha; <keyword arguments>)

Compute the probability of drawing a halo with the given virial mass and Einasto shape parameter using the alpha-peak height relation of Gao et al. 2008.

Arguments

  • mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")
  • cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensity
  • z::Real = 0.0: redshift at which to evaluate the overdensity
  • sigma0::Real = 0.16: scatter (in dex) of Mvir at redshift z = 0
  • sigmaz::Real = 0.03: slope with redshift in scatter (in dex) of Mvir
source
Slomo.Halos.drhodr_ABGMethod
drhodr_ABG(r, rs, rhos, alpha, beta, gamma)

Compute the linear density derivative for αβγ model.

  • r: radii in kpc
  • rs: scale radius in kpc
  • rhos: scale density in Msun / kpc3
  • alpha: transition sharpness
  • beta: negative outer log slope
  • gamma: negative inner log slope
source
Slomo.Halos.drhodr_NFWMethod
drhodr_NFW(r, rs, rhos)

Compute the linear density derivative for NFW model.

  • r: radii in kpc
  • rs: scale radius in kpc
  • rhos: scale density in Msun / kpc3
source
Slomo.Halos.drhodr_solMethod
drhodr_sol(r, rsol, rhosol)

Compute the linear density derivative for the soliton model.

  • rsol: soliton scale radius in kpc
  • rhosol: soliton scale density in kpc
source
Slomo.Halos.m22_from_solMethod
m22_from_sol(rsol, rhosol; <keyword arguments>)

Calculate the axion mass from the soliton scale parameters. See Marsh & Pop 2015, equation 8. Axion mass has units of 1e-22 eV.

  • rsol: soliton scale radius in kpc
  • rhosol: soliton scale density in kpc

Arguments

  • cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensity
  • z::Real = 0.0: redshift at which to evaluate the overdensity
  • alpha_mp::Real = 0.23: the alpha fitting parameter from Marsh & Pop 2015
source
Slomo.Halos.matching_radiusMethod
matching_radius(ρ1, ρ2; <keyword arguments>)

Compute the matching radius from the specified density profiles.

Arguments

  • x_start::Real = 2.0: factor of radius scale to start the search for the root
  • radius_scale::Real = 1.0: radius scale
  • maxevals:Int = 100: maximum number of function evalutions for the root finding
source
Slomo.Halos.peak_heightMethod
peak_height(Mvir; <keyword arguments>)

Computes the dimensionless peak height,

\[\nu = \delta_\mathrm{crit}(z) / \sigma(M, z)\]

using the approximation of Dutton & Maccio 2014 for a Planck 2013 cosmology.

Arguments

  • mdef::AbstractString = default_mdef: halo mass definition (e.g., "200c", "vir")
  • cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensity
  • z::Real = 0.0: redshift at which to evaluate the overdensity
source
Slomo.Halos.rho_ABGMethod
rho_ABG(r, rs, rhos, alpha, beta, gamma)

Compute the local volume density for an αβγ model.

  • r: radii in kpc
  • rs: scale radius in kpc
  • rhos: scale density in Msun / kpc3
  • alpha: transition sharpness
  • beta: negative outer log slope
  • gamma: negative inner log slope
source
Slomo.Halos.rho_EinastoMethod
rho_Einasto(r, rs, rhos, alpha)

Compute the local volume density for Einasto model.

  • r: radii in kpc
  • rs: scale radius in kpc
  • rhos: scale density in Msun / kpc3
  • alpha: Einasto shape index, lower is more cored
source
Slomo.Halos.rho_GNFWMethod
rho_GNFW(r, rs, rhos, gamma)

Compute the local volume density for NFW model.

  • r: radii in kpc
  • rs: scale radius in kpc
  • rhos: scale density in Msun / kpc3
  • gamma: negative of inner density log slope
source
Slomo.Halos.rho_NFWMethod
rho_nfw(r, rs, rhos)

Compute the local volume density for NFW model.

  • r: radii in kpc
  • rs: scale radius in kpc
  • rhos: scale density in Msun / kpc3
source
Slomo.Halos.rho_solMethod
rho_sol(r, rsol, rhosol)

Compute the local volume density for the soliton model.

  • rsol: soliton scale radius in kpc
  • rhosol: soliton scale density in kpc
source
Slomo.Halos.rhosol_from_rsolMethod
rhosol_from_rsol(m22, rsol; <keyword arguments>)

Calculate the scale density of the soliton core from m22 and the scale radius.

  • m22: ultralight axion mass in units of 1e-22 eV
  • rsol: soliton scale radius in kpc

Arguments

  • cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensity
  • z::Real = 0.0: redshift at which to evaluate the overdensity
  • alpha_mp::Real = 0.23: the alpha fitting parameter from Marsh & Pop 2015
source
Slomo.Halos.rsol_from_MvirMethod
rsol_from_Mvir(m22, Mvir)

Calculate the scale radius of the soliton core from the halo mass scaling relation (see Robles+2019).

source
Slomo.Halos.Δρ_from_mdefMethod
Δρ_from_mdef(mdef; <keyword arguments>)

Parse the halo mass definition to get the spherical overdensity in Msun / kpc3. If mdef is "vir", then use the spherical overdensity relation of Bryan & Norman 1998.

mdef should be either "vir" or of the form "[0-9]*[c|m]" where "c" refers to a multiple of the critical density and "m" refers to a multiple of the mass density, e.g.,

mdef in ["vir", "200c", "500c", "200m", "500m", ...]

Arguments

  • cosmo::AbstractCosmology = default_cosmo: cosmology under which to evaluate the overdensity
  • z::Real = 0.0: redshift at which to evaluate the overdensity
source

Jeans modeling

Slomo.Jeans.integral_from_kernelMethod
integral_from_kernel(Rgrid, M, ρ, K; <keyword arguments>)

Compute the outer Jeans integral for a model which has an analytic form for the Jeans projection kernel.

source
Slomo.Jeans.integral_from_vr2Method
integral_from_vr2(Rgrid, M, ρ, K; <keyword arguments>)

Compute the double Jeans integral for a model which has no analytic form for the Jeans projection kernel.

source

Cosmology

Slomo.CosmologyToolsModule

Copy of Cosmology.jl to remove Unitful dependence.

Original license:

The MIT License (MIT)

Copyright (c) 2013 Mike Nolta <mike@nolta.net>

Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:

The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.

THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

source
Slomo.CosmologyTools.cosmologyMethod
cosmology(;h = 0.69,
           Neff = 3.04,
           OmegaK = 0,
           OmegaM = 0.29,
           OmegaR = nothing,
           Tcmb = 2.7255,
           w0 = -1,
           wa = 0)

Parameters

  • h - Dimensionless Hubble constant
  • OmegaK - Curvature density (Ω_k)
  • OmegaM - Matter density (Ω_m)
  • OmegaR - Radiation density (Ω_r)
  • Tcmb - CMB temperature in Kelvin; used to compute Ω_γ
  • Neff - Effective number of massless neutrino species; used to compute Ω_ν
  • w0 - CPL dark energy equation of state; w = w0 + wa(1-a)
  • wa - CPL dark energy equation of state; w = w0 + wa(1-a)

Examples

julia> c = cosmology()
Cosmology.FlatLCDM{Float64}(0.69, 0.7099122024007928, 0.29, 8.77975992071536e-5)
julia> c = cosmology(OmegaK=0.1)
Cosmology.OpenLCDM{Float64}(0.69, 0.1, 0.6099122024007929, 0.29, 8.77975992071536e-5)
julia> c = cosmology(w0=-0.9, OmegaK=-0.1)
Cosmology.ClosedWCDM{Float64}(0.69, -0.1, 0.8099122024007929, 0.29, 8.77975992071536e-5, -0.9, 0.0)
source

Integrate

Utils

Constants