pypeit.core.wavecal.wv_fitting module

Module for finding patterns in arc line spectra

class pypeit.core.wavecal.wv_fitting.WaveFit(spat_id, ech_order=None, pypeitfit=None, pixel_fit=None, wave_fit=None, ion_bits=None, cen_wave=None, cen_disp=None, spec=None, wave_soln=None, sigrej=None, shift=None, tcent=None, rms=None, xnorm=None, fwhm=None)[source]

Bases: DataContainer

DataContainer for the output from BuildWaveCalib

All of the items in the datamodel are required for instantiation, although they can be None (but shouldn’t be).

The datamodel attributes are:

Version: 1.1.1

Attribute

Type

Array Type

Description

cen_disp

float

Approximate wavelength dispersion

cen_wave

float

Central wavelength

ech_order

int, numpy.integer

Echelle order number.

fwhm

float

Estimate FWHM of arc lines in binned pixels of the input arc frame

ion_bits

numpy.ndarray

numpy.integer

Ion bit values for the Ion names

pixel_fit

numpy.ndarray

numpy.floating

Pixel values of arc lines

pypeitfit

PypeItFit

Fit to 1D wavelength solutions

rms

float

RMS of the solution

shift

float

Shift applied

sigrej

float

Final sigma rejection applied

spat_id

int, numpy.integer

Spatial position of slit/order for this fit. Required for I/O

spec

numpy.ndarray

numpy.floating

Arc spectrum

tcent

numpy.ndarray

numpy.floating

Pixel centroids of all arc lines found

wave_fit

numpy.ndarray

numpy.floating

Wavelength IDs assigned

wave_soln

numpy.ndarray

numpy.floating

Evaluated wavelengths at pixel_fit

xnorm

float

Normalization for fit

When written to an output-file HDU, all numpy.ndarray elements are bundled into an astropy.io.fits.BinTableHDU, the pypeitfit attribute is written to a separate extension (see PypeItFit), and the other elements are written as header keywords. Any datamodel elements that are None are not included in the output. The two HDU extensions are given names according to their spatial ID; see hduext_prefix_from_spatid().

_bundle(**kwargs)[source]

Over-ride DataContainer._bundle() to deal with PYPEITFIT

Parameters:

kwargs – Passed to DataContainer._bundle()

Return type:

list

bitmask = <pypeit.core.wavecal.defs.LinesBitMask object>
datamodel = {'cen_disp': {'descr': 'Approximate wavelength dispersion', 'otype': <class 'float'>}, 'cen_wave': {'descr': 'Central wavelength', 'otype': <class 'float'>}, 'ech_order': {'descr': 'Echelle order number.', 'otype': (<class 'int'>, <class 'numpy.integer'>)}, 'fwhm': {'descr': 'Estimate FWHM of arc lines in binned pixels of the input arc frame', 'otype': <class 'float'>}, 'ion_bits': {'atype': <class 'numpy.integer'>, 'descr': 'Ion bit values for the Ion names', 'otype': <class 'numpy.ndarray'>}, 'pixel_fit': {'atype': <class 'numpy.floating'>, 'descr': 'Pixel values of arc lines', 'otype': <class 'numpy.ndarray'>}, 'pypeitfit': {'descr': 'Fit to 1D wavelength solutions', 'otype': <class 'pypeit.core.fitting.PypeItFit'>}, 'rms': {'descr': 'RMS of the solution', 'otype': <class 'float'>}, 'shift': {'descr': 'Shift applied', 'otype': <class 'float'>}, 'sigrej': {'descr': 'Final sigma rejection applied', 'otype': <class 'float'>}, 'spat_id': {'descr': 'Spatial position of slit/order for this fit. Required for I/O', 'otype': (<class 'int'>, <class 'numpy.integer'>)}, 'spec': {'atype': <class 'numpy.floating'>, 'descr': 'Arc spectrum', 'otype': <class 'numpy.ndarray'>}, 'tcent': {'atype': <class 'numpy.floating'>, 'descr': 'Pixel centroids of all arc lines found', 'otype': <class 'numpy.ndarray'>}, 'wave_fit': {'atype': <class 'numpy.floating'>, 'descr': 'Wavelength IDs assigned', 'otype': <class 'numpy.ndarray'>}, 'wave_soln': {'atype': <class 'numpy.floating'>, 'descr': 'Evaluated wavelengths at pixel_fit', 'otype': <class 'numpy.ndarray'>}, 'xnorm': {'descr': 'Normalization for fit', 'otype': <class 'float'>}}

Provides the class data model. This is undefined in the abstract class and should be overwritten in the derived classes.

The format of the datamodel needed for each implementation of a DataContainer derived class is as follows.

The datamodel itself is a class attribute (i.e., it is a member of the class, not just of an instance of the class). The datamodel is a dictionary of dictionaries: Each key of the datamodel dictionary provides the name of a given datamodel element, and the associated item (dictionary) for the datamodel element provides the type and description information for that datamodel element. For each datamodel element, the dictionary item must provide:

  • otype: This is the type of the object for this datamodel item. E.g., for a float or a numpy.ndarray, you would set otype=float and otype=np.ndarray, respectively.

  • descr: This provides a text description of the datamodel element. This is used to construct the datamodel tables in the pypeit documentation.

If the object type is a numpy.ndarray, you should also provide the atype keyword that sets the type of the data contained within the array. E.g., for a floating point array containing an image, your datamodel could be simply:

datamodel = {'image' : dict(otype=np.ndarray, atype=float, descr='My image')}

More advanced examples are given in the top-level module documentation.

Currently, datamodel components are restricted to have otype that are tuple, int, float, numpy.integer, numpy.floating, numpy.ndarray, or astropy.table.Table objects. E.g., datamodel values for otype cannot be dict.

