pypeit.sensfunc module
Implements the objects used to construct sensitivity functions.
- class pypeit.sensfunc.IRSensFunc(spec1dfiles, sensfile, par, par_fluxcalib=None, debug=False, write_qa=True, chk_version=True)[source]
Bases:
SensFuncDetermine a sensitivity functions from standard-star spectra. Should only be used with NIR spectra (\(\lambda > 7000\) angstrom).
- Parameters:
spec1dfile (
str) – PypeIt spec1d file for the standard file.sensfile (
str) – File name for the sensitivity function data.par (
SensFuncPar, optional) – The parameters required for the sensitivity function computation.debug (
bool, optional) – Run in debug mode.
- _algorithm = 'IR'
Algorithm used for the sensitivity calculation.
- compute_zeropoint()[source]
Calls routine to compute the sensitivity function.
- Returns:
TelObj – Best-fitting telluric model
- Return type:
- eval_zeropoint(wave, iorddet)[source]
Evaluate the sensitivity function zero-points at the input wavelength
- Parameters:
wave (numpy.ndarray, shape is (nspec)) – Wavelength array
iorddet (
int) – Order or detector (0-indexed)
- Returns:
zeropoint – Zeropoint array evaluated at the input wavelength grid and with the gpm applied.
- Return type:
numpy.ndarray, shape is (nspec,)
- class pypeit.sensfunc.SensFunc(spec1dfiles, sensfile, par, par_fluxcalib=None, debug=False, write_qa=True, chk_version=True)[source]
Bases:
DataContainerBase class for generating sensitivity functions from a standard-star spectrum.
This class should not be instantated by itself; instead instantiate either
UVISSensFuncorIRSensFunc, depending on the wavelength range of your data (UVIS for \(\lambda < 7000\) angstrom, IR for \(\lambda > 7000\) angstrom.)The datamodel attributes are:
Version: 1.0.2
Attribute
Type
Array Type
Description
PYP_SPECstr
PypeIt spectrograph name
airmassfloat
Airmass of the observation
algorithmstr
Algorithm used for the sensitivity calculation.
exptimefloat
Exposure time
extrstr
Extraction method used for the standard star (OPT or BOX)
pypelinestr
PypeIt pipeline reduction path
sensTable with the sensitivity function
spec1dfstr
PypeIt spec1D file(s) used to for sensitivity function
std_calstr
File name (or shorthand) with the standard flux data
std_decfloat
DEC of the standard source
std_namestr
Type of standard source
std_rafloat
RA of the standard source
telluricTelluric model; see
Telluricthroughputfloat
Spectrograph throughput measurements
wavefloat
Wavelength vectors
zeropointfloat
Sensitivity function zeropoints
- Parameters:
spec1dfiles (
list) – PypeIt spec1d file(s) for the standard file.sensfile (
str) – File name for the sensitivity function data.par (
SensFuncPar) – The parameters required for the sensitivity function computation.par_fluxcalib (
FluxCalibratePar, optional) – The parameters required for flux calibration. These are only used for flux calibration of the standard star spectrum for the QA plot. If None, defaults will be used.debug (
bool, optional) – Run in debug mode, sending diagnostic information to the screen.chk_version (
bool, optional) – Check the version of the data model.write_qa (
bool, optional) – If True, write QA plots to PDF files.
- _algorithm = None
Algorithm used for the sensitivity calculation.
- compute_throughput()[source]
Compute the spectroscopic throughput
- Returns:
throughput (numpy.ndarray,
float, shape is (nspec, norders)) – Throughput measurementsthroughput_splice (numpy.ndarray,
float, shape is (nspec_splice, norders)) – Throughput measurements for spliced spectra
- datamodel = {'PYP_SPEC': {'descr': 'PypeIt spectrograph name', 'otype': <class 'str'>}, 'airmass': {'descr': 'Airmass of the observation', 'otype': <class 'float'>}, 'algorithm': {'descr': 'Algorithm used for the sensitivity calculation.', 'otype': <class 'str'>}, 'exptime': {'descr': 'Exposure time', 'otype': <class 'float'>}, 'extr': {'descr': 'Extraction method used for the standard star (OPT or BOX)', 'otype': <class 'str'>}, 'pypeline': {'descr': 'PypeIt pipeline reduction path', 'otype': <class 'str'>}, 'sens': {'descr': 'Table with the sensitivity function', 'otype': <class 'astropy.table.table.Table'>}, 'spec1df': {'descr': 'PypeIt spec1D file(s) used to for sensitivity function', 'otype': <class 'str'>}, 'std_cal': {'descr': 'File name (or shorthand) with the standard flux data', 'otype': <class 'str'>}, 'std_dec': {'descr': 'DEC of the standard source', 'otype': <class 'float'>}, 'std_name': {'descr': 'Type of standard source', 'otype': <class 'str'>}, 'std_ra': {'descr': 'RA of the standard source', 'otype': <class 'float'>}, 'telluric': {'descr': 'Telluric model; see :class:`~pypeit.core.telluric.Telluric`', 'otype': <class 'pypeit.core.telluric.Telluric'>}, 'throughput': {'atype': <class 'float'>, 'descr': 'Spectrograph throughput measurements', 'otype': <class 'numpy.ndarray'>}, 'wave': {'atype': <class 'float'>, 'descr': 'Wavelength vectors', 'otype': <class 'numpy.ndarray'>}, 'zeropoint': {'atype': <class 'float'>, 'descr': 'Sensitivity function zeropoints', 'otype': <class 'numpy.ndarray'>}}
DataContainer datamodel.
- static empty_sensfunc_table(norders, nspec, nspec_in, ncoeff=1)[source]
Construct an empty astropy.table.Table for the sensitivity function.
- Parameters:
norders (
int) – The number of slits/orders on the detector.nspec (
int) – The number of spectral pixels for the zeropoint arrays.nspec_in (
int) – The number of spectral pixels on the detector for the input standard star spectrum.ncoeff (
int, optional) – Number of coefficients for smooth model fit to zeropoints
- Returns:
Instance of the empty sensitivity function table.
- Return type:
- eval_zeropoint(wave, iorddet)[source]
Evaluate at a given wavelength the sensitivity function zeropoint for a given order/detector. This is a dummy method, overloaded by subclasses.
- Parameters:
wave (numpy.ndarray) – Wavelength array at which to evaluate the zeropoint
iorddet (
int) – The order/detector for which to evaluate the zeropoint
- Returns:
zeropoint – The zeropoint evaluated given wavelength array and order/detector
- Return type:
- extrapolate(samp_fact=1.5)[source]
Extrapolates the sensitivity function to cover an extra wavelength range set by the
extrapl_bluandextrap_redparameters. This is important for making sure that the sensitivity function can be applied to data with slightly different wavelength coverage etc.- Parameters:
samp_fact (
float) – Parameter governing the sampling of the wavelength grid used for the extrapolation.- Returns:
wave_extrap (numpy.ndarray) – Extrapolated wavelength array
zeropoint_extrap (numpy.ndarray) – Extrapolated sensitivity function
- classmethod from_hdu(hdu, hdu_prefix=None, chk_version=True)[source]
Instantiate the object from an HDU extension.
This overrides the base-class method, essentially just to handle the fact that the ‘TELLURIC’ extension is not called ‘MODEL’.
- Parameters:
hdu (astropy.io.fits.HDUList, astropy.io.fits.ImageHDU, astropy.io.fits.BinTableHDU) – The HDU(s) with the data to use for instantiation.
hdu_prefix (
str, optional) – Maintained for consistency with the base class but is not used by this method.chk_version (
bool, optional) – If True, raise an error if the datamodel version or type check failed. If False, throw a warning only.
- classmethod get_instance(spec1dfile, sensfile, par, par_fluxcalib=None, debug=False, write_qa=True, chk_version=True)[source]
Instantiate the relevant subclass based on the algorithm provided in
par.
- internals = ['sensfile', 'spec1d_arr', 'spectrograph', 'par', 'par_fluxcalib', 'qafile', 'thrufile', 'fstdfile', 'debug', 'sobjs_std', 'wave_cnts', 'counts', 'counts_ivar', 'counts_mask', 'log10_blaze_function', 'nspec_in', 'norderdet', 'wave_splice', 'zeropoint_splice', 'throughput_splice', 'steps', 'splice_multi_det', 'meta_spec', 'std_spec', 'write_qa', 'chk_version']
A list of strings identifying a set of internal attributes that are not part of the datamodel.
- classmethod sensfunc_weights(sensfile, waves, ech_order_vec=None, debug=False, extrap_sens=True, chk_version=True)[source]
Get the weights based on the sensfunc
- Parameters:
sensfile (str) – the name of your fits format sensfile
waves (numpy.ndarray) – wavelength grid for your output weights. Shape is (nspec, norders, nexp) or (nspec, norders).
ech_order_vec (numpy.ndarray, optional) – Vector of echelle orders. Only used for echelle data.
debug (bool) – default=False show the weights QA
extrap_sens (bool) – default=True Extrapolate the sensitivity function
chk_version (
bool, optional) – When reading in existing files written by PypeIt, perform strict version checking to ensure a valid file. If False, the code will try to keep going, but this may lead to faults and quiet failures. User beware!
- Returns:
sensfunc weights evaluated on the input waves wavelength grid, shape = same as waves
- Return type:
- splice()[source]
Routine to splice together sensitivity functions into one global sensitivity function for spectrographs with multiple detectors extending across the wavelength direction.
- Returns:
wave_splice (numpy.ndarray, shape is (nspec_splice, 1)) – wavelength array
zeropoint_splice (numpy.ndarray, shape is (nspec_splice, 1)) – zero-point array
- unpack_std()[source]
Unpack the standard star data from a 1D spectrum file(s) with a SpecObj or OneSpec class.
- Returns:
sobjs_std – The SpecObjs class with the standard star spectrum.
- Return type:
- version = '1.0.2'
Datamodel version.
- class pypeit.sensfunc.UVISSensFunc(spec1dfile, sensfile, par, par_fluxcalib=None, debug=False, write_qa=True, chk_version=True)[source]
Bases:
SensFuncDetermine a sensitivity functions from standard-star spectra. Should only be used with UVIS spectra (\(\lambda < 7000\) angstrom).
- Parameters:
spec1dfile (
str) – PypeIt spec1d file for the standard file.sensfile (
str) – File name for the sensitivity function data.par (
SensFuncPar, optional) – The parameters required for the sensitivity function computation.debug (
bool, optional) – Run in debug mode.
- _algorithm = 'UVIS'
Algorithm used for the sensitivity calculation.
- eval_zeropoint(wave, iorddet)[source]
Evaluate the sensitivity function zero-points at the input wavelength
- Parameters:
wave (numpy.ndarray, shape is (nspec)) – Wavelength array
iorddet (
int) – Order or detector (0-indexed)
- Returns:
zeropoint – Zeropoint array evaluated at the input wavelength grid and with the gpm applied.
- Return type:
numpy.ndarray, shape is (nspec,)