PhotometricFilters
This package provides access to, and operations on, photometric filter curves. Such filter curves are defined by a filter's transmission as a function of wavelength. Transmission and wavelength vectors are therefore the foundation of a filter curve, but it is also important to note whether the filter is used for photon counter or energy counter detectors, as the integrals used to calculate statistics over a filter curve are different between these two types of detectors.
Types
All photometric filter types should be subtypes of the AbstractFilter type. We define a minimal API that should be followed so that new types can make use of the generic filter operations we define.
The simplest concrete filter type we provide to represent photometric filters is PhotometricFilter.
PhotometricFilters.PhotometricFilter — TypePhotometricFilter(wavelength::AbstractVector, throughput::AbstractVector{T};
detector::DetectorType=Photon(), filtername::Union{String, Nothing}=nothing) where TStruct representing a photometric filter, defined by vectors of wavelengths (wavelength) and filter throughputs (throughput). wavelength can have Unitful units attached, otherwise they are assumed to be Å. Optional keyword arguments define the detector type for which the filter is valid and a name to identify the filter.
julia> using PhotometricFilters: PhotometricFilter, Photon, wavelength, throughput
julia> using Unitful
julia> f = PhotometricFilter(1000:2000, vcat(fill(0.25, 250), fill(0.5, 500), fill(0.25, 251))) # Specify only wavelength and throughput
1001-element PhotometricFilter{Float64}: nothing
reference wave.: 1478.1028279485677 Å
min. wave.: 1000 Å
max. wave.: 2000 Å
effective wave.: 1602.7669435459648 Å
mean wave.: 1499.8333333333333 Å
central wave.: 1499.5 Å
pivot wave.: 1478.1028279485677 Å
eff. width: 750.0 Å
fwhm: 501.0 Å
julia> f == PhotometricFilter(uconvert.(Unitful.nm, wavelength(f)), throughput(f)) # Can also specify wavelength argument with Unitful units
true
julia> f[10] # Indexing into the filter as `f[i]` returns `(wavelength(f)[i], throughput(f)[i])`
(1009 Å, 0.25)
julia> f(1001.1) # Calling `f` like a function interpolates the throughput
0.25
julia> f(100.11 * Unitful.nm) # Can also specify wavelength with units
0.25The data contained in the struct can be accessed with the following methods:
PhotometricFilters.filtername — Functionfiltername(f::AbstractFilter)Returns a string indicating a human-readable name for the filter (e.g., "SDSS_u").
julia> using PhotometricFilters: SDSS_u, filtername
julia> filtername(SDSS_u())
"SDSS_u"PhotometricFilters.wavelength — Functionwavelength(f::AbstractFilter)Returns the wavelength vector of the filter transmission curve with proper Unitful.jl units.
julia> using PhotometricFilters: SDSS_u, wavelength
julia> using Unitful: Quantity
julia> wavelength(SDSS_u()) isa Vector{<:Quantity}
truePhotometricFilters.throughput — Functionthroughput(f::AbstractFilter)Returns the throughput vector of the filter transmission curve (no units).
julia> using PhotometricFilters: SDSS_u, throughput
julia> throughput(SDSS_u()) isa Vector{<:Number}
truePhotometricFilters.detector_type — Functiondetector_type(f::AbstractFilter)Return an instance of PhotometricFilters.Energy if the filter is defined for energy-counting detectors or PhotometricFilters.Photon for photon-counting detectors.
julia> using PhotometricFilters: SDSS_u, detector_type, Photon
julia> detector_type(SDSS_u()) === Photon()
trueUsers can construct their own filter curvers from raw data using this type.
Accessing Filter Curves
We provide a modest collection of filter curves through a data dependency. The available filter curves are accessible via the FILTER_NAMES module constant,
using PhotometricFilters
PhotometricFilters.FILTER_NAMES |> println["2MASS_H", "2MASS_J", "2MASS_Ks", "CFHT_CFH12K_CFH7406", "CFHT_CFH12K_CFH7504", "CFHT_MEGAPRIME_CFH7605", "CFHT_MEGAPRIME_CFH7701", "CFHT_MEGAPRIME_CFH7803", "CFHT_MEGAPRIME_CFH9301", "CFHT_MEGAPRIME_CFH9401", "CFHT_MEGAPRIME_CFH9601", "CFHT_MEGAPRIME_CFH9701", "CFHT_MEGAPRIME_CFH9702", "CFHT_MEGAPRIME_CFH9801", "CFHT_WIRCAM_CFH8002", "CFHT_WIRCAM_CFH8101", "CFHT_WIRCAM_CFH8102", "CFHT_WIRCAM_CFH8103", "CFHT_WIRCAM_CFH8104", "CFHT_WIRCAM_CFH8201", "CFHT_WIRCAM_CFH8202", "CFHT_WIRCAM_CFH8203", "CFHT_WIRCAM_CFH8204", "CFHT_WIRCAM_CFH8301", "CFHT_WIRCAM_CFH8302", "CFHT_WIRCAM_CFH8303", "CFHT_WIRCAM_CFH8304", "CFHT_WIRCAM_CFH8305", "GALEX_FUV", "GALEX_NUV", "GROUND_BESSELL_H", "GROUND_BESSELL_J", "GROUND_BESSELL_K", "GROUND_COUSINS_I", "GROUND_COUSINS_R", "GROUND_JOHNSON_B", "GROUND_JOHNSON_U", "GROUND_JOHNSON_V", "GaiaDR2_BP", "GaiaDR2_G", "GaiaDR2_RP", "GaiaDR2_weiler_BPbright", "GaiaDR2_weiler_BPfaint", "GaiaDR2_weiler_G", "GaiaDR2_weiler_RP", "GaiaDR2v2_BP", "GaiaDR2v2_G", "GaiaDR2v2_RP", "Gaia_BP", "Gaia_G", "Gaia_MAW_BP_bright", "Gaia_MAW_BP_faint", "Gaia_MAW_G", "Gaia_MAW_RP", "Gaia_RP", "Gaia_rvs", "HERSCHEL_PACS_BLUE", "HERSCHEL_PACS_GREEN", "HERSCHEL_PACS_RED", "HERSCHEL_SPIRE_PLW", "HERSCHEL_SPIRE_PLW_EXT", "HERSCHEL_SPIRE_PMW", "HERSCHEL_SPIRE_PSW", "HERSCHEL_SPIRE_PSW_EXT", "HST_ACS_HRC_F220W", "HST_ACS_HRC_F250W", "HST_ACS_HRC_F330W", "HST_ACS_HRC_F344N", "HST_ACS_HRC_F435W", "HST_ACS_HRC_F475W", "HST_ACS_HRC_F502N", "HST_ACS_HRC_F550M", "HST_ACS_HRC_F555W", "HST_ACS_HRC_F606W", "HST_ACS_HRC_F625W", "HST_ACS_HRC_F658N", "HST_ACS_HRC_F660N", "HST_ACS_HRC_F775W", "HST_ACS_HRC_F814W", "HST_ACS_HRC_F850LP", "HST_ACS_HRC_F892N", "HST_ACS_WFC_F435W", "HST_ACS_WFC_F475W", "HST_ACS_WFC_F502N", "HST_ACS_WFC_F550M", "HST_ACS_WFC_F555W", "HST_ACS_WFC_F606W", "HST_ACS_WFC_F625W", "HST_ACS_WFC_F658N", "HST_ACS_WFC_F660N", "HST_ACS_WFC_F775W", "HST_ACS_WFC_F814W", "HST_ACS_WFC_F850LP", "HST_ACS_WFC_F892N", "HST_NIC2_F110W", "HST_NIC2_F160W", "HST_NIC2_F205W", "HST_NIC3_F108N", "HST_NIC3_F110W", "HST_NIC3_F113N", "HST_NIC3_F150W", "HST_NIC3_F160W", "HST_NIC3_F164N", "HST_NIC3_F166N", "HST_NIC3_F175W", "HST_NIC3_F187N", "HST_NIC3_F190N", "HST_NIC3_F196N", "HST_NIC3_F200N", "HST_NIC3_F205M", "HST_NIC3_F212N", "HST_NIC3_F215N", "HST_NIC3_F222M", "HST_NIC3_F240M", "HST_WFC3_F098M", "HST_WFC3_F105W", "HST_WFC3_F110W", "HST_WFC3_F125W", "HST_WFC3_F126N", "HST_WFC3_F127M", "HST_WFC3_F128N", "HST_WFC3_F130N", "HST_WFC3_F132N", "HST_WFC3_F139M", "HST_WFC3_F140W", "HST_WFC3_F153M", "HST_WFC3_F160W", "HST_WFC3_F164N", "HST_WFC3_F167N", "HST_WFC3_F200LP", "HST_WFC3_F218W", "HST_WFC3_F225W", "HST_WFC3_F275W", "HST_WFC3_F280N", "HST_WFC3_F300X", "HST_WFC3_F336W", "HST_WFC3_F343N", "HST_WFC3_F350LP", "HST_WFC3_F373N", "HST_WFC3_F390M", "HST_WFC3_F390W", "HST_WFC3_F395N", "HST_WFC3_F410M", "HST_WFC3_F438W", "HST_WFC3_F467M", "HST_WFC3_F469N", "HST_WFC3_F475W", "HST_WFC3_F475X", "HST_WFC3_F487N", "HST_WFC3_F502N", "HST_WFC3_F547M", "HST_WFC3_F555W", "HST_WFC3_F600LP", "HST_WFC3_F606W", "HST_WFC3_F621M", "HST_WFC3_F625W", "HST_WFC3_F631N", "HST_WFC3_F645N", "HST_WFC3_F656N", "HST_WFC3_F657N", "HST_WFC3_F658N", "HST_WFC3_F665N", "HST_WFC3_F673N", "HST_WFC3_F680N", "HST_WFC3_F689M", "HST_WFC3_F763M", "HST_WFC3_F775W", "HST_WFC3_F814W", "HST_WFC3_F845M", "HST_WFC3_F850LP", "HST_WFC3_F953N", "HST_WFC3_FQ232N", "HST_WFC3_FQ243N", "HST_WFC3_FQ378N", "HST_WFC3_FQ387N", "HST_WFC3_FQ422M", "HST_WFC3_FQ436N", "HST_WFC3_FQ437N", "HST_WFC3_FQ492N", "HST_WFC3_FQ508N", "HST_WFC3_FQ575N", "HST_WFC3_FQ619N", "HST_WFC3_FQ634N", "HST_WFC3_FQ672N", "HST_WFC3_FQ674N", "HST_WFC3_FQ727N", "HST_WFC3_FQ750N", "HST_WFC3_FQ889N", "HST_WFC3_FQ906N", "HST_WFC3_FQ924N", "HST_WFC3_FQ937N", "HST_WFPC2_F170W", "HST_WFPC2_F218W", "HST_WFPC2_F255W", "HST_WFPC2_F300W", "HST_WFPC2_F336W", "HST_WFPC2_F439W", "HST_WFPC2_F450W", "HST_WFPC2_F555W", "HST_WFPC2_F606W", "HST_WFPC2_F622W", "HST_WFPC2_F675W", "HST_WFPC2_F791W", "HST_WFPC2_F814W", "HST_WFPC2_F850LP", "JWST_NIRCAM_F070W", "JWST_NIRCAM_F090W", "JWST_NIRCAM_F115W", "JWST_NIRCAM_F140M", "JWST_NIRCAM_F150W", "JWST_NIRCAM_F150W2", "JWST_NIRCAM_F162M", "JWST_NIRCAM_F164N", "JWST_NIRCAM_F182M", "JWST_NIRCAM_F187N", "JWST_NIRCAM_F200W", "JWST_NIRCAM_F210M", "JWST_NIRCAM_F212N", "JWST_NIRCAM_F250M", "JWST_NIRCAM_F277W", "JWST_NIRCAM_F300M", "JWST_NIRCAM_F322W2", "JWST_NIRCAM_F323N", "JWST_NIRCAM_F335M", "JWST_NIRCAM_F356W", "JWST_NIRCAM_F360M", "JWST_NIRCAM_F405N", "JWST_NIRCAM_F410M", "JWST_NIRCAM_F430M", "JWST_NIRCAM_F444W", "JWST_NIRCAM_F460M", "JWST_NIRCAM_F466N", "JWST_NIRCAM_F470N", "JWST_NIRCAM_F480M", "KEPLER_Kp", "NGTS_I", "PS1_g", "PS1_i", "PS1_r", "PS1_w", "PS1_y", "PS1_z", "SDSS_g", "SDSS_i", "SDSS_r", "SDSS_u", "SDSS_z", "SPITZER_IRAC_36", "SPITZER_IRAC_45", "SPITZER_IRAC_58", "SPITZER_IRAC_80", "STROMGREN_b", "STROMGREN_u", "STROMGREN_v", "STROMGREN_y", "TESS", "TYCHO_B_MvB", "TYCHO_V_MvB", "WISE_RSR_W1", "WISE_RSR_W2", "WISE_RSR_W3", "WISE_RSR_W4", "ZTF_g", "ZTF_i", "ZTF_r"]These included filter curves can be accessed like so,
using PhotometricFilters: SDSS_u, SDSS_g, SDSS_r, SDSS_i, SDSS_z, fwhm
filts = [SDSS_u(), SDSS_g(), SDSS_r(), SDSS_i(), SDSS_z()]5-element Vector{PhotometricFilter{Float64, Vector{Unitful.Quantity{Float64, 𝐋, Unitful.FreeUnits{(Å,), 𝐋, nothing}}}, PhotometricFilters.Photon, Vector{Float64}, String, Interpolations.FilledExtrapolation{Float64, 1, Interpolations.GriddedInterpolation{Float64, 1, Vector{Float64}, Interpolations.Gridded{Interpolations.Linear{Interpolations.Throw{Interpolations.OnGrid}}}, Tuple{Vector{Unitful.Quantity{Float64, 𝐋, Unitful.FreeUnits{(Å,), 𝐋, nothing}}}}}, Interpolations.Gridded{Interpolations.Linear{Interpolations.Throw{Interpolations.OnGrid}}}, Float64}}}:
SDSS_u
SDSS_g
SDSS_r
SDSS_i
SDSS_zNOTE THAT THESE INCLUDED FILTER CURVES ARE NOT GUARANTEED TO BE UP-TO-DATE. If you are using a filter/instrument that may have recently had its filter curves updated (e.g., JWST/NIRCAM), you should use our SVO query interface to make sure you get the most up-to-date filter curves. SVO also provides additional metadata that is useful for some applications (e.g., filter zeropoints). If you know the SVO-designated name of the filter you want, you can use get_filter to retrieve its transmission data, which returns an instance of SVOFilter.
PhotometricFilters.get_filter — Functionget_filter(filtername::AbstractString, magsys::Symbol=:Vega)Query the online SVO filter service for data on a photometric filter.
Arguments
filtername::AbstractString: The desired filter ID, in the correct SVO specification (e.g.,"2MASS/2MASS.J").magsys::Symbol: Desired magnitude system for associated metadata (e.g.,"ZeroPoint"). Can be any of(:AB, :Vega, :ST). SVO uses Vega by default, so we mirror that choice here.
Returns
An SVOFilter instance containing the results of the query.
Examples
julia> using PhotometricFilters: get_filter
julia> filt = get_filter("2MASS/2MASS.J", :Vega)
107-element PhotometricFilters.SVOFilter{PhotometricFilter{Float64}}: 2MASS/2MASS.J
reference wave.: 12350.0 Å
min. wave.: 10806.470589792389 Å
max. wave.: 14067.974683578484 Å
effective wave.: 12284.994608629975 Å
mean wave.: 12410.5170694321 Å
central wave.: 12390.584132888223 Å
pivot wave.: 12393.093155655275 Å
eff. width: 1624.3245065600008 Å
fwhm: 2149.1445403830403 ÅPhotometricFilters.SVOFilter — TypeSVOFilter(filter::PhotometricFilter, metadata) <: AbstractFilterType for containing the photometric filter information returned by the SVO filter service. A result of this type is returned by get_filter. Contains two fields:
filteris aPhotometricFiltertype that is used to support common operations.metadatais a dictionary (currently anOrderedCollections.OrderedDict) that contains the metadata returned by SVO.
These fields are considered internal (subject to change) and users should interact with instances of this type via the public accessor methods instead. Example usage is below.
julia> using PhotometricFilters: get_filter, SVOFilter, PhotometricFilter
julia> filt = get_filter("2MASS/2MASS.J", :Vega);
julia> filt isa SVOFilter
true
julia> PhotometricFilter(filt) isa PhotometricFilter # Access simpler PhotometricFilter type
true
julia> Dict(filt) isa Dict # Retrieve full metadata dictionary
true
julia> filt["ZeroPoint"] # Can retrieve metadata directly
1594.0 Jy
julia> filtername(filt) # `filtername`, `detector_type`, `wavelength`, `throughput` all work
"2MASS/2MASS.J"If you'd like to perform a search on the filters available through the SVO filter service, you can use query_filters.
PhotometricFilters.query_filters — Functionquery_filters(; queries...)Queries the filters available from the SVO filter service with search parameters and returns a table of the filters found.
The available search parameters can be found with PhotometricFilters.get_metadata. The following should be available in general:
WavelengthRef: Tuple of NumbersWavelengthMean: Tuple of NumbersWavelengthEff: Tuple of NumbersWavelengthMin: Tuple of NumbersWavelengthMax: Tuple of NumbersWidthEff: Tuple of NumbersFWHM: Tuple of NumbersInstrument: StringFacility: StringPhotSystem: String
The returned table is a DataFrame from the DataFrames package with all the columns of the response VOTable.
The filter information and transmission data can be obtained by calling get_filter with the ID obtained from the filterID column.
Examples
julia> using PhotometricFilters: query_filters, SVOFilter
julia> using DataFrames: DataFrame
julia> df = query_filters(; Facility="SLOAN", WavelengthEff=(1000, 5000));
julia> df isa DataFrame
true
julia> id = df.filterID[3]
"SLOAN/SDSS.g"
julia> get_filter(id) isa SVOFilter
true# Other examples for querying
query_filters(; Facility="SLOAN") # all filters from a given facility
query_filters(; Instrument="BUSCA", WavelengthEff=(1000u"angstrom", 5000u"angstrom")) # Unitful wavelengthsInteracting with the Filter Cache
After you first access a filter with get_filter, it is cached to disk for future use. It is expected that users should not typically have to manually interact with the cache. As such, the cache-related utilities discussed here are not exported from the package and must be explicitly imported (e.g., using PhotometricFilters: update_filters; update_filters()) or used via the qualified syntax (e.g., using PhotometricFilters; PhotometricFilters.update_filters()).
You can list the currently cached filters with PhotometricFilters.cached_filters. To update cached filters, ensuring you have the most up-to-date data, you can use PhotometricFilters.update_filter. You can delete filters from the cache with PhotometricFilters.clear_filter.
PhotometricFilters.cached_filters — Functioncached_filters()Returns a Vector{Tuple{String, Symbol}} containing the filter identifier and magnitude system for each SVO filter in the cache.
julia> using PhotometricFilters: cached_filters, get_filter
julia> get_filter("2MASS/2MASS.J", :Vega); # Load SVO filter, will be cached if not already
julia> ("2MASS/2MASS.J", :Vega) in cached_filters() # Check that filter is in the cache
truePhotometricFilters.update_filter — Functionupdate_filter(f::AbstractString, magsys::Symbol)Reacquires filter with SVO identifier f in the magnitude system magsys from SVO and saves it into the cache. Returns the updated filter.
julia> using PhotometricFilters: cached_filters, update_filter, SVOFilter
julia> get_filter("2MASS/2MASS.J", :Vega); # Load SVO filter, will be cached if not already
julia> update_filter("2MASS/2MASS.J", :Vega) isa SVOFilter # Updates cached file, returns filter
trueupdate_filter()When called with no arguments, updates all filters in the cache. Returns nothing.
PhotometricFilters.clear_filter — Functionclear_filter(f::AbstractString, magsys::Symbol)Deletes filter f in magnitude system magsys from the cache.
julia> using PhotometricFilters: get_filter, clear_filter, cached_filters
julia> get_filter("2MASS/2MASS.J", :Vega); # Ensure filter in cache
julia> clear_filter("2MASS/2MASS.J", :Vega) # Remove filter from cache
julia> ("2MASS/2MASS.J", :Vega) in cached_filters() # Check that filter was removed from cache
falseclear_filter()When called with no arguments, deletes all filters from the cache.
We include functions for performing many common operations on photometric filters, summarized below.
Applying Filter Curves to Spectra
PhotometricFilters.apply_throughput — Functionapply_throughput(f::AbstractFilter, wavelengths, flux)Use linear interpolation to map the wavelengths of the photometric filter f to the given wavelengths and apply the filter throughput to the flux. The provided wavelengths and those of the filter must be compatible. This means if one has units, the other one needs units, too.
julia> using PhotometricFilters: SDSS_u, wave_unit, apply_throughput
julia> f = SDSS_u();
julia> λ = 3000:4000
3000:4000
julia> flux = fill(1.0, length(λ)); # If `flux` is all `1`, `apply_throughput` reduces to `f` interpolated at `λ`
julia> apply_throughput(f, λ, flux) == f(λ)
true
julia> λ_u = λ .* wave_unit # Can also put units on λ
(3000:4000) Å
julia> apply_throughput(f, λ_u, flux) == f.(λ_u)
truePhotometricFilters.apply_throughput! — Functionapply_throughput!(f::AbstractFilter, wavelengths, flux, out)In-place version of apply_throughput which modifies out. It should have a compatible element type with flux.
PhotometricFilters.mean_flux_density — Functionmean_flux_density(filt::AbstractFilter, wavelengths, flux)Returns the mean flux density of a spectrum (defined by wavelengths wavelengths and fluxes flux) when integrated over the provided filter filt.
For photon counting detectors, this is
\[\overline{f_\lambda} = \frac{\int_\lambda \lambda \, f_\lambda \, T(\lambda) \, d\lambda}{\int_\lambda \lambda \, T(\lambda) \, d\lambda}\]
which can also be interpreted as the mean photon rate density, while for energy counting detectors, this is
\[\overline{f_\lambda} = \frac{\int_\lambda f_\lambda \, T(\lambda) \, d\lambda}{\int_\lambda T(\lambda) \, d\lambda}\]
which is essentially just the mean flux weighted by the filter throughput.
Below we show example usage that can be compared against this example from pyphot.
julia> using PhotometricFilters: mean_flux_density, HST_WFC3_F110W, Vega
julia> using Unitful, UnitfulAstro
julia> v = Vega("alpha_lyr_stis_006");
julia> mfd = mean_flux_density(HST_WFC3_F110W(), v.wave, v.flux);
julia> isapprox(mfd, 4.082289e-10 * u"erg/s/cm^2/angstrom"; rtol=1e-3)
truePhotometricFilters.F_nu — FunctionF_nu(F_lambda::SpectralFluxDensity, λref)
F_nu(F_lambda::SpectralFluxDensity, f::AbstractFilter)Convert a spectral flux density F_lambda into a spectral energy density. Assuming F_lambda in erg / s / cm^2 / Angstrom, and F_nu in Jy, this conversion is
\[F_\nu = \frac{10^5}{10^{-8} \, c} \, \lambda^2_r \, F_\lambda\]
where $c$ is the speed of light in m/s and $\lambda_r$ is the reference wavelength (λref) in Angstroms. If providing an AbstractFilter as the second argument, the reference wavelength will be automatically computed with reference_wavelength.
PhotometricFilters.F_lambda — FunctionF_lambda(F_nu::SpectralEnergyDensity, λref)
F_lambda(F_nu::SpectralEnergyDensity, f::AbstractFilter)Convert a spectral energy density F_nu into a spectral flux density. Assuming F_nu in Jy and F_lambda in erg / s / cm^2 / Angstrom, this conversion is
\[F_\lambda = \frac{10^{-8} \, c}{10^5} \, \lambda^{-2}_r \, F_\nu\]
where $c$ is the speed of light in m/s and $\lambda_r$ is the reference wavelength (λref) in Angstroms. If providing an AbstractFilter as the second argument, the reference wavelength will be automatically computed with reference_wavelength.
Zero Points
We utilize multiple dispatch to support dynamic calculation of zeropoints in the magnitude systems below.
PhotometricFilters.MagnitudeSystem — TypeAbstract supertype for magnitude systems like AB, ST, and Vega. Subtypes should implement zeropoint_mag, zeropoint_flux, and zeropoint_Jy.
PhotometricFilters.AB — TypeAB() <: MagnitudeSystemSingleton struct representing the AB magnitude system. This system is defined such that
\[m_\text{AB} = -2.5 \ \text{log} \left( \bar{f_ν} \right) - 48.60\]
when $f_ν$ is in units of erg / s / cm^2 / Hz. This corresponds to a constant zeropoint_Jy value in all filters of exp10(48.6 / -2.5 + 23) ≈ 3631 Jy.
When passed to methods such as zeropoint_flux, indicates that you wish to have the AB zeropoint flux returned.
PhotometricFilters.ST — TypeST() <: MagnitudeSystemSingleton struct representing the ST magnitude system. This system is defined so that a source with uniform $f_\lambda$ has identical magnitude in every filter.
\[m_\text{ST} = -2.5 \ \text{log} \left( \bar{f_\lambda} \right) - 21.1\]
When passed to methods such as zeropoint_flux, indicates that you wish to have the ST zeropoint flux returned.
PhotometricFilters.Vega — TypeVega(wave, flux, name::String) <: MagnitudeSystemStruct for containing a Vega reference spectrum with wavelengths wave and flux values flux. wave should be provided in units of Å and flux should be provided in $f_λ$ units of erg / s / cm^2 / Angstrom.
The Vega magnitude system is defined so that the star Alpha Lyr (i.e., Vega) has magnitude 0 in every filter.
\[m_\text{Vega} = -2.5 \ \text{log} \left( \frac{\bar{f_\lambda}}{\bar{f_\lambda} \left( \text{Vega} \right)} \right)\]
When passed to methods such as zeropoint_flux, indicates that you wish to have the Vega zeropoint flux returned.
Any reference spectrum hosted by CALSPEC can be used to construct an instance of Vega using the method below.
PhotometricFilters.Vega — MethodVega(name::String = "alpha_lyr_stis_011")Loads the reference spectrum with filename name and returns an appropriate instance of Vega that can be used to compute zeropoints and magnitudes in the Vega magnitude system.
If the provided name is the full path to the existing file on disk, the spectrum is loaded from that file. Otherwise, it is downloaded from the CALSPEC database of standard stars maintained by STScI. Files downloaded this way are cached for future use. Specifically, we draw from the extended catalog here. The CALSPEC catalog with the most recent reference spectrum for each star is located here. Standard Vega spectra start with "alpha_lyr". names like "alpha_lyr_stis_XXX" are based on composites of calibrated stellar models and HST/STIS data, while "alpha_lyr_mod_XXX" are based on stellar models only.
Sometimes Vega is not used as the standard star for photometric systems even when the system follows the Vega magnitude convention. For example, in the near-IR Sirius is often used as the standard reference spectrum rather than Vega. This is the case for the definition of the JWST zeropoints, which presently use the "sirius_stis_005.fits" CALSPEC spectrum as their standard. To load a different standard, simply provide the corresponding name to Vega. For example, to load the Sirius spectrum, use Vega("sirius_stis_005") (the ".fits" extension is optional).
The list of available spectral standards can be retrieved with PhotometricFilters.get_calspec_names.
PhotometricFilters.get_calspec_names — Functionget_calspec_names([substring::AbstractString])Returns a list of the names of the available spectral standards that can be download from CALSPEC and used as a standard in the Vega magnitude system.
If the optional substring::AbstractString argument is provided, then the list of names is filtered to only include those that contain the provided substring.
This method is not currently exported.
julia> using PhotometricFilters: Vega, get_calspec_names
julia> names = get_calspec_names();
julia> names isa Vector{String}
true
julia> Vega(names[1]) isa Vega
true
julia> vega_standards = get_calspec_names("alpha_lyr");
julia> all(map(x -> occursin("alpha_lyr", x), vega_standards))
trueZeropoints can be computed with methods below.
PhotometricFilters.zeropoint_flux — Functionzeropoint_flux(f::AbstractFilter, T::MagnitudeSystem)Returns the flux zero point of the filter f in magnitude system T in units of erg / s / cm^2 / Angstrom.
julia> using PhotometricFilters: zeropoint_flux, AB, ST, Vega, HST_WFC3_F110W
julia> using Unitful
julia> isapprox(zeropoint_flux(HST_WFC3_F110W(), AB()), 8.159816925e-10 * u"erg/s/cm^2/angstrom"; rtol=1e-3)
true
julia> isapprox(zeropoint_flux(HST_WFC3_F110W(), ST()), 3.6307805e-9 * u"erg/s/cm^2/angstrom"; rtol=1e-3)
true
julia> isapprox(zeropoint_flux(HST_WFC3_F110W(), Vega("alpha_lyr_stis_006")), 4.082289e-10 * u"erg/s/cm^2/angstrom"; rtol=1e-3)
truePhotometricFilters.zeropoint_Jy — Functionzeropoint_Jy(f::AbstractFilter, T::MagnitudeSystem)Returns the flux zeropoint in Jansky in magnitude system T.
Note that for the AB system, this is often approximated as 3631 Jy, following from the definition $m_\text{AB} = -2.5 \, \text{log} \left( \bar{f_\nu} \right) - 48.6$ where $\bar{f_\nu}$ is in units of erg / s / cm^2 / Hz. This can be solved for $m_\text{AB} = 0$ to give $\bar{f}_{\nu, 0} = 10^{\frac{48.6}{-2.5}}$ which is approximately $3.631 \times 10^{-20}$ erg / s / cm^2 / Hz, or ≈ 3631 Jy. This function returns the exact value.
julia> using PhotometricFilters: zeropoint_Jy, AB, ST, Vega, HST_WFC3_F110W
julia> using Unitful, UnitfulAstro
julia> isapprox(zeropoint_Jy(HST_WFC3_F110W(), AB()), 3630.78054 * u"Jy"; rtol=1e-3)
true
julia> isapprox(zeropoint_Jy(HST_WFC3_F110W(), ST()), 16155.46954* u"Jy"; rtol=1e-3)
true
julia> isapprox(zeropoint_Jy(HST_WFC3_F110W(), Vega("alpha_lyr_stis_006")), 1816.43597 * u"Jy"; rtol=1e-3)
truePhotometricFilters.zeropoint_mag — Functionzeropoint_mag(f::AbstractFilter, T::MagnitudeSystem)Returns the magnitude zero point of the filter f in the magnitude system T. This is used by the magnitude method to calculate magnitudes from spectra in units of F_lambda as
\[m_{\text{AB}} = -2.5 * \text{log} \left( \bar{f_\lambda} \right) - \text{Zpt}\]
For the ST magnitude system, this is always equal to 21.1 by definition.
julia> using PhotometricFilters: zeropoint_mag, AB, ST, Vega, HST_WFC3_F110W
julia> isapprox(zeropoint_mag(HST_WFC3_F110W(), AB()), 22.7207989; rtol=1e-3)
true
julia> isapprox(float(zeropoint_mag(HST_WFC3_F110W(), ST())), 21.1; rtol=1e-3)
true
julia> isapprox(zeropoint_mag(HST_WFC3_F110W(), Vega("alpha_lyr_stis_006")), 23.4727487; rtol=1e-3)
trueSynthetic Photometry
PhotometricFilters.magnitude — Functionmagnitude(f::AbstractFilter, T::MagnitudeSystem, wavelengths, flux)Calculates the magnitude in the given filter f in the magnitude system T from a spectrum defined by arrays wavelengths and flux, both of which must have valid Unitful units.
julia> using PhotometricFilters: magnitude, Vega, ST, AB, HST_WFC3_F110W
julia> v = Vega("alpha_lyr_stis_006");
julia> isapprox(magnitude(HST_WFC3_F110W(), AB(), v.wave, v.flux), 0.7519497; rtol=1e-3)
true
julia> isapprox(magnitude(HST_WFC3_F110W(), ST(), v.wave, v.flux), 2.372748728; rtol=1e-3)
true
julia> isapprox(magnitude(HST_WFC3_F110W(), v, v.wave, v.flux), 0; rtol=1e-3)
trueStatistics
PhotometricFilters.reference_wavelength — Functionreference_wavelength(f::AbstractFilter)Returns the reference wavelength of the filter f, used for conversions of the flux and for determination of magnitudes.
By default the pivot wavelength is returned (pivot_wavelength), but filter providers sometimes provide their own specified values.
PhotometricFilters.mean_wavelength — Functionmean_wavelength(f::AbstractFilter)Returns the mean wavelength of the filter f, defined as
\[\frac{\int \lambda \, T(\lambda) \, d\lambda}{\int T(\lambda) \, d\lambda}\]
PhotometricFilters.central_wavelength — Functioncentral_wavelength(f::AbstractFilter)Returns the central wavelength of the filter f, defined as the central wavelength between the two wavelengths used for the FWHM (fwhm).
PhotometricFilters.effective_wavelength — Functioneffective_wavelength(f::AbstractFilter, v::Vega = Vega())Returns the effective wavelength of the filter f using the Vega spectrum contained in v as a standard. Defined as
\[\frac{\int \lambda \, T(\lambda) \text{Vg}(\lambda) \, d\lambda}{\int T(\lambda) \text{Vg}(\lambda) \, d\lambda}\]
where $T(\lambda)$ is the filter transmission at wavelength $\lambda$ and $\text{Vg}(\lambda)$ is the spectrum of Vega.
PhotometricFilters.pivot_wavelength — Functionpivot_wavelength(f::AbstractFilter)Returns the pivot wavelength of the filter f, defined for filters with Energy detector types as
\[\sqrt{ \frac{\int T(\lambda) \, d\lambda}{\int T(\lambda) / \lambda^2 \, d\lambda} }\]
For filters with Photon detector types, $\lambda \, T(\lambda)$ is substituted for $T(\lambda)$ in the above expression.
Internally integration is carried out using trapezoidal integration. It can be convenient to think of this as the "center of mass" of the filter.
PhotometricFilters.min_wave — Functionmin_wave(f::AbstractFilter; level=0.01)Returns the shortest wavelength at which the filter transmission is equal to level * maximum(transmission).
PhotometricFilters.max_wave — Functionmax_wave(f::AbstractFilter; level=0.01)Returns the longest wavelength at which the filter transmission is equal to level * maximum(transmission).
PhotometricFilters.fwhm — Functionfwhm(f::AbstractFilter)Returns the difference between the furthest two wavelengths for which the filter transmission is equal to half its maximum value.
PhotometricFilters.width — Functionwidth(f::AbstractFilter)Returns the effective width of the filter, defined as the horizontal size of a rectangle with height equal to the maximum transmission of the filter such that the area of the rectangle is equal to the area under the filter transmission curve. This is calculated as
\[\frac{\int T(\lambda) \, d\lambda}{\text{max}(T(\lambda))}\]
Internals
PhotometricFilters.AbstractFilter — TypeAbstractFilter{T}Abstract supertype for representing photometric filters. Most functions provided by this package (e.g., effective_wavelength and similar methods) are designed to work with any subtype of AbstractFilter so long as a minimal API is defined for new subtypes. The methods that should be implemented for new types to conform to this API are summarized below:
filtername(f::NewType)should return a string indicating a human-readable name for the filter (e.g., "SDSS_u").wavelength(f::NewType)should return the wavelength vector of the filter transmission curve with properUnitful.jlunits.throughput(f::NewType)should return the throughput vector of the filter transmission curve (no units).detector_type(f::NewType)should return an instance ofPhotometricFilters.Energyif the filter is defined for energy-counting detectors orPhotometricFilters.Photonfor photon-counting detectors.
Additionally, all subtypes should support filter interpolation at user-defined wavelengths with a call signature (f::NewType)(wavelengths). To support this, new types should implement a method like (f::PhotometricFilter)(wavelength::Q) where Q <: Unitful.Length. A generic fallback for inputs without units is already defined.
PhotometricFilters.get_metadata — Functionget_metadata()Returns a table of the available parameters that can be used to query the SVO filter service from the FORMAT=metadata VOTable they provide.
The table is a DataFrame from the DataFrames package with the following columns:
parameter: parameter name that can be used for queries usingquery_filtersunit:Unitfulunit of the parameterdatatype:Typeof the parameterdescription: description of the parametervalues: vector of the possible values that the respective parameter can take on (e.g. for Instrument), or a vector of the minimum and maximum values that the parameter can assume (e.g. for WavelengthEff)
Example
julia> using DataFrames: DataFrame
julia> df = PhotometricFilters.get_metadata();
julia> df isa DataFrame
true
julia> facilities = df[findfirst(==("Facility"), df.parameter), :].values;
julia> facilities isa Vector{String}
trueThis is not exported.
Index
PhotometricFilters.ABPhotometricFilters.AbstractFilterPhotometricFilters.MagnitudeSystemPhotometricFilters.PhotometricFilterPhotometricFilters.STPhotometricFilters.SVOFilterPhotometricFilters.VegaPhotometricFilters.VegaPhotometricFilters.F_lambdaPhotometricFilters.F_nuPhotometricFilters.apply_throughputPhotometricFilters.apply_throughput!PhotometricFilters.cached_filtersPhotometricFilters.central_wavelengthPhotometricFilters.clear_filterPhotometricFilters.detector_typePhotometricFilters.effective_wavelengthPhotometricFilters.filternamePhotometricFilters.fwhmPhotometricFilters.get_calspec_namesPhotometricFilters.get_filterPhotometricFilters.get_metadataPhotometricFilters.magnitudePhotometricFilters.max_wavePhotometricFilters.mean_flux_densityPhotometricFilters.mean_wavelengthPhotometricFilters.min_wavePhotometricFilters.pivot_wavelengthPhotometricFilters.query_filtersPhotometricFilters.reference_wavelengthPhotometricFilters.throughputPhotometricFilters.update_filterPhotometricFilters.wavelengthPhotometricFilters.widthPhotometricFilters.zeropoint_JyPhotometricFilters.zeropoint_fluxPhotometricFilters.zeropoint_mag