pypeit.wavetilts module

Module for guiding Arc/Sky line tracing

class pypeit.wavetilts.BuildWaveTilts(mstilt, slits, spectrograph, par, wavepar, det=1, qa_path=None, spat_flexure=None)[source]

Bases: object

Class to guide arc/sky tracing

Parameters:
  • mstilt (TiltImage) – Tilt image. QA file naming inherits the calibration key (calib_key) from this object.

  • slits (SlitTraceSet) – Slit edges

  • spectrograph (Spectrograph) – Spectrograph object

  • par (WaveTiltsPar or None) – The parameters used to fuss with the tilts

  • wavepar (WavelengthSolutionPar or None) – The parameters used for the wavelength solution

  • det (int) – Detector index

  • qa_path (str, optional) – Directory for QA output.

  • spat_flexure (float, optional) – If input, the slitmask and slit edges are shifted prior to tilt analysis.

spectrograph

??

Type:

Spectrograph

steps

??

Type:

list

mask

boolean array; True = Ignore this slit

Type:

numpy.ndarray

all_trcdict

All trace dict’s

Type:

list

tilts

Tilts for a single slit/order

Type:

numpy.ndarray

all_tilts

Tuple of tilts numpy.ndarray objects

Type:

list

final_tilts

Final tilts image

Type:

numpy.ndarray

gpm

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

Type:

numpy.ndarray

_parse_param(par, key, slit)[source]

Grab a parameter for a given slit

Parameters:
Returns:

Value of the parameter

Return type:

object

extract_arcs()[source]

Extract the arcs down each slit/order

Wrapper to arc.get_censpec()

Returns:

Extracted arcs in two numpy.ndarray objects

Return type:

tuple

find_lines(arcspec, slit_cen, slit_idx, bpm=None, debug=False)[source]

Find the lines for tracing

Wrapper to tracewave.tilts_find_lines()

Parameters:
  • () (slit_cen) –

    ??

  • ()

    ??

  • slit_idx (int) – Slit index, zero-based

  • bpm (numpy.ndarray, optional) –

    ??

  • debug (bool, optional) –

    ??

Returns:

2 objectcs

Return type:

tuple

fit_tilts(trc_tilt_dict, thismask, slit_cen, spat_order, spec_order, slit_idx, show_QA=False, doqa=True)[source]

Fit the tilts

all_fit_dict and all_trace_dict are filled in place

Parameters:
  • trc_tilt_dict (dict) – Contains information from tilt tracing

  • slit_cen (numpy.ndarray) – (nspec,) Central trace for this slit

  • spat_order (int) – Order of the 2d polynomial fit for the spatial direction

  • spec_order (int) – Order of the 2d polytnomial fit for the spectral direction

  • slit_idx (int) – zero-based, integer index for the slit in question

  • show_QA (bool, optional) – show the QA instead of writing it out to the outfile

  • doqa (bool, optional) – Construct the QA plot

Returns:

Array containing the coefficients for the 2d legendre polynomial fit. Shape is (spat_order + 1, spec_order+1).

Return type:

numpy.ndarray

make_tbl_tilt_traces()[source]

Make a Table of the traced and fitted tilts from self.all_trace_dict to save into the WaveTilts object

Returns:

Table including the traced and fitted tilts for each slit. Columns are:

  • slit_ids: Slit IDs for which the tilts were traced and fit

  • goodpix_spat: Good spatial pixels of the traced tilts

  • goodpix_tilt: Good spectral pixels of the traced tilts

  • goodpix_lid: Line IDs for each goodpix. This is needed to associate goodpix to each line

  • badpix_spat: Masked spatial pixels of the traced tilts

  • badpix_tilt: Masked spectral pixels of the traced tilts

  • badpix_lid: Line IDs for each badpix. This is needed to associate badpix to each line

  • good2dfit_spat: Good spatial pixels of the 2D fit of the tilts

  • good2dfit_tilt: Good spectral pixels of the 2D fit of the tilts

  • good2dfit_lid: Line IDs for each good2dfit. This is needed to associate good2dfit to each line

  • bad2dfit_spat: Rejected spatial pixels of the 2D fit of the tilts

  • bad2dfit_tilt: Rejected spectral pixels of the 2D fit of the tilts

  • bad2dfit_lid: Line IDs for each bad2dfit. This is needed to associate bad2dfit to each line

