pyphot.astropy package#
Submodules#
pyphot.astropy.config module#
pyphot.astropy.sandbox module#
Sandbox of new developments
Use at your own risks
Photometric package using Astropy Units#
Defines a Filter class and associated functions to extract photometry.
This also include functions to keep libraries up to date
Note
integrations are done using trapezoid()
Why not Simpsons? Simpsons principle is to take sequence of 3 points to
make a quadratic interpolation. Which in the end, when filters have sharp
edges, the error due to this “interpolation” are extremely large in
comparison to the uncertainties induced by trapeze integration.
- class Constants[source]#
Bases:
object
A namespace for constants
- c = <Quantity 2.99792458e+18 Angstrom / s>#
- h = <Quantity 6.62607015e-27 erg s>#
- class UncertainFilter[source]#
Bases:
UnitFilter
What could be a filter with uncertainties
- wavelength#
wavelength sequence defining the filter transmission curve
- Type:
ndarray
- name#
name of the passband
- Type:
string
- dtype#
detector type, either “photon” or “energy” counter
- Type:
str
- unit#
wavelength units
- Type:
str
- property AB_zero_Jy#
AB flux zero point in Jansky (Jy)
- property AB_zero_flux#
AB flux zero point in erg/s/cm2/AA
- property AB_zero_mag#
AB magnitude zero point
- ABmag = -2.5 * log10(f_nu) - 48.60
= -2.5 * log10(f_lamb) - 2.5 * log10(lpivot ** 2 / c) - 48.60 = -2.5 * log10(f_lamb) - zpts
- property ST_zero_Jy#
ST flux zero point in Jansky (Jy)
- property ST_zero_flux#
ST flux zero point in erg/s/cm2/AA
- property ST_zero_mag#
ST magnitude zero point STmag = -2.5 * log10(f_lamb) -21.1
- property Vega_zero_Jy#
Vega flux zero point in Jansky (Jy)
- property Vega_zero_flux#
Vega flux zero point in erg/s/cm2/AA
- property Vega_zero_mag#
Vega magnitude zero point Vegamag = -2.5 * log10(f_lamb) + 2.5 * log10(f_vega) Vegamag = -2.5 * log10(f_lamb) - zpts
- property Vega_zero_photons#
Vega number of photons per wavelength unit
Note
see self.get_Nphotons
- __init__(wavelength, mean_transmit, samples, name='', dtype='photon', unit=None)[source]#
Constructor
- apply_transmission(slamb, sflux)[source]#
Apply filter transmission to a spectrum (with reinterpolation of the filter)
- Parameters:
slamb (ndarray) – spectrum wavelength definition domain
sflux (ndarray) – associated flux
- Returns:
flux – new spectrum values accounting for the filter
- Return type:
float
- property cl#
Unitwise wavelength definition
- classmethod from_gp_model(model, xprime=None, n_samples=10, **kwargs)[source]#
Generate a filter object from a sklearn GP model
- Parameters:
model (sklearn.gaussian_process.GaussianProcessRegressor) – model of the passband
xprime (ndarray) – wavelength to express the model in addition to the training points
n_samples (int) – number of samples to generate from the model.
kwargs (dict) – UncertainFilter keywords
- property fwhm#
the difference between the two wavelengths for which filter transmission is half maximum
- ..note::
This calculation is not exact but rounded to the nearest passband data points
- getFlux(slamb, sflux, axis=-1)[source]#
Integrate the flux within the filter and return the integrated energy If you consider applying the filter to many spectra, you might want to consider extractSEDs.
- Parameters:
slamb (ndarray(dtype=float, ndim=1)) – spectrum wavelength definition domain
sflux (ndarray(dtype=float, ndim=1)) – associated flux
- Returns:
flux – Energy of the spectrum within the filter
- Return type:
float
- get_Nphotons(slamb, sflux, axis=-1)[source]#
getNphot the number of photons through the filter (Ntot / width in the documentation)
getflux() * leff / hc
- Parameters:
slamb (ndarray(dtype=float, ndim=1)) – spectrum wavelength definition domain
sflux (ndarray(dtype=float, ndim=1)) – associated flux in erg/s/cm2/AA
- Returns:
N – Number of photons of the spectrum within the filter
- Return type:
float
- property leff#
Unitwise Effective wavelength leff = int (lamb * T * Vega dlamb) / int(T * Vega dlamb)
- property lmax#
Calculated as the last value with a transmission at least 1% of maximum transmission
- property lmin#
Calculate das the first value with a transmission at least 1% of maximum transmission
- property lphot#
Photon distribution based effective wavelength. Defined as
lphot = int(lamb ** 2 * T * Vega dlamb) / int(lamb * T * Vega dlamb)
which we calculate as
lphot = get_flux(lamb * vega) / get_flux(vega)
- property lpivot#
Unitwise wavelength definition
- to_Table(**kwargs)[source]#
Export filter to a SimpleTable object
- Parameters:
fname (str) – filename
parameters (Uses SimpleTable)
- property transmit#
Transmission curves
- property wavelength#
Unitwise wavelength definition
- property wavelength_unit#
Unit wavelength definition
- property width#
Effective width Equivalent to the horizontal size of a rectangle with height equal to maximum transmission and with the same area that the one covered by the filter transmission curve.
W = int(T dlamb) / max(T)
- class UnitAscii_Library[source]#
Bases:
UnitLibrary
Interface one or multiple directory or many files as a filter library
>>> lib = Ascii_Library(['ground', 'hst', 'myfilter.csv'])
- add_filters(filter_object, fmt='%.6f', **kwargs)[source]#
Add a filter to the library permanently
- Parameters:
filter_object (Filter object) – filter to add
- load_filters(names, interp=True, lamb=None, filterLib=None)[source]#
load a limited set of filters
- Parameters:
names (list[str]) – normalized names according to filtersLib
interp (bool) – reinterpolate the filters over given lambda points
lamb (ndarray[float, ndim=1]) – desired wavelength definition of the filter
filterLib (path) – path to the filter library hd5 file
- Returns:
filters – list of filter objects
- Return type:
list[filter]
- class UnitFilter[source]#
Bases:
object
Evolution of Filter that makes sure the input spectra and output fluxes have units to avoid mis-interpretation.
- Note the usual (non SI) units of flux definitions:
flam = erg/s/cm**2/AA fnu = erg/s/cm**2/Hz photflam = photon/s/cm**2/AA photnu = photon/s/cm**2/Hz
Define a filter by its name, wavelength and transmission The type of detector (energy or photon counter) can be specified for adapting calculations. (default: photon)
- name#
name of the filter
- Type:
str
- cl#
central wavelength of the filter
- Type:
float
- norm#
normalization factor of the filter
- Type:
float
- lpivot#
pivot wavelength of the filter
- Type:
float
- wavelength#
wavelength sequence defining the filter transmission curve
- Type:
ndarray
- transmit#
transmission curve of the filter
- Type:
ndarray
- dtype#
detector type, either “photon” or “energy” counter
- Type:
str
- unit#
wavelength units
- Type:
str
- property AB_zero_Jy#
AB flux zero point in Jansky (Jy)
- property AB_zero_flux#
AB flux zero point in erg/s/cm2/AA
- property AB_zero_mag#
AB magnitude zero point
- ABmag = -2.5 * log10(f_nu) - 48.60
= -2.5 * log10(f_lamb) - 2.5 * log10(lpivot ** 2 / c) - 48.60 = -2.5 * log10(f_lamb) - zpts
- property ST_zero_Jy#
ST flux zero point in Jansky (Jy)
- property ST_zero_flux#
ST flux zero point in erg/s/cm2/AA
- property ST_zero_mag#
ST magnitude zero point STmag = -2.5 * log10(f_lamb) -21.1
- property Vega_zero_Jy#
Vega flux zero point in Jansky (Jy)
- property Vega_zero_flux#
Vega flux zero point in erg/s/cm2/AA
- property Vega_zero_mag#
vega magnitude zero point vegamag = -2.5 * log10(f_lamb) + 2.5 * log10(f_vega) vegamag = -2.5 * log10(f_lamb) - zpts
- property Vega_zero_photons#
Vega number of photons per wavelength unit
Note
see self.get_Nphotons
- apply_transmission(slamb, sflux)[source]#
Apply filter transmission to a spectrum (with reinterpolation of the filter)
- Parameters:
slamb (ndarray) – spectrum wavelength definition domain
sflux (ndarray) – associated flux
- Returns:
flux – new spectrum values accounting for the filter
- Return type:
float
- property cl#
Unitwise wavelength definition
- property fwhm#
the difference between the two wavelengths for which filter transmission is half maximum
- ..note::
This calculation is not exact but rounded to the nearest passband data points
- getFlux(slamb, sflux, axis=-1)[source]#
Integrate the flux within the filter and return the integrated energy If you consider applying the filter to many spectra, you might want to consider extractSEDs.
- Parameters:
slamb (ndarray(dtype=float, ndim=1)) – spectrum wavelength definition domain
sflux (ndarray(dtype=float, ndim=1)) – associated flux
- Returns:
flux – Energy of the spectrum within the filter
- Return type:
float
- get_Nphotons(slamb, sflux, axis=-1)[source]#
getNphot the number of photons through the filter (Ntot / width in the documentation)
getflux() * leff / hc
- Parameters:
slamb (ndarray(dtype=float, ndim=1)) – spectrum wavelength definition domain
sflux (ndarray(dtype=float, ndim=1)) – associated flux in erg/s/cm2/AA
- Returns:
N – Number of photons of the spectrum within the filter
- Return type:
float
- get_flux(slamb, sflux, axis=-1)[source]#
getFlux Integrate the flux within the filter and return the integrated energy If you consider applying the filter to many spectra, you might want to consider extractSEDs.
- Parameters:
slamb (ndarray(dtype=float, ndim=1)) – spectrum wavelength definition domain
sflux (ndarray(dtype=float, ndim=1)) – associated flux
- Returns:
flux – Energy of the spectrum within the filter
- Return type:
float
- property leff#
Unitwise Effective wavelength leff = int (lamb * T * Vega dlamb) / int(T * Vega dlamb)
- property lmax#
Calculated as the last value with a transmission at least 1% of maximum transmission
- property lmin#
Calculate das the first value with a transmission at least 1% of maximum transmission
- property lphot#
Photon distribution based effective wavelength. Defined as
lphot = int(lamb ** 2 * T * Vega dlamb) / int(lamb * T * Vega dlamb)
which we calculate as
lphot = get_flux(lamb * vega) / get_flux(vega)
- property lpivot#
Unitwise wavelength definition
- classmethod make_integration_filter(lmin, lmax, name='', dtype='photon', unit=None)[source]#
Generate an heavyside filter between lmin and lmax
- to_Table(**kwargs)[source]#
Export filter to a SimpleTable object
- Parameters:
fname (str) – filename
parameters (Uses SimpleTable)
- property wavelength#
Unitwise wavelength definition
- property width#
Effective width Equivalent to the horizontal size of a rectangle with height equal to maximum transmission and with the same area that the one covered by the filter transmission curve.
W = int(T dlamb) / max(T)
- class UnitHDF_Library[source]#
Bases:
UnitLibrary
Storage based on HDF
- __init__(source=PosixPath('/home/runner/work/pyphot/pyphot/pyphot/libs/new_filters.hd5'), mode='r')[source]#
Construct the library
- add_filter(f, **kwargs)[source]#
Add a filter to the library permanently
- Parameters:
f (Filter object) – filter to add
- load_all_filters(interp=True, lamb=None)[source]#
load all filters from the library
- Parameters:
interp (bool) – reinterpolate the filters over given lambda points
lamb (ndarray[float, ndim=1]) – desired wavelength definition of the filter
- Returns:
filters – list of filter objects
- Return type:
list[filter]
- load_filters(names, interp=True, lamb=None, filterLib=None)[source]#
load a limited set of filters
- Parameters:
names (list[str]) – normalized names according to filtersLib
interp (bool) – reinterpolate the filters over given lambda points
lamb (ndarray[float, ndim=1]) – desired wavelength definition of the filter
filterLib (path) – path to the filter library hd5 file
- Returns:
filters – list of filter objects
- Return type:
list[filter]
- class UnitLibrary[source]#
Bases:
object
Common grounds for filter libraries
- __init__(source=PosixPath('/home/runner/work/pyphot/pyphot/pyphot/libs/new_filters.hd5'), *args, **kwargs)[source]#
Construct the library
- property content#
Get the content list
- class UnitLickIndex[source]#
Bases:
object
Define a Lick Index similarily to a Filter object
- __init__(name, lick, unit='AA')[source]#
Constructor
- Parameters:
name (str) – name of the index
lick (dict) – expecting ‘blue’, ‘red’, ‘band’, and ‘unit’ definitions blue and red are used to continuum normalize the spectra band covers the index itself. unit gives the index measurement units, either magnitudes (mag) or equivalent width (ew)
unit (str) – wavelength unit of the intervals
- property band#
Unitwise band definition
- property blue#
Unitwise band definition
- classmethod continuum_normalized_region_around_line(wi, fi, blue, red, band=None, degree=1)[source]#
cut out and normalize flux around a line
- Parameters:
wi (ndarray (nw, )) – array of wavelengths in AA
fi (ndarray (N, nw)) – array of flux values for different spectra in the series
blue (tuple(2)) – selection for blue continuum estimate
red (tuple(2)) – selection for red continuum estimate
band (tuple(2), optional) – select region in this band only. default is band = (min(blue), max(red))
degree (int) – degree of the polynomial fit to the continuum
- Returns:
wnew (ndarray (nw1, )) – wavelength of the selection in AA
f (ndarray (N, len(wnew))) – normalized flux in the selection region
Example
# indice of CaII # wavelength are always supposed in AA w, f = region_around_line( wavelength, flux, [3925, 3930],[3938, 3945]] )
- get(wave, flux, **kwargs)[source]#
compute spectral index after continuum subtraction
- Parameters:
w (ndarray (nw, )) – array of wavelengths in AA
flux (ndarray (N, nw)) – array of flux values for different spectra in the series
degree (int (default 1)) – degree of the polynomial fit to the continuum
nocheck (bool) – set to silently pass on spectral domain mismatch. otherwise raises an error when index is not covered
- Returns:
ew – equivalent width or magnitude array
- Return type:
ndarray (N,)
- Raises:
ValueError – when the spectral coverage wave does not cover the index:
range –
- property index_unit#
- property red#
Unitwise band definition
- class UnitLickLibrary[source]#
Bases:
object
Collection of Lick indices
- __init__(fname=PosixPath('/home/runner/work/pyphot/pyphot/pyphot/libs/licks.dat'), comment='#')[source]#
- property content#
- property description#
any comment in the input file
- get_library(fname=PosixPath('/home/runner/work/pyphot/pyphot/pyphot/libs/new_filters.hd5'), **kwargs)[source]#
Finds the appropriate class to load the library
- reduce_resolution(wi, fi, fwhm0=<Quantity 0.55 Angstrom>, sigma_floor=<Quantity 0.2 Angstrom>)[source]#
Adapt the resolution of the spectra to match the lick definitions
Lick definitions have different resolution elements as function of wavelength. These definition are hard-coded in this function
- Parameters:
wi (ndarray (n, )) – wavelength definition
fi (ndarray (nspec, n) or (n, )) – spectra to convert
fwhm0 (float) – initial broadening in the spectra fi
sigma_floor (float) – minimal dispersion to consider
- Returns:
flux_red – reduced spectra
- Return type:
ndarray (nspec, n) or (n, )
pyphot.astropy.sun module#
Handle the Sun Spectrum
- class Sun[source]#
Bases:
object
Class that handles the Sun’s spectrum and references.
Observed solar spectrum comes from: ftp://ftp.stsci.edu/cdbs/current_calspec/sun_reference_stis_001.fits
and theoretical spectrum comes from: ftp://ftp.stsci.edu/cdbs/grid/k93models/standards/sun_kurucz93.fits
The theoretical spectrum is scaled to match the observed spectrum from 1.5 - 2.5 microns, and then it is used where the observed spectrum ends. The theoretical model of the Sun from Kurucz’93 atlas using the following parameters when the Sun is at 1 au.
log_Z T_eff log_g V_{Johnson} +0.0 5777 +4.44 -26.75
- source#
filename of the sun library
- Type:
str
- data#
data table
- Type:
- units#
detected units from file header
- Type:
tuple
- wavelength#
wavelength (with units when found)
- Type:
array
- flux#
flux(wavelength) values (with units when provided)
- Type:
array
- distance#
distance to the observed Sun (default, 1 au)
- Type:
float
- flavor#
either ‘observed’ using the stis reference, or ‘theoretical’ for the Kurucz model.
- Type:
str, (default theoretical)
- property flux#
- property wavelength#
pyphot.astropy.vega module#
Handle vega spec/mags/fluxes manipulations
Works with both ascii and hd5 files for back-compatibility
Vega.wavelength and Vega.flux have now units!
- class Vega[source]#
Bases:
object
Class that handles vega spectrum and references. This class know where to find the Vega synthetic spectrum (Bohlin 2007) in order to compute fluxes and magnitudes in given filters
- source#
filename of the vega library
- Type:
str
- data#
data table
- Type:
- units#
detected units from file header
- Type:
tuple
- wavelength#
wavelength (with units when found)
- Type:
array
- flux#
flux(wavelength) values (with units when provided)
- Type:
array
- An instance can be used as a context manager as
- >>> filters = ['HST_WFC3_F275W', 'HST_WFC3_F336W', 'HST_WFC3_F475W', 'HST_WFC3_F814W', 'HST_WFC3_F110W', 'HST_WFC3_F160W']
- with Vega() as v:
vega_f, vega_mag, flamb = v.getSed(filters)
print vega_f, vega_mag, flamb
- property flux#
flux(wavelength) values (with units when provided)
- property wavelength#
wavelength (with units when found)