classmethod from_hdu(hdu, hdu_prefix=None, spat_id=None, **kwargs)[source]

Parse the data from one or more fits objects.

Note

spat_id and hdu_prefix are mutually exclusive, with hdu_prefix taking precedence. If both are None, the code attempts to read the spatial ID number from the HDU extension name(s) to construct the appropriate prefix. In the latter case, only the first WaveFit object will be read if the HDU object contains many.

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) – Only parse HDUs with extension names matched to this prefix. If None, spat_id will be used to construct the prefix if it is provided, otherwise all extensions are parsed.

  • spat_id (int, optional) – Spatial ID number for the relevant slit/order for this wavelength fit. If provided, this is used to construct the hdu_prefix using hduext_prefix_from_spatid().

  • kwargs – Passed directly to the base-class from_hdu.

Returns:

Object instantiated from data in the provided HDU.

Return type:

WaveFit

static hduext_prefix_from_spatid(spat_id)[source]

Use the slit spatial ID number to construct the prefix for an output HDU.

Parameters:

spat_id (int) – Slit/order spatial ID number

Returns:

The prefix to use for HDU extensions associated with this data container.

Return type:

str

property ions

Returns an array of ion labels

Returns:

Array of the ion label for each line as recorded in ion_bits

Return type:

numpy.ndarray

static parse_spatid_from_hduext(ext)[source]

Parse the slit spatial ID integer from a full HDU extension name.

Parameters:

ext (str) – Full extension name.

Returns:

The parsed slit spatial ID number. Returns None if the extension name does not start with 'SPAT_ID-'.

Return type:

int

to_hdu(**kwargs)[source]

Over-ride base-class function to force force_to_bintbl to always be true.

See base-class to_hdu for arguments.

Parameters:

kwargs – Passed directly to the base-class to_hdu. If force_to_bintbl is present, it is forced to be True.

Returns:

A list of HDUs, where the type depends on the value of add_primary.

Return type:

list, astropy.io.fits.HDUList

version = '1.1.1'

Provides the string representation of the class version.

This is currently put to minimal use so far, but will used for I/O verification in the future.

Each derived class should provide a version to guard against data model changes during development.

pypeit.core.wavecal.wv_fitting.fit_slit(spec, patt_dict, tcent, line_lists, vel_tol=1.0, outroot=None, slittxt='Slit', thar=False, match_toler=3.0, func='legendre', n_first=2, sigrej_first=2.0, n_final=4, sigrej_final=3.0, verbose=False)[source]

Perform a fit to the wavelength solution. Wrapper for iterative fitting code.

Parameters:
  • spec (ndarray) – arc spectrum

  • patt_dict (dict) – dictionary of patterns

  • tcent (ndarray) – List of the detections in this slit to be fit using the patt_dict

  • line_lists (astropy.table.Table) – Table containing the line list

  • vel_tol (float, default = 1.0) – Tolerance in km/s for matching lines in the IDs to lines in the NIST database. The default is 1.0 km/s

  • outroot (str) – Path for QA file.

  • slittxt (str) – Label used for QA

  • thar (bool, default = False) – True if this is a ThAr fit

  • match_toler (float, default = 3.0) – Matching tolerance when searching for new lines. This is the difference in pixels between the wavlength assigned to an arc line by an iteration of the wavelength solution to the wavelength in the line list.

  • func (str, default = 'legendre') – Name of function used for the wavelength solution

  • n_first (int, default = 2) – Order of first guess to the wavelength solution.

  • sigrej_first (float, default = 2.0) – Number of sigma for rejection for the first guess to the wavelength solution.

  • n_final (int, default = 4) – Order of the final wavelength solution fit

  • sigrej_final (float, default = 3.0) – Number of sigma for rejection for the final fit to the wavelength solution.

  • verbose (bool) – If True, print out more information.

  • plot_fil – Filename for plotting some QA?

Returns:

final_fit – A dictionary containing all of the information about the fit

Return type:

dict

pypeit.core.wavecal.wv_fitting.iterative_fitting(spec, tcent, ifit, IDs, llist, dispersion, match_toler=2.0, func='legendre', n_first=2, sigrej_first=2.0, n_final=4, sigrej_final=3.0, input_only=False, weights=None, plot_fil=None, verbose=False)[source]

Routine for iteratively fitting wavelength solutions.

Parameters:
  • spec (ndarray, shape = (nspec,)) – arcline spectrum

  • tcent (ndarray) – Centroids in pixels of lines identified in spec

  • ifit (ndarray) – Indices of the lines that will be fit

  • IDs (ndarray) – wavelength IDs of the lines that will be fit (I think?)

  • llist (dict) – Linelist dictionary

  • dispersion (float) – dispersion

  • match_toler (float, default = 3.0) – Matching tolerance when searching for new lines. This is the difference in pixels between the wavlength assigned to an arc line by an iteration of the wavelength solution to the wavelength in the line list.

  • func (str, default = 'legendre') – Name of function used for the wavelength solution

  • n_first (int, default = 2) – Order of first guess to the wavelength solution.

  • sigrej_first (float, default = 2.0) – Number of sigma for rejection for the first guess to the wavelength solution.

  • n_final (int, default = 4) – Order of the final wavelength solution fit

  • sigrej_final (float, default = 3.0) – Number of sigma for rejection for the final fit to the wavelength solution.

  • input_only (bool) – If True, the routine will only perform a robust polyfit to the input IDs. If False, the routine will fit the input IDs, and then include additional lines in the linelist that are a satisfactory fit.

  • weights (ndarray) – Weights to be used?

  • verbose (bool) – If True, print out more information.

  • plot_fil – Filename for plotting some QA?

Returns:

final_fit – Fit result

Return type:

WaveFit