Helper Functions

Distances and Sizes

StarFormationHistories.arcsec_to_pcFunction
arcsec_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)\]

source

Magnitudes and Luminosities

StarFormationHistories.mag2fluxFunction
mag2flux(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)
true
source
StarFormationHistories.flux2magFunction
flux2mag(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
true
source
StarFormationHistories.magerrFunction
magerr(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)
true
source
StarFormationHistories.fluxerrFunction
fluxerr(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)
true
source
StarFormationHistories.snr_magerrFunction
snr_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
true
source

Metallicities

StarFormationHistories.Y_from_ZFunction
Y_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.

source
StarFormationHistories.X_from_ZFunction
X_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.

source
StarFormationHistories.MH_from_ZFunction
MH_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.

The present-day solar Z is measured to be 0.01524 (Caffau et al. 2011), but for PARSEC isochrones an [M/H] of 0 corresponds to Z=0.01471. This is because of a difference between the Sun's initial and present helium content caused by diffusion. If you want to reproduce PARSEC's scaling, you should set solZ=0.01471.

This function is an approximation and may not be suitable for precision calculations.

source
StarFormationHistories.Z_from_MHFunction
Z_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.

source
StarFormationHistories.mdf_amrFunction
(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])
source
(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])
source

Photometric Error and Completeness Models

StarFormationHistories.Martin2016_completeFunction
η(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.

source
StarFormationHistories.exp_photerrFunction
exp_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.

source
StarFormationHistories.process_ASTsFunction
process_ASTs(ASTs::Union{DataFrames.DataFrame,
                         TypedTables.Table},
             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

  • ASTs is the table of artificial stars to be analyzed.
  • inmag is the column name in symbol format (e.g., :F606Wi) that corresponds to the intrinsic (input) magnitudes of the artificial stars.
  • outmag is the column name in symbol format (e.g., :F606Wo) that corresponds to the measured (output) magnitude of the artificial stars.
  • bins give the bin edges to be used when computing the binned statistics.
  • selectfunc is a method that takes a single row from ASTs, corresponding to a single artificial star, and returns a boolean that is true if the star is considered successfully measured.

Keyword Arguments

  • statistic is the method that will be used to determine the bias and error, i.e., bias = statistic(out .- in) and error = statistic(abs.(out .- in)). By default we use StatsBase.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 which selectfunc returned true.
  • result[3] contains the photometric bias measured for each bin, defined as statistic(out .- in), where out are the measured (output) magnitudes and in are the intrinsic (input) magnitudes.
  • result[4] contains the photometric error measured for each bin, defined as statistic(abs.(out .- in)), with out and in defined as above.

Examples

Let

  • F606Wi be a vector containing the input magnitudes of your artificial stars
  • F606Wo be a vector containing the measured magnitudes of the artificial stars, where a value of 99.999 indicates a non-detection.
  • flag be 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.

source