pypeit.sensfunc module

Implements the objects used to construct sensitivity functions.

class pypeit.sensfunc.IRSensFunc(spec1dfile, sensfile, par=None, debug=False)[source]

Bases: SensFunc

Determine 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

Telluric

eval_zeropoint(wave, iorddet)[source]

Evaluate the sensitivity function zero-points

Parameters
  • wave (numpy.ndarray, shape is (nspec)) – Wavelength array

  • iorddet (int) – Order or detector (0-indexed)

Returns

zeropoint

Return type

numpy.ndarray, shape is (nspec,)

class pypeit.sensfunc.SensFunc(spec1dfile, sensfile, par=None, debug=False)[source]

Bases: DataContainer

Base class for generating sensitivity functions from a standard-star spectrum.

This class should not be instantated by itself; instead instantiate either UVISSensFunc or IRSensFunc, 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.1

Attribute

Type

Array Type

Description

PYP_SPEC

str

PypeIt spectrograph name

airmass

float

Airmass of the observation

algorithm

str

Algorithm used for the sensitivity calculation.

exptime

float

Exposure time

pypeline

str

PypeIt pipeline reduction path

sens

astropy.table.table.Table

Table with the sensitivity function

spec1df

str

PypeIt spec1D file used to for sensitivity function

std_cal

str

File name (or shorthand) with the standard flux data

std_dec

float

DEC of the standard source

std_name

str

Type of standard source

std_ra

float

RA of the standard source

telluric

Telluric

Telluric model; see Telluric

throughput

numpy.ndarray

float

Spectrograph throughput measurements

wave

numpy.ndarray

float

Wavelength vectors

zeropoint

numpy.ndarray

float

Sensitivity function zeropoints

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, sending diagnostic information to the screen.

_algorithm = None

Algorithm used for the sensitivity calculation.

_bundle()[source]

Bundle the object for writing using to_hdu().

compute_throughput()[source]

Compute the spectroscopic throughput

Returns

  • throughput (numpy.ndarray, float, shape is (nspec, norders)) – Throughput measurements

  • throughput_splice (numpy.ndarray, float, shape is (nspec_splice, norders)) – Throughput measurements for spliced spectra

compute_zeropoint()[source]

Dummy method overloaded by subclasses

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'>}, '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 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, 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 on the detector.

Returns

Instance of the empty sensitivity function table.

Return type

astropy.table.Table

eval_zeropoint(wave, iorddet)[source]

Dummy method, overloaded by subclasses

extrapolate(samp_fact=1.5)[source]

Extrapolates the sensitivity function to cover an extra wavelength range set by the extrapl_blu and extrap_red parameters. 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
classmethod get_instance(spec1dfile, sensfile, par, debug=False)[source]

Instantiate the relevant subclass based on the algorithm provided in par.

internals = ['sensfile', 'spectrograph', 'par', 'qafile', 'thrufile', 'debug', 'wave_cnts', 'counts', 'counts_ivar', 'counts_mask', 'nspec_in', 'norderdet', 'wave_splice', 'zeropoint_splice', 'throughput_splice', 'steps', 'splice_multi_det', 'meta_spec', 'std_dict']

A list of strings identifying a set of internal attributes that are not part of the datamodel.

run()[source]

Execute the sensitivity function calculations.

classmethod sensfunc_weights(sensfile, waves, debug=False, extrap_sens=True)[source]

Get the weights based on the sensfunc

Parameters
  • sensfile (str) – the name of your fits format sensfile

  • waves (ndarray) – (nspec, norders, nexp) or (nspec, norders) wavelength grid for your output weights

  • debug (bool) – default=False show the weights QA

Returns

sensfunc weights evaluated on the input waves wavelength grid

Return type

ndarray

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

version = '1.0.1'

Datamodel version.

write_QA()[source]

Write out zeropoint QA files

class pypeit.sensfunc.UVISSensFunc(spec1dfile, sensfile, par=None, debug=False)[source]

Bases: SensFunc

Determine 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.

compute_zeropoint()[source]

Calls routine to compute the sensitivity function.

eval_zeropoint(wave, iorddet)[source]

Evaluate the sensitivity function zero-points

Parameters
  • wave (numpy.ndarray, shape is (nspec)) – Wavelength array

  • iorddet (int) – Order or detector (0-indexed)

Returns

zeropoint

Return type

numpy.ndarray, shape is (nspec,)