Return type:

astropy.table.Table

model_arc_continuum(debug=False)[source]

Model the continuum of the arc image.

The method uses the arc spectra extracted using extract_arcs and fits a characteristic low-order continuum for each slit/order using robust_fit() and the parameters cont_function, cont_order, and cont_rej from par. The characteristic continuum is then rescaled to match the continuum at each spatial position in the slit/order.

Note

The approach used here may be too simplistic (in the robustness of the continuum fit and then how the continuum is rescaled and projected for each slit/order). Tests should be performed to make sure that the approach is good enough to trace the centroid of the arc/sky lines without biasing the centroids due to a non-zero continuum level.

Parameters:

debug (bool, optional) – Run the method in debug mode.

Returns:

Returns a 2D image with the same shape as mstilt with the model continuum.

Return type:

numpy.ndarray

run(doqa=True, debug=False, show=False)[source]

Main driver for tracing arc lines

Code flow:

  1. Extract an arc spectrum down the center of each slit/order

  2. Loop on slits/orders
    1. Trace and fit the arc lines (This is done twice, once with trace_crude as the tracing crutch, then again with a PCA model fit as the crutch).

    2. Repeat trace.

    3. 2D Fit to the offset from slitcen

    4. Save

Parameters:
Return type:

WaveTilts

trace_tilts(arcimg, lines_spec, lines_spat, thismask, slit_cen, debug_pca=False, show_tracefits=False)[source]

Trace the tilts

Parameters:
  • arcimg (numpy.ndarray) – Arc image. Shape is (nspec, nspat).

  • lines_spec (numpy.ndarray) – Array containing the spectral pixel location of each line found for this slit. Shape is (nlines,).

  • lines_spat (numpy.ndarray) – Array containing the spatial pixel location of each line, which is the slitcen evaluate at the spectral position position of the line stored in lines_spec. Shape is (nlines,).

  • thismask (numpy.ndarray) – Image indicating which pixels lie on the slit in equation. True = on the slit. False = not on slit. Shape is (nspec, nspat) with dtype=bool.

  • slit_cen (int) – Integer index indicating the slit in question.

Returns:

Dictionary containing information on the traced tilts required to fit the filts.

Return type:

dict

class pypeit.wavetilts.WaveTilts(coeffs, nslit, spat_id, spat_order, spec_order, func2d, bpmtilts=None, spat_flexure=None, PYP_SPEC=None, slits_filename=None, tiltimg_filename=None, tilt_traces=None)[source]

Bases: CalibFrame

Calibration frame containing the wavelength tilt 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.2.0

Attribute

Type

Array Type

Description

PYP_SPEC

str

PypeIt spectrograph name

bpmtilts

numpy.ndarray

numpy.integer

Bad pixel mask for tilt solutions. Keys are taken from SlitTraceSetBitmask

coeffs

numpy.ndarray

numpy.floating

2D coefficents for the fit on the initial slits. One set per slit/order (3D array).

func2d

str

Function used for the 2D fit

nslit

int

Total number of slits. This can include masked slits

slits_filename

str

Path to SlitTraceSet file. This helps to find the Slits calibration file when running pypeit_chk_tilts()

spat_flexure

float

Flexure shift from the input TiltImage

spat_id

numpy.ndarray

numpy.integer

Slit spat_id

spat_order

numpy.ndarray

numpy.integer

Order for spatial fit (nslit)

spec_order

numpy.ndarray

numpy.integer

Order for spectral fit (nslit)

tilt_traces

astropy.table.table.Table

Table with the positions of the traced and fitted tilts for all the slits. see make_tbl_tilt_traces() for more details.

tiltimg_filename

str

Path to Tiltimg file. This helps to find Tiltimg file when running pypeit_chk_tilts()

When written to an output-file HDU, all numpy.ndarray elements are bundled into an astropy.io.fits.BinTableHDU, and the other elements are written as header keywords. Any datamodel elements that are None are not included in the output.

_bundle()[source]

Bundle the data in preparation for writing to a fits file.

See _bundle(). Data is always written to a ‘TILTS’ extension.

