pypeit.wavecalib module

Module for guiding 1D Wavelength Calibration

class pypeit.wavecalib.BuildWaveCalib(msarc, slits, spectrograph, par, lamps, meta_dict=None, det=1, qa_path=None, msbpm=None)[source]

Bases: object

Class to guide wavelength calibration

Parameters:
  • msarc (PypeItImage) – Arc image, created by the ArcImage class

  • slits (SlitTraceSet) – Slit edges

  • spectrograph (Spectrograph) – The Spectrograph instance that sets the instrument used to take the observations. Used to set spectrograph.

  • par (WavelengthSolutionPar) – The parameters used for the wavelength solution. Uses ['calibrations']['wavelengths'].

  • meta_dict (dict, optional) –

    Dictionary containing meta information required for wavelength calibration. Specifically for non-fixed format echelles this dict must contain the following keys:

    • 'echangle': the echelle angle

    • 'xdangle': the cross-disperser angle

    • 'dispmame': the disperser name

  • det (int, optional) – Detector number

  • msbpm (numpy.ndarray, optional) – Bad pixel mask image

  • qa_path (str, optional) – For QA

steps

List of the processing steps performed

Type:

list

wv_calib

Primary output. Keys 0, 1, 2, 3 are solution for individual previously slit steps

Type:

dict

arccen

(nwave, nslit) Extracted arc(s) down the center of the slit(s)

Type:

numpy.ndarray

maskslits

Slits to ignore because they were not extracted. WARNING: Outside of this Class, it is best to regenerate the mask using make_maskslits()

Type:

numpy.ndarray

gpm

Good pixel mask Eventually, we might attach this to self.msarc although that would then require that we write it to disk with self.msarc.image

Type:

numpy.ndarray

nonlinear_counts

Specifies saturation level for the arc lines

Type:

float

wvc_bpm

Mask for slits attempted to have a wv_calib solution

Type:

numpy.ndarray

build_wv_calib(arccen, method, skip_QA=False, prev_wvcalib=None)[source]

Main routine to generate the wavelength solutions in a loop over slits Wrapper to arc.simple_calib or arc.calib_with_arclines

self.maskslits is updated for slits that fail

Parameters:
  • method – str ‘simple’ – arc.simple_calib ‘arclines’ – arc.calib_with_arclines ‘holy-grail’ – wavecal.autoid.HolyGrail ‘reidentify’ – wavecal.auotid.ArchiveReid ‘identify’ – wavecal.identify.Identify ‘full_template’ – wavecal.auotid.full_template

  • skip_QA (bool, optional)

  • prev_wvcalib (WaveCalib, optional) – Previous wavelength calibration

Returns:

self.wv_calib

Return type:

dict

echelle_2dfit(wv_calib, debug=False, skip_QA=False)[source]

Fit a two-dimensional wavelength solution for echelle data.

Primarily a wrapper for pypeit.core.arc.fit2darc(), using data unpacked from the wv_calib dictionary.

Parameters:
  • wv_calib (pypeit.wavecalib.WaveCalib) – Wavelength calibration object

  • debug (bool, optional) – Show debugging info

  • skip_QA (bool, optional) – Flag to skip construction of the nominal QA plots.

Returns:

  • fit2ds (list of pypeit.fitting.PypeItFit) – Contains information from 2-d fit. Frequently a list of 1 fit. The main exception is for a mosaic when one sets echelle_separate_2d=True.

  • dets (list) – List of integers for the detector numbers.

  • save_order_dets (list) – List of integer lists providing list of the orders.

extract_arcs(slitIDs=None)[source]

Extract the arcs down each slit/order

Wrapper to arc.get_censpec()

Parameters:

slitIDs (list, optional) – A list of the slit IDs to extract (if None, all slits will be extracted)

Returns:

  • arccen (numpy.ndarray) – arc spectrum for all slits, shape is (nspec, nslit):

  • wvc_bpm (numpy.ndarray) – boolean array containing a mask indicating which slits are good. True = masked (bad).

