Helper Functions
Distances and Sizes
StarFormationHistories.arcsec_to_pc — Functionarcsec_to_pc(arcsec, dist_mod)Converts on-sky angle in arcseconds to physical separation based on distance modulus under the small-angle approximation.
\[r ≈ 10^{μ/5 + 1} \times \text{atan}(θ/3600)\]
StarFormationHistories.pc_to_arcsec — Functionpc_to_arcsec(pc, dist_mod)Inverse of arcsec_to_pc.
\[θ ≈ \text{tan}\left( r / 10^{μ/5 + 1} \right) \times 3600\]
StarFormationHistories.distance_modulus — Functiondistance_modulus(distance)Finds distance modulus for distance in parsecs.
\[μ = 5 \times \log_{10}(d) - 5\]
StarFormationHistories.distance_modulus_to_distance — Functiondistance_modulus_to_distance(dist_mod)Converts distance modulus to distance in parsecs.
\[d = 10^{μ/5 + 1}\]
Magnitudes and Luminosities
StarFormationHistories.mag2flux — Functionmag2flux(m, zpt=0)Convert a magnitude m to a flux assuming a photometric zeropoint of zpt, defined as the magnitude of an object that produces one count (or data number, DN) per second.
julia> mag2flux(15.0, 25.0) ≈ exp10(4 * (25.0 - 15.0) / 10)
trueStarFormationHistories.flux2mag — Functionflux2mag(f, zpt=0)Convert a flux f to a magnitude assuming a photometric zeropoint of zpt, defined as the magnitude of an object that produces one count (or data number, DN) per second.
julia> flux2mag(10000.0, 25.0) ≈ 25.0 - 5 * log10(10000.0) / 2
trueStarFormationHistories.magerr — Functionmagerr(f, σf)Returns an error in magnitudes given a flux and a flux uncertainty.
julia> magerr(100.0, 1.0) ≈ 2.5 / log(10) * (1.0 / 100.0)
trueStarFormationHistories.fluxerr — Functionfluxerr(f, σm)Returns an error in flux given a flux and a magnitude uncertainty.
julia> fluxerr(100.0, 0.01) ≈ (0.01 * 100.0) / 2.5 * log(10)
trueStarFormationHistories.snr_magerr — Functionsnr_magerr(σm)Returns a signal-to-noise ratio $(f/σf)$ given an uncertainty in magnitudes.
julia> snr_magerr(0.01) ≈ 2.5 / log(10) / 0.01
trueStarFormationHistories.magerr_snr — Functionmagerr_snr(snr)Returns a magnitude uncertainty given a signal-to-noise ratio $(f/σf)$.
julia> magerr_snr(100.0) ≈ 2.5 / log(10) / 100.0
trueMetallicities
StarFormationHistories.Y_from_Z — FunctionY_from_Z(Z, Y_p=0.2485, γ=1.78)Calculates the helium mass fraction (Y) for a star given its metal mass fraction (Z) using the approximation Y = Y_p + γ * Z, with Y_p being the primordial helium abundance Y_p=0.2485 as assumed for PARSEC isochrones and γ=1.78 matching the input scaling for PARSEC as well.
StarFormationHistories.X_from_Z — FunctionX_from_Z(Z[, Yp, γ])Calculates the hydrogen mass fraction (X) for a star given its metal mass fraction (Z) via X = 1 - (Z + Y), with the helium mass fraction Y approximated via StarFormationHistories.Y_from_Z. You may also provide the primordial helium abundance Y_p and γ such that Y = Y_p + γ * Z; these are passed through to X_from_Z.
StarFormationHistories.MH_from_Z — FunctionMH_from_Z(Z, solZ=0.01524; Y_p = 0.2485, γ = 1.78)Calculates [M/H] = log(Z/X) - log(Z/X)⊙. Given the provided solar metal mass fraction solZ, it calculates the hydrogen mass fraction X for both the Sun and the provided Z with StarFormationHistories.X_from_Z. You may also provide the primordial helium abundance Y_p and γ such that Y = Y_p + γ * Z; these are passed through to X_from_Z.
PARSEC Normalization
For PARSEC tracks, [M/H]=0 will result in isochrones with Z=0.014711, instead of the Z=0.01524 expected for the present Sun (Caffau et al. 2011). This small offset is caused by the intrinsic difference between the helium content in the model for the present Sun, and the initial helium content for which the PARSEC evolutionary tracks are computed. That is:
- The model for the present Sun is forced to reproduce the observed ratio of Z⊙/X⊙=0.0207. This ratio defines the zero-point of the [M/H] scale, [M/H] = log(Z/X)-log(Z⊙/X⊙).
- The same solar model has an initial surface composition of Zinitial=0.01774, Yinitial=0.28 (see Table 3 in Bressan et al. 2012). Together with the primordial helium content of Yp=0.2485, these numbers define the relation Y = Yp + 1.78Z, which was used to build the set of initial (Z, Y) values for which we compute the tracks.
- Therefore, the Y(Z) relation used to compute the grids, is respected by the initial Sun, but not by the present Sun (which is affected by diffusion). This creates a 0.015 dex offset between the metallicity scale of the present Sun, and the initial metallicity scale of the computed grids of tracks. See Table 4 in Bressan et al. (2012) for a more complete list of [M/H] values that follow from this approximation.
The above note was taken from the CMD webform FAQ, copyright Leo Girardi.
This function is an approximation and may not be suitable for precision calculations.
StarFormationHistories.dMH_dZ — FunctiondMH_dZ(Z, solZ=0.01524; Y_p = 0.2485, γ = 1.78)Partial derivative of MH_from_Z with respect to the input metal mass fraction Z. Used for LogarithmicAMR.
StarFormationHistories.Z_from_MH — FunctionZ_from_MH(MH, solZ=0.01524; Y_p = 0.2485, γ = 1.78)Calculates metal mass fraction Z assuming that the solar metal mass fraction is solZ and using the PARSEC relation for the helium mass fraction Y = Y_p + γ * Z with primordial helium abundance Y_p = 0.2485, and γ = 1.78.
StarFormationHistories.dZ_dMH — FunctiondZ_dMH(MH, solZ=0.01524; Y_p = 0.2485, γ = 1.78)Partial derivative of Z_from_MH with respect to the input metal abundance MH. Used for LogarithmicAMR.
StarFormationHistories.mdf_amr — Function(unique_MH, mass_mdf) =
mdf_amr(coeffs::AbstractVector{<:Number},
logAge::AbstractVector{<:Number},
metallicities::AbstractVector{<:Number})Calculates the mass-weighted metallicity distribution function given a set of stellar mass coefficients coeffs for stellar populations with logarithmic ages logAge=log10(age [yr]) and metallicities given by metallicities. This is calculated as
\[P_j = \frac{ \sum_k r_{j,k} \, [\text{M} / \text{H}]_k}{\sum_{j,k} r_{j,k} \, [\text{M} / \text{H}]_k}\]
where $r_{j,k}$ are the elements of coeffs where $j$ indexes over unique entries in logAge and $k$ indexes over unique entries in metallicities. This is the same nomenclature used in the the documentation on constrained metallicity evolutions. The return values are sorted so that unique_MH is in increasing order.
Examples
julia> mdf_amr([1.0, 2.0, 1.0], [10, 10, 10], [-2, -1.5, -1])
([-2.0, -1.5, -1.0], [0.25, 0.5, 0.25])(unique_MH, mass_mdf) =
mdf_amr(coeffs::AbstractVector{<:Number},
logAge::AbstractVector{<:Number},
metallicities::AbstractVector{<:Number},
models::Union{AbstractVector{<:AbstractMatrix{<:Number}},
AbstractMatrix{<:Number}})Calculates the number-weighted metallicity distribution function given a set of coefficients coeffs for stellar populations with logarithmic ages logAge=log10(age [yr]), metallicities given by metallicities, and Hess diagram templates models. This function constructs a composite Hess diagram (see composite! for definition) out of the coeffs and models that match each unique entry in metallicites. The sums over the bins of these composite Hess diagrams then give the total predicted number of stars in the Hess diagram viewport for each metallicity. The raw number counts for each unique metallicity are returned as mass_mdf – these can be renormalized afterward to sum to one by calculating mass_mdf ./ sum(mass_mdf).
Examples
julia> mdf_amr([1.0, 2.0, 1.0], [10, 10, 10], [-2, -1.5, -1],
map(Base.Fix2(fill, (100, 100)), [1.0, 2.0, 1.0]))
([-2.0, -1.5, -1.0], [10000.0, 40000.0, 10000.0])Photometric Error and Completeness Models
StarFormationHistories.Martin2016_complete — Functionη(m) = Martin2016_complete(m, A, m50, ρ)Completeness model of Martin et al. 2016 implemented as their Equation 7:
\[\eta(m) = \frac{A}{1 + \text{exp} \left( \frac{m - m_{50}}{\rho} \right)}\]
m is the magnitude of interest, A is the maximum completeness, m50 is the magnitude at which the data are 50% complete, and ρ is an effective slope modifier.
StarFormationHistories.exp_photerr — Functionexp_photerr(m, a, b, c, d)Exponential model for photometric errors of the form
\[\sigma(m) = a^{b \times \left( m-c \right)} + d\]
Reported values for some HST data were a=1.05, b=10.0, c=32.0, d=0.01.
StarFormationHistories.process_ASTs — Functionprocess_ASTs(ASTs,
inmag::Symbol,
outmag::Symbol,
bins::AbstractVector{<:Real},
selectfunc;
statistic=StatsBase.median)Processes a table of artificial stars to calculate photometric completeness, bias, and error across the provided bins. This method has no default implementation and is implemented in package extensions that rely on either DataFrames.jl or TypedTables.jl being loaded into your Julia session to load the relevant method. This method therefore requires Julia 1.9 or greater to use.
Arguments
ASTsis the table of artificial stars to be analyzed. This function should support most table types that follow the Tables.jl interface, but we specifically test to ensureDataFrames.DataFrameandTypedTables.Tablewill work.inmagis the column name in symbol format (e.g., :F606Wi) that corresponds to the intrinsic (input) magnitudes of the artificial stars.outmagis the column name in symbol format (e.g., :F606Wo) that corresponds to the measured (output) magnitude of the artificial stars.binsgive the bin edges to be used when computing the binned statistics.selectfuncis a method that takes a single row fromASTs, corresponding to a single artificial star, and returns a boolean that istrueif the star is considered successfully measured.
Keyword Arguments
statisticis the method that will be used to determine the bias and error, i.e.,bias = statistic(out .- in)anderror = statistic(abs.(out .- in)). By default we useStatsBase.median, but you could instead use a simple or sigma-clipped mean if so desired.
Returns
This method returns a result of type NTuple{4,Vector{Float64}}. Each vector is of length length(bins)-1. result contains the following elements, each of which are computed over the provided bins considering only artificial stars for which selectfunc returned true:
result[1]contains the mean input magnitude of the stars in each bin.result[2]contains the completeness value measured for each bin, defined as the fraction of input stars in each bin for whichselectfuncreturnedtrue.result[3]contains the photometric bias measured for each bin, defined asstatistic(out .- in), whereoutare the measured (output) magnitudes andinare the intrinsic (input) magnitudes.result[4]contains the photometric error measured for each bin, with the effect of photometric bias removed. This is defined asstatistic(abs.(out .- in .- bias)), withoutandindefined as above, andbiasdefined as inresult[3]above. This choice is made to simplify the modeling of the probability distribution for the observed magnitude given an intrinsic magnitude, which we may write asP(O|I). For Gaussian errors, this may be written asP(O|I) = Normal(I + bias, σ)wherebiasis the photometric bias as returned byresult[3]. This definition requires the error σ to be measured with the effect of photometric bias removed, as we implement here.
Examples
Let
F606Wibe a vector containing the input magnitudes of your artificial starsF606Wobe a vector containing the measured magnitudes of the artificial stars, where a value of 99.999 indicates a non-detection.flagbe a vector of booleans that indicates whether the artificial star passed additional quality cuts (star-galaxy separation, etc.)
You could call this method as
import TypedTables: Table
process_ASTs(Table(input=F606Wi, output=F606Wo, good=flag),
:input, :output, minimum(F606Wi):0.1:maximum(F606Wi),
x -> (x.good==true) & (x.output != 99.999))See also the tests in test/utilities/process_ASTs_test.jl.
Miscellaneous
StarFormationHistories.integrate_sfr — Functionintegrate_sfr(la_l, la_h, sfr)Returns the total integrated stellar mass given log(age [yr]) bins and the corresponding star formation rates in solar masses per year. la_l and la_h must be the lower and upper limits of log(age [yr]) bins and sfr are the corresponding star formation rates in solar masses per yr sfr. Assumes constant SFR in each time bin.
julia> using StarFormationHistories: integrate_sfr
julia> la_l = 6.6:0.1:10.0; # Lower log(age [yr]) bins
julia> la_h = vcat(la_l[2:end], 10.13); # Upper log(age [yr]) bins
julia> integrate_sfr(la_l, la_h, fill(1e-4, length(la_l))) ≈ (exp10(10.13) - exp10(6.6)) * 1e-4
true