calib_file_format = 'fits'

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

calib_type = 'Tilts'

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

datamodel = {'PYP_SPEC': {'descr': 'PypeIt spectrograph name', 'otype': <class 'str'>}, 'bpmtilts': {'atype': <class 'numpy.integer'>, 'descr': 'Bad pixel mask for tilt solutions. Keys are taken from SlitTraceSetBitmask', 'otype': <class 'numpy.ndarray'>}, 'coeffs': {'atype': <class 'numpy.floating'>, 'descr': '2D coefficents for the fit on the initial slits.  One set per slit/order (3D array).', 'otype': <class 'numpy.ndarray'>}, 'func2d': {'descr': 'Function used for the 2D fit', 'otype': <class 'str'>}, 'nslit': {'descr': 'Total number of slits.  This can include masked slits', 'otype': <class 'int'>}, 'slits_filename': {'descr': 'Path to SlitTraceSet file. This helps to find the Slits calibration file when running pypeit_chk_tilts()', 'otype': <class 'str'>}, 'spat_flexure': {'descr': 'Flexure shift from the input TiltImage', 'otype': <class 'float'>}, 'spat_id': {'atype': <class 'numpy.integer'>, 'descr': 'Slit spat_id ', 'otype': <class 'numpy.ndarray'>}, 'spat_order': {'atype': <class 'numpy.integer'>, 'descr': 'Order for spatial fit (nslit)', 'otype': <class 'numpy.ndarray'>}, 'spec_order': {'atype': <class 'numpy.integer'>, 'descr': 'Order for spectral fit (nslit)', 'otype': <class 'numpy.ndarray'>}, 'tilt_traces': {'descr': 'Table with the positions of the traced and fitted tilts for all the slits. see :func:`~pypeit.wavetilts.BuildWaveTilts.make_tbl_tilt_traces` for more details. ', 'otype': <class 'astropy.table.table.Table'>}, 'tiltimg_filename': {'descr': 'Path to Tiltimg file. This helps to find Tiltimg file when running pypeit_chk_tilts()', 'otype': <class 'str'>}}

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, ...}
fit2tiltimg(slitmask, flexure=None)[source]

Generate a tilt image from the fit parameters

Mainly to allow for flexure

Parameters:
  • slitmask (numpy.ndarray) –

    ??

  • flexure (float, optional) – Spatial shift of the tilt image onto the desired frame (typically a science image)

Returns:

New tilt image

Return type:

numpy.ndarray

is_synced(slits)[source]

Confirm the slits in WaveTilts are aligned to that in SlitTraceSet

Barfs if not

Parameters:

slits (SlitTraceSet)

show(waveimg=None, wcs_match=True, in_ginga=True, show_traces=False, chk_version=True)[source]

Show in ginga or mpl Tiltimg with the tilts traced and fitted overlaid

Parameters:
  • waveimg (numpy.ndarray, optional) – Image with the wavelength solution.

  • wcs_match (bool, optional) – If True, use this image as a reference image for the WCS and match all image in other channels to it.

  • in_ginga (bool, optional) – If True, show the image in ginga. Otherwise, use matplotlib.

  • show_traces (bool, optional) – If True, show the traces of the tilts on the image.

  • 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!

spatid_to_zero(spat_id)[source]

Convert slit spat_id to zero-based Mainly for coeffs

Parameters:

spat_id (int)

Return type:

int

version = '1.2.0'

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.wavetilts.show_tilts_mpl(tilt_img, tilt_traces, show_traces=False, left_edges=None, right_edges=None, slit_ids=None, cut=None)[source]

Show the TiltImage with the traced and 2D fitted tilts overlaid

Parameters:
  • tilt_img (numpy.ndarray) – TiltImage

  • tilt_traces (astropy.table.Table) – Table containing the traced and fitted tilts. See make_tbl_tilt_traces() for information on the table columns.

  • show_traces (bool, optional) – Show the traced tilts

  • left_edges (numpy.ndarray, optional) – Left edges of the slits

  • right_edges (numpy.ndarray, optional) – Right edges of the slits

  • slit_ids (numpy.ndarray, optional) – Slit IDs

  • cut (tuple, optional) – Lower and upper levels cut for the image display