frametype = 'wv_calib'
redo_echelle_orders(bad_orders: ndarray, dets: ndarray, order_dets: ndarray, bad_orders_maxfrac: float = 0.1, frac_rms_thresh: float = 1.1)[source]

Attempt to redo the wavelength calibration for a set of bad echelle orders

Parameters:
  • bad_orders (np.ndarray) – Array of bad order numbers

  • dets (np.ndarray) – detectors of the spectrograph Multiple numbers for mosaic (typically)

  • order_dets (np.ndarray) – Orders on each detector

  • bad_orders_maxfrac (float) – Maximum fraction of bad orders in a detector for attempting a refit

  • frac_rms_thresh (float) – Fractional change in the RMS threshold for accepting a refit

Returns:

True if any of the echelle orders were successfully redone

Return type:

bool

run(skip_QA=False, debug=False, prev_wvcalib=None)[source]

Main method for wavelength calibration

Code flow:
  1. Extract 1D arc spectra down the center of each unmasked slit/order

  2. Generate the 1D wavelength fits

  3. If echelle, perform 2D fit(s).

  4. Deal with masking

  5. Return a WaveCalib object

Parameters:
  • skip_QA (bool, optional) – Skip QA?

  • prev_wvcalib (WaveCalib, optional) – Previous wavelength calibration object (from disk, typically)

Returns:

wavelength calibration object

Return type:

WaveCalib

show(item, slit=None)[source]

Show one of the class internals

Parameters:
  • item (str) – ‘spec’ – Show the fitted points and solution; requires slit ‘fit’ – Show fit QA; requires slit

  • slit (int, optional)

Returns:

update_wvmask()[source]

(re)Generate the mask for wv_calib based on its contents This is the safest way to go…

Parameters:

nslit (int) – Number of slits/orders

Returns:

self.wvc_bpm, boolean array – True = masked, i.e. do not use

Return type:

numpy.ndarray

class pypeit.wavecalib.WaveCalib(wv_fits=None, fwhm_map=None, nslits=None, spat_ids=None, ech_orders=None, PYP_SPEC=None, strpar=None, wv_fit2d=None, arc_spectra=None, lamps=None, det_img=None)[source]

Bases: CalibFrame

Calibration frame containing the wavelength calibration.

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

Attribute

Type

Array Type

Description

PYP_SPEC

str

PypeIt spectrograph name

arc_spectra

numpy.ndarray

numpy.floating

2D array: 1D extracted spectra, slit by slit (nspec, nslits)

det_img

numpy.ndarray

numpy.integer

Detector image which indicates which pixel in the mosaic corresponds to which detector; used occasionally by echelle. Currently only saved if ech_separate_2d=True

ech_orders

numpy.ndarray

numpy.integer

Echelle order ID numbers. Defined only for echelle.

fwhm_map

numpy.ndarray

PypeItFit

A fit that determines the spectral FWHM at every location of every slit

lamps

str

List of arc lamps used for the wavelength calibration

nslits

int

Total number of slits. This can include masked slits

spat_ids

numpy.ndarray

numpy.integer

Slit spat_ids. Named distinctly from that in WaveFit

strpar

str

Parameters as a string

wv_fit2d

numpy.ndarray

PypeItFit

2D wavelength solution(s) (echelle). If there is more than one, they must be aligned to the separate detectors analyzed

wv_fits

numpy.ndarray

WaveFit

WaveFit to each 1D wavelength solution

_bundle()[source]

Override base class function to write one HDU per image. Any extras are in the HDU header of the primary image.

Returns:

A list of dictionaries, each list element is written to its own fits extension.

Return type:

list

build_fwhmimg(tilts, slits, initial=False, spat_flexure=None)[source]

Generates an image of the instrument spectral FWHM (units=pixels) at every pixel on the detector.

Parameters:
  • tilts (numpy.ndarray) – Image holding tilts

  • slits (pypeit.slittrace.SlitTraceSet) – Properties of the slits

  • initial (bool, optional) – If True, the initial slit locations will be used. Otherwise, the tweaked edges will be used.

  • spat_flexure (float, optional) – Spatial flexure correction in pixels.

Returns:

The spectral FWHM image.

