Index
XGPaint.AbstractForegroundModelXGPaint.CIB_Planck2013XGPaint.CO_CROWNEDXGPaint.LRG_Yuan23XGPaint.LogInterpolatorProfileXGPaint.Radio_Sehgal2009XGPaint.P_e_losXGPaint.P_th_losXGPaint.SZpackXGPaint.SZpack_rkszXGPaint.T_mass_calcXGPaint.T_vir_calcXGPaint.bl2beamXGPaint.build_c_lnm2r_interpolatorXGPaint.build_interpolatorXGPaint.build_jiang_interpolatorXGPaint.build_r2z_interpolatorXGPaint.build_shang_interpolatorXGPaint.build_z2r_interpolatorXGPaint.build_zhengcen_interpolatorXGPaint.build_zhengsat_interpolatorXGPaint.chunkXGPaint.cleanup_negatives!XGPaint.compton_yXGPaint.ellpadXGPaint.fquench_UMXGPaint.generate_sourcesXGPaint.generate_sourcesXGPaint.generate_sourcesXGPaint.generate_subhalo_offsetsXGPaint.get_anglesXGPaint.get_basic_halo_propertiesXGPaint.get_cosmologyXGPaint.get_interpolatorsXGPaint.get_interpolatorsXGPaint.get_ring_disc_rangesXGPaint.hod_sehgalXGPaint.l2fXGPaint.load_example_halosXGPaint.load_precomputed_battagliaXGPaint.m2rXGPaint.mz2cXGPaint.paint!XGPaint.paint!XGPaint.paint!XGPaint.rSZ_perturbativeXGPaint.ra_dec_redshift_to_xyzXGPaint.read_halo_catalog_hdf5XGPaint.realspacegaussbeamXGPaint.rsz_szpack_table_filenameXGPaint.sehgal_LF!XGPaint.transform_profile_grid!XGPaint.vectorhealpixmapXGPaint.z_evo
XGPaint.AbstractForegroundModel — TypeAll foreground models inherit from this type.
XGPaint.CIB_Planck2013 — TypeCIB_Planck2013{T}(; kwargs...)Define CIB model parameters. Defaults are from Viero et al. 2013. All numbers not typed are converted to type T. This model has the following parameters and default values:
nside::Int64 = 4096min_redshift = 0.0max_redshift = 5.0min_mass = 1e12box_size = 40000shang_zplat = 2.0shang_Td = 20.7shang_betan = 1.6shang_eta = 2.4shang_alpha = 0.2shang_Mpeak = 10^12.3shang_sigmaM = 0.3shang_Msmin = 1e11shang_Mmin = 1e10shang_I0 = 92jiang_gamma_1 = 0.13jiang_alpha_1 = -0.83jiang_gamma_2 = 1.33jiang_alpha_2 = -0.02jiang_beta_2 = 5.67jiang_zeta = 1.19
XGPaint.CO_CROWNED — TypeCO_CROWNED{T}(; kwargs...)Define broadband CO model parameters. All numbers not typed are converted to type T. This model has the following parameters and default values:
nside::Int64 = 4096hod::String = "shang"
XGPaint.LRG_Yuan23 — TypeLRG_Yuan23{T}(; kwargs...)Define LRG model parameters. Defaults from Yuan et al. 2023 [arXiv:2306.06314]. All numbers not typed are converted to type T. This model has the following parameters and default values:
nside::Int64 = 4096hod::String = "zheng"min_redshift = 0.0max_redshift = 5.0min_mass = 1e12box_size = 40000shang_Msmin = 1e11zheng_Mcut = 10^12.7zheng_M1 = 10^13.6zheng_sigma = 0.2zheng_kappa = 0.08zheng_alpha = 1.15yuan_ic = 0.8jiang_gamma_1 = 0.13jiang_alpha_1 = -0.83jiang_gamma_2 = 1.33jiang_alpha_2 = -0.02jiang_beta_2 = 5.67jiang_zeta = 1.19
XGPaint.LogInterpolatorProfile — TypeLogInterpolatorProfile{T, P, I1}A profile that interpolates over a positive-definite function (θ, z, M_halo), but internally interpolates over log(θ) and log10(M) using a given interpolator. Evaluation of this profile is then done by exponentiating the result of the interpolator.
f(θ, z, M) = exp(itp(log(θ), z, log10(M)))This is useful for interpolating over a large range of scales and masses, where the profile is expected to be smooth in log-log space. It wraps the original model and also the interpolator object itself.
XGPaint.Radio_Sehgal2009 — TypeRadio_Sehgal2009{T}(model parameters...)Define CIB model parameters. Defaults are from Viero et al. 2013.
model = CIBModel{Float32}(a_0=0.4)XGPaint.P_e_los — MethodLine-of-sight integrated electron pressure
XGPaint.P_th_los — MethodLine-of-sight integrated thermal pressure
XGPaint.SZpack — MethodSZpack(model_szp, θ, M_200, z; showT=false)Outputs the integrated compton-y signal calculated using SZpack along the line of sight. Note: M_200 requires units.
XGPaint.SZpack_rksz — MethodSZpack_rksz(model, r, M_200, z, vel, mu)Computes the relativistic kSZ signal using the SZpack tables.
XGPaint.T_mass_calc — MethodT_mass_calc(model,M,z::T; scale_type="Ty", sim_type="combination") where TCalculates the temperature for a given halo using https://arxiv.org/pdf/2207.05834.pdf.
XGPaint.T_vir_calc — MethodT_vir_calc(model, M, z)Calculates the virial temperature for a given halo using Wang et al. 2007.
XGPaint.bl2beam — MethodCompute the real-space beam from a harmonic-space beam.
XGPaint.build_c_lnm2r_interpolator — Methodbuild_c_lnm2r_interpolator(;cmin::T=1f-3, cmax::T=25.0f0,
mmin::T=-7.1f0, mmax::T=0.0f0, nbin=100) where TGenerate a LinearInterpolation object that turns concentration and ln(M_halo) into satellite radius.
XGPaint.build_interpolator — MethodHelper function to build a (θ, z, Mh) interpolator
XGPaint.build_jiang_interpolator — MethodBuild a linear interpolation function which maps log(Mh) to Nsh.
XGPaint.build_r2z_interpolator — Methodbuild_r2z_interpolator(min_z::T, max_z::T,
cosmo::Cosmology.AbstractCosmology; n_bins=2000) where TConstruct a fast r2z linear interpolator.
XGPaint.build_shang_interpolator — MethodBuild a linear interpolation function which maps log(Mh) to Nsat.
XGPaint.build_z2r_interpolator — MethodConstruct a z2r linear interpolator.
XGPaint.build_zhengcen_interpolator — MethodBuild a linear interpolation function which maps log(Mh) to Ncen.
XGPaint.build_zhengsat_interpolator — MethodBuild a linear interpolation function which maps log(Mh) to Nsat.
XGPaint.chunk — MethodGenerates a list of tuples which describe starting and ending chunk indices. Useful for parallelizing an array operation.
XGPaint.cleanup_negatives! — Methodprune a profile grid for negative values, extrapolate instead
XGPaint.compton_y — Methodcompton_y(model, r, M_200c, z)Calculate the Compton y parameter for a given model at a given radius, mass, and redshift. Mass needs to have units!
XGPaint.ellpad — MethodUtility function which prepends some zeros to an array. It makes a copy instead of modifying the input.
XGPaint.fquench_UM — MethodQuiescent fraction recipe from UniverseMachine
XGPaint.generate_sources — Methodgenerate_sources(model, cosmo, halo_pos_inp, halo_mass_inp; verbose=true)Produce a source catalog from a model and halo catalog. This converts the halo arrays into the type specified by model.
Arguments:
model::AbstractCIBModel{T}: source model parameterscosmo::Cosmology.FlatLCDM{T}: background cosmologyHealpix_res::Resolution: Healpix map resolutionhalo_pos_inp::AbstractArray{TH,2}: halo positions with dims (3, nhalos)halo_mass_inp::AbstractArray{TH,1}: halo masses
Keywords
verbose::Bool=true: print out progress details
XGPaint.generate_sources — Methodgenerate_sources(model, cosmo, halo_pos_inp, halo_mass_inp; verbose=true)Produce a source catalog from a model and halo catalog. This converts the halo arrays into the type specified by model.
Arguments:
model::AbstractCIBModel{T}: source model parameterscosmo::Cosmology.FlatLCDM{T}: background cosmologyHealpix_res::Resolution: Healpix map resolutionhalo_pos_inp::AbstractArray{TH,2}: halo positions with dims (3, nhalos)halo_mass_inp::AbstractArray{TH,1}: halo masses
Keywords
verbose::Bool=true: print out progress details
XGPaint.generate_sources — MethodProduce a source catalog from a model and halo catalog.
XGPaint.generate_subhalo_offsets — MethodGenerate an array where the value at index i corresponds to the index of the first source of halo i. Takes an array where the value at index i corresponds to the number of subhalos that halo i has.
XGPaint.get_angles — MethodCompute angles of halos
XGPaint.get_basic_halo_properties — MethodFill in basic halo properties.
XGPaint.get_cosmology — Methodget_cosmology(::Type{T}; h=0.69, Neff=3.04, OmegaK=0.0,
OmegaM=0.29, OmegaR=nothing, Tcmb=2.7255, w0=-1, wa=0)Construct a background cosmology. This function duplicates the cosmology() function in Cosmology.jl, but with the numeric type specified.
Arguments:
::Type{T}: numerical type to use for calculations
Keywords
h- Dimensionless Hubble constantOmegaK- 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)
Example
julia> get_cosmology(Float32; h=0.7)
Cosmology.FlatLCDM{Float32}(0.7f0, 0.7099147f0, 0.29f0, 8.5307016f-5)XGPaint.get_interpolators — MethodConstruct the necessary interpolator set.
XGPaint.get_interpolators — MethodConstruct the necessary interpolator set.
XGPaint.get_ring_disc_ranges — Methodget_ring_disc_ranges(workspace, ring_idx, θ_center, ϕ_center, radius)Returns (range1, range2) of pixel indices on a ring that lie within the disc.
The spherical distance between points (θ₁,φ₁) and (θ₂,φ₂) is given by: cos(d) = cos(θ₁)cos(θ₂) + sin(θ₁)sin(θ₂)cos(φ₁ - φ₂)
For a point to be included in the disc: d ≤ radius.
Returns two ranges to handle φ wraparound at the 0/2π boundary:
- range1: Main contiguous range of pixels
- range2: Additional range when disc crosses φ=0/2π seam (empty if no crossing)
XGPaint.hod_sehgal — MethodPopulate halos with radio sources according to the HOD in Sehgal et al. 2009.
The optional rng parameter provides an array of random number generators, one for each thread.
XGPaint.l2f — MethodInverse square law with redshift dependence.
XGPaint.load_example_halos — MethodReads a collection of example halos out of the WebSky halo catalogs, and returns RA (rad), DEC (rad), redshift, and halo mass (M200c).
XGPaint.load_precomputed_battaglia — MethodReads in a standard Battaglia 2016 model from disk for tSZ.
XGPaint.m2r — Methodm2r(m::T, cosmo::Cosmology.FlatLCDM{T}) where TConvert virial mass to virial radius.
XGPaint.mz2c — Methodmz2c(m::T, z::T, cosmo::Cosmology.FlatLCDM{T}) where TCompute concentration factor from Duffy et al. 2008.
XGPaint.paint! — MethodPaint a source catalog onto a map.
This function creates the arrays for you.
XGPaint.paint! — Methodpaint!(result_map, nu_obs, model, sources, fluxes_cen, fluxes_sat)Paint a source catalog onto a map, recording the fluxes in fluxes_cen and fluxes_sat.
Arguments:
result_map::HealpixMap{T_map, RingOrder}: Healpix map to paintnu_obs: frequency in Hzmodel::AbstractCIBModel{T}: source model parameterssources: NamedTuple containing source information from generate_sourcesfluxes_cen::AbstractArray: buffer for writing fluxes of centralsfluxes_sat::AbstractArray: buffer for writing fluxes of satellites
XGPaint.paint! — Methodpaint!(result_map, model, sources, min_redshift, max_redshift)'Paint' the LRG catalog onto a map.
Arguments:
result_map::HealpixMap{T_map, RingOrder}: Healpix map to paintmodel::AbstractCIBModel{T}: source model parameterssources: NamedTuple containing source information from generate_sourcesfluxes_cen::AbstractArray: buffer for writing fluxes of centralsfluxes_sat::AbstractArray: buffer for writing fluxes of satellites
XGPaint.rSZ_perturbative — MethodrSZ_perturbative(model, r, M_200, z; T_scale="virial", sim_type="combination", showT=true)Calculates the integrated relativistic compton-y signal along the line of sight. Mass requires units.
XGPaint.ra_dec_redshift_to_xyz — MethodConvert RA (rad), DEC (rad), and redshift to xyz comoving radial dist.
XGPaint.read_halo_catalog_hdf5 — MethodUtility function to read an HDF5 table with x, y, z, M_h as the four rows. The hdf5 record is "halos".
XGPaint.realspacegaussbeam — MethodComputes a real-space beam interpolator and a maximum
XGPaint.rsz_szpack_table_filename — MethodReturns a stored artifact: a precomputed SZpack table
XGPaint.sehgal_LF! — MethodFills the result array with draws from the luminosity function.
XGPaint.transform_profile_grid! — MethodApply a beam to a profile grid
XGPaint.vectorhealpixmap — Methodvecmap([T=Float64], nside::Int)Generate a map of 3-tuples, showing the (x,y,z) values of the points on the unit 2-sphere for a Healpix map.
XGPaint.z_evo — MethodCompute redshift evolution factor for LF.