Return type:

numpy.ndarray

build_waveimg(tilts, slits, spat_flexure=None, spec_flexure=None)[source]

Main algorithm to build the wavelength image

Only applied to good slits, which means any non-flagged or flagged

in the exclude_for_reducing list

Parameters:
  • tilts (numpy.ndarray) – Image holding tilts

  • slits (pypeit.slittrace.SlitTraceSet) – Properties of the slits

  • spat_flexure (float, optional) – Spatial flexure correction in pixels.

  • spec_flexure (float, numpy.ndarray, optional) – Spectral flexure correction in pixels. If a float, the same spectral flexure correction will be applied to all slits. If a numpy array, the length of the array should be the same as the number of slits. The value of each element is the spectral shift in pixels to be applied to each slit.

Returns:

The wavelength image.

Return type:

numpy.ndarray

calib_file_format = 'fits'

The extension and file format of the output file. Should be 'fits' or 'fits.gz' (for gzipped output).

calib_type = 'WaveCalib'

The type of the calibration frame, primarily used to set the name of the output file.

chk_synced(slits)[source]

Confirm the slits in WaveCalib are aligned to that in SlitTraceSet

Barfs if not

Parameters:

slits (pypeit.slittrace.SlitTraceSet)

datamodel = {'PYP_SPEC': {'descr': 'PypeIt spectrograph name', 'otype': <class 'str'>}, 'arc_spectra': {'atype': <class 'numpy.floating'>, 'descr': '2D array: 1D extracted spectra, slit by slit (nspec, nslits)', 'otype': <class 'numpy.ndarray'>}, 'det_img': {'atype': <class 'numpy.integer'>, 'descr': 'Detector image which indicates which pixel in the mosaic corresponds to which detector; used occasionally by echelle.  Currently only saved if ech_separate_2d=True', 'otype': <class 'numpy.ndarray'>}, 'ech_orders': {'atype': <class 'numpy.integer'>, 'descr': 'Echelle order ID numbers.  Defined only for echelle.', 'otype': <class 'numpy.ndarray'>}, 'fwhm_map': {'atype': <class 'pypeit.core.fitting.PypeItFit'>, 'descr': 'A fit that determines the spectral FWHM at every location of every slit', 'otype': <class 'numpy.ndarray'>}, 'lamps': {'descr': 'List of arc lamps used for the wavelength calibration', 'otype': <class 'str'>}, 'nslits': {'descr': 'Total number of slits.  This can include masked slits', 'otype': <class 'int'>}, 'spat_ids': {'atype': <class 'numpy.integer'>, 'descr': 'Slit spat_ids. Named distinctly from that in WaveFit ', 'otype': <class 'numpy.ndarray'>}, 'strpar': {'descr': 'Parameters as a string', 'otype': <class 'str'>}, 'wv_fit2d': {'atype': <class 'pypeit.core.fitting.PypeItFit'>, 'descr': '2D wavelength solution(s) (echelle).  If there is more than one, they must be aligned to the separate detectors analyzed', 'otype': <class 'numpy.ndarray'>}, 'wv_fits': {'atype': <class 'pypeit.core.wavecal.wv_fitting.WaveFit'>, 'descr': 'WaveFit to each 1D wavelength solution', 'otype': <class 'numpy.ndarray'>}}

Default datamodel for any CalibFrame. Derived classes should instantiate their datamodels by first inheriting from the base class. E.g.:

class ArcFrame(CalibFrame):
    datamodel = {**CalibFrame.datamodel, ...}
classmethod from_hdu(hdu, chk_version=True, **kwargs)[source]

Instantiate the object from an HDU extension.

This overrides the base-class method. Overriding this method is preferrable to overriding the _parse method because it makes it easier to deal with the DataContainer nesting of this object.

Parameters:
property par
version = '1.1.2'

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.

wave_diagnostics(print_diag=False)[source]

Create a table with wavecalib diagnostics

Parameters:

print_diag (bool, optional) – If True, the diagnostic table is printed to screen

Returns:

wavecalib diagnostics table

Return type:

astropy.table.Table