pypeit.slittrace module
Implements the objects used to hold slit edge data.
- class pypeit.slittrace.SlitTraceBitMask[source]
Bases:
BitMask
Mask bits used during slit tracing.
- property exclude_for_flexure
- property exclude_for_reducing
- version = '1.0.1'
- class pypeit.slittrace.SlitTraceSet(left_init, right_init, pypeline, detname=None, nspec=None, nspat=None, PYP_SPEC=None, mask_init=None, specmin=None, specmax=None, binspec=1, binspat=1, pad=0, spat_id=None, maskdef_id=None, maskdef_designtab=None, maskfile=None, maskdef_posx_pa=None, maskdef_offset=None, maskdef_objpos=None, maskdef_slitcen=None, ech_order=None, nslits=None, left_tweak=None, right_tweak=None, center=None, mask=None)[source]
Bases:
CalibFrame
Defines a generic class for holding and manipulating image traces organized into left-right slit pairs.
The datamodel attributes are:
Version: 1.1.5
Attribute
Type
Array Type
Description
PYP_SPEC
str
PypeIt spectrograph name
binspat
int
Number of pixels binned in the spatial direction.
binspec
int
Number of pixels binned in the spectral direction.
center
Spatial coordinates of the slit centers from left_init, right_init. Shape is Nspec by Nslits.
detname
str
Identifier for detector or mosaic
ech_order
int, numpy.integer
Slit ID number echelle order
left_init
Spatial coordinates (pixel indices) of all left edges, one per slit. Derived from the TraceImage. Shape is Nspec by Nslits.
left_tweak
Spatial coordinates (pixel indices) of all left edges, one per slit. These traces have been adjusted by the flat-field. Shape is Nspec by Nslits.
mask
Bit mask for slits (fully good slits have 0 value). Shape is Nslits.
mask_init
Bit mask for slits at instantiation. Used to reset
maskdef_designtab
Table with slitmask design and object info
maskdef_id
int, numpy.integer
Slit ID number slitmask
maskdef_objpos
Object positions expected by the slitmask design [relative pixels]
maskdef_offset
float
Slitmask offset (pixels) from position expected by the slitmask design
maskdef_posx_pa
float
PA that aligns with spatial dimension of the detector
maskdef_slitcen
Slit centers expected by the slitmask design
maskfile
str
Data file that yielded the slitmask info
nslits
int
Total number of slits, derived from shape of left_init.
nspat
int
Number of pixels in the image spatial direction.
nspec
int
Number of pixels in the image spectral direction.
pad
int
Integer number of pixels to consider beyond the slit edges.
pypeline
str
PypeIt pypeline name
right_init
Spatial coordinates (pixel indices) of all right edges, one per slit. Derived from the TraceImage. Shape is Nspec by Nslits.
right_tweak
Spatial coordinates (pixel indices) of all right edges, one per slit. These traces have been adjusted by the flat-field. Shape is Nspec by Nslits.
spat_id
int, numpy.integer
Slit ID number from SPAT measured at half way point.
specmax
Maximum spectral position (pixel units) allowed for each slit/order. Shape is Nslits.
specmin
Minimum spectral position (pixel units) allowed for each slit/order. Shape is Nslits.
- left_flexure
Convenient spot to hold flexure corrected left
- Type:
- right_flexure
Convenient spot to hold flexure corrected right
- Type:
- _base_header(hdr=None)[source]
Construct the baseline header for all HDU extensions.
This appends the
SlitTraceBitMask
data to the supplied header.- Parameters:
hdr (astropy.io.fits.Header, optional) – Baseline header for additional data. If None, set by
pypeit.io.initialize_header()
.- Returns:
Header object to include in all HDU extensions.
- Return type:
- _bundle()[source]
Bundle the data in preparation for writing to a fits file.
See
pypeit.datamodel.DataContainer._bundle()
. Data is always written to a ‘SLITS’ extension.
- classmethod _parse(hdu, hdu_prefix=None, **kwargs)[source]
Parse the data that was previously written to a fits file.
See
pypeit.datamodel.DataContainer._parse()
. Data is always read from the ‘SLITS’ extension.
- assign_maskinfo(sobjs, plate_scale, spat_flexure, TOLER=1.0)[source]
Assign RA, DEC, Name to objects. Modified in place.
- Parameters:
sobjs (
pypeit.specobjs.SpecObjs
) – List of SpecObj that have been found and traced.plate_scale (
float
) – platescale for the current detector.spat_flexure (
float
) – Shifts, in spatial pixels, between this image and SlitTrace.det_buffer (
int
) – Minimum separation between detector edges and a slit edge.TOLER (
float
, optional) – Matching tolerance in arcsec.
- Returns:
Updated list of SpecObj that have been found and traced.
- Return type:
- Returns:
Updated list of SpecObj that have been found and traced
- Return type:
- bitmask = <pypeit.slittrace.SlitTraceBitMask object>
Bit interpreter for slit masks.
- calib_file_format = 'fits.gz'
File format for the calibration frame file.
- calib_type = 'Slits'
Name for type of calibration frame.
- datamodel = {'PYP_SPEC': {'descr': 'PypeIt spectrograph name', 'otype': <class 'str'>}, 'binspat': {'descr': 'Number of pixels binned in the spatial direction.', 'otype': <class 'int'>}, 'binspec': {'descr': 'Number of pixels binned in the spectral direction.', 'otype': <class 'int'>}, 'center': {'atype': <class 'numpy.floating'>, 'descr': 'Spatial coordinates of the slit centers from left_init, right_init. Shape is Nspec by Nslits.', 'otype': <class 'numpy.ndarray'>}, 'detname': {'descr': 'Identifier for detector or mosaic', 'otype': <class 'str'>}, 'ech_order': {'atype': (<class 'int'>, <class 'numpy.integer'>), 'descr': 'Slit ID number echelle order', 'otype': <class 'numpy.ndarray'>}, 'left_init': {'atype': <class 'numpy.floating'>, 'descr': 'Spatial coordinates (pixel indices) of all left edges, one per slit. Derived from the TraceImage. Shape is Nspec by Nslits.', 'otype': <class 'numpy.ndarray'>}, 'left_tweak': {'atype': <class 'numpy.floating'>, 'descr': 'Spatial coordinates (pixel indices) of all left edges, one per slit. These traces have been adjusted by the flat-field. Shape is Nspec by Nslits.', 'otype': <class 'numpy.ndarray'>}, 'mask': {'atype': <class 'numpy.integer'>, 'descr': 'Bit mask for slits (fully good slits have 0 value). Shape is Nslits.', 'otype': <class 'numpy.ndarray'>}, 'mask_init': {'atype': <class 'numpy.integer'>, 'descr': 'Bit mask for slits at instantiation. Used to reset', 'otype': <class 'numpy.ndarray'>}, 'maskdef_designtab': {'descr': 'Table with slitmask design and object info', 'otype': <class 'astropy.table.table.Table'>}, 'maskdef_id': {'atype': (<class 'int'>, <class 'numpy.integer'>), 'descr': 'Slit ID number slitmask', 'otype': <class 'numpy.ndarray'>}, 'maskdef_objpos': {'atype': <class 'numpy.floating'>, 'descr': 'Object positions expected by the slitmask design [relative pixels]', 'otype': <class 'numpy.ndarray'>}, 'maskdef_offset': {'descr': 'Slitmask offset (pixels) from position expected by the slitmask design', 'otype': <class 'float'>}, 'maskdef_posx_pa': {'descr': 'PA that aligns with spatial dimension of the detector', 'otype': <class 'float'>}, 'maskdef_slitcen': {'atype': <class 'numpy.floating'>, 'descr': 'Slit centers expected by the slitmask design', 'otype': <class 'numpy.ndarray'>}, 'maskfile': {'descr': 'Data file that yielded the slitmask info', 'otype': <class 'str'>}, 'nslits': {'descr': 'Total number of slits, derived from shape of left_init.', 'otype': <class 'int'>}, 'nspat': {'descr': 'Number of pixels in the image spatial direction.', 'otype': <class 'int'>}, 'nspec': {'descr': 'Number of pixels in the image spectral direction.', 'otype': <class 'int'>}, 'pad': {'descr': 'Integer number of pixels to consider beyond the slit edges.', 'otype': <class 'int'>}, 'pypeline': {'descr': 'PypeIt pypeline name', 'otype': <class 'str'>}, 'right_init': {'atype': <class 'numpy.floating'>, 'descr': 'Spatial coordinates (pixel indices) of all right edges, one per slit. Derived from the TraceImage. Shape is Nspec by Nslits.', 'otype': <class 'numpy.ndarray'>}, 'right_tweak': {'atype': <class 'numpy.floating'>, 'descr': 'Spatial coordinates (pixel indices) of all right edges, one per slit. These traces have been adjusted by the flat-field. Shape is Nspec by Nslits.', 'otype': <class 'numpy.ndarray'>}, 'spat_id': {'atype': (<class 'int'>, <class 'numpy.integer'>), 'descr': 'Slit ID number from SPAT measured at half way point.', 'otype': <class 'numpy.ndarray'>}, 'specmax': {'atype': <class 'numpy.floating'>, 'descr': 'Maximum spectral position (pixel units) allowed for each slit/order. Shape is Nslits.', 'otype': <class 'numpy.ndarray'>}, 'specmin': {'atype': <class 'numpy.floating'>, 'descr': 'Minimum spectral position (pixel units) allowed for each slit/order. Shape is Nslits.', 'otype': <class 'numpy.ndarray'>}}
Provides the class data model.
- det_of_slit(spat_id: int, det_img: ndarray, slit_img: ndarray = None)[source]
Identify the ‘best’ detector for this slit/order The metric is the detector where the slit appears the most frequently.
Only sensibly used for mosaic images
- Parameters:
spat_id (int) – spat_id value for the slit of interest
det_img (numpy.ndarray) –
int
image specifying the detector number (1-based) for each pixel in the mosaicslit_img (numpy.ndarray, optional) – image identifying each pixel with its associated slit.
- Returns:
Detector number for the slit (1-based)
- Return type:
- classmethod from_hdu(hdu, chk_version=True, **kwargs)[source]
Instantiate the object from an HDU extension.
This overrides the base-class method, only to add checks (or not) for the bitmask.
- Parameters:
hdu (astropy.io.fits.HDUList, astropy.io.fits.ImageHDU, astropy.io.fits.BinTableHDU) – The HDU(s) with the data to use for instantiation.
chk_version (
bool
, optional) – If True, raise an error if the datamodel version or type check failed. If False, throw a warning only.**kwargs – Passed directly to
_parse()
.
- get_maskdef_extract_fwhm(sobjs, platescale, fwhm_parset, find_fwhm)[source]
This method determines the fwhm to use for the optimal extraction of maskdef_extract (i.e., undetected) objects. If the user provides a fwhm, it would be used. Otherwise fwhm will be computed using the average fwhm of the detected objects.
- Parameters:
sobjs (
SpecObjs
) – List of SpecObj that have been found and traced.platescale (
float
) – Platescale.fwhm_parset (
float
, optional) –Parset
that guides the determination of the fwhm of the maskdef_extract objects. If None (default) the fwhm are computed as the averaged from the detected objects, if it is a number it will be adopted as the fwhm.find_fwhm (
float
) – Initial guess of the objects fwhm in pixels (used in object finding)
- Returns:
FWHM in pixels to be used in the optimal extraction
- Return type:
- get_maskdef_objpos(plate_scale, det_buffer)[source]
Determine the object positions expected by the slitmask design
- get_maskdef_offset(sobjs, platescale, spat_flexure, slitmask_off, bright_maskdefid, snr_thrshd, use_alignbox, dither_off=None)[source]
Determine the Slitmask offset (pixels) from position expected by the slitmask design
- Parameters:
sobjs (
pypeit.specobjs.SpecObjs
) – List of SpecObj that have been found and tracedplatescale (
float
) – Platescalespat_flexure (
float
) – Shifts, in spatial pixels, between this image and SlitTraceslitmask_off (
float
) – User provided slitmask offset in pixelsbright_maskdefid (
str
) – User provided maskdef_id of a bright object to be used to measure offsetsnr_thrshd (
float
) – Objects detected above this S/N ratio threshold will be use to compute the slitmask offsetuse_alignbox (
bool
) – Flag that determines if the alignment boxes are used to measure the offsetdither_off (
float
, optional) – dither offset recorded in the header of the observations
- get_radec_image(wcs, alignSplines, tilts, slit_compute=None, slice_offset=None, initial=True, flexure=None)[source]
Generate an RA and DEC image for every pixel in the frame NOTE: This function is currently only used for SlicerIFU reductions.
- Parameters:
wcs (astropy.wcs.WCS) – The World Coordinate system of a science frame
alignSplines (
pypeit.alignframe.AlignmentSplines
) – An instance of the AlignmentSplines class that allows one to build and transform between detector pixel coordinates and WCS pixel coordinates.tilts (numpy.ndarray) – Spectral tilts.
slit_compute (int, list, numpy.ndarray, optional) – The slit indices to compute the RA and DEC image for. If None, all slits are computed. If an integer, the RA and DEC image for the slit with this index is computed. If a list or numpy.ndarray, the RA and DEC image for the slits with these indices are computed.
slice_offset (float, optional) – Offset to apply to the slice positions. A value of 0.0 means that the slice positions are the centre of the slits. A value of +/-0.5 means that the slice positions are at the edges of the slits. If None, the slice_offset is set to 0.0.
initial (bool) – Select the initial slit edges?
flexure (float, optional) – If provided, offset each slit by this amount.
- Returns:
raimg (numpy.ndarray) – Image with the RA coordinates of each pixel in degrees. Shape is (nspec, nspat).
decimg (numpy.ndarray) – Image with the DEC coordinates of each pixel in degrees. Shape is (nspec, nspat).
minmax (numpy.ndarray) – The minimum and maximum difference (in pixels) between the WCS reference (usually the centre of the slit) and the edges of the slits. Shape is (nslits, 2).
- get_slitlengths(initial=False, median=False)[source]
Get the length of each slit in pixels.
By default, the method will return the tweaked slit lengths if they have been defined. If they haven’t been defined the nominal edges (
left
andright
) are returned. Useinitial=True
to return the nominal edges regardless of the presence of the tweaked edges.- Parameters:
initial (
bool
, optional) – To use the initial edges regardless of the presence of the tweaked edges, set this to True.median (
bool
, optional) – The default is to return the slit length as a function of the spectral coordinate. If median is set to true, the median slit length of each slit is returned.
- Returns:
Slit lengths.
- Return type:
- internals = ['calib_id', 'calib_key', 'calib_dir', 'left_flexure', 'right_flexure']
Attributes kept separate from the datamodel.
- mask_add_missing_obj(sobjs, spat_flexure, fwhm, boxcar_rad)[source]
Generate new SpecObj and add them into the SpecObjs object for any slits missing the targeted source.
- Parameters:
- Returns:
Updated list of SpecObj that have been found and traced
- Return type:
- mask_flats(flatImages)[source]
Mask based on a
FlatImages
object.- Parameters:
flatImages (
FlatImages
)
- mask_wavetilts(waveTilts)[source]
Mask from a
pypeit.wavetilts.WaveTilts
object- Parameters:
waveTilts (
pypeit.wavetilts.WaveTilts
)
- select_edges(initial=False, flexure=None)[source]
Select between the initial or tweaked slit edges and allow for flexure correction.
By default, the method will return the tweaked slits if they have been defined. If they haven’t been defined the nominal edges (
left
andright
) are returned. Useoriginal=True
to return the nominal edges regardless of the presence of the tweaked edges.- Parameters:
- Returns:
Returns the full arrays containing the left and right edge coordinates and the mask, respectively. These are returned as copies.
- Return type:
- slit_img(pad=None, slitidx=None, initial=False, flexure=None, exclude_flag=None, use_spatial=True)[source]
Construct an image identifying each pixel with its associated slit.
The output image has the same shape as the original trace image. Each pixel in the image is set to the index of its associated slit (i.e, the pixel value is \(0..N_{\rm slit}-1\)). Pixels not associated with any slit are given values of -1.
The width of the slit is extended at either edge by a fixed number of pixels using the pad parameter in
par
. This value can be overridden using the method keyword argument.Warning
The function does not check that pixels end up in multiple pixels or that the padding is sensible given the separation between slit edges!
All slits are identified in the image, even if they are masked with
mask
.
- Parameters:
pad (
float
,int
,tuple
, optional) – The number of pixels used to pad (extend) the edge of each slit. This can be a single scale to pad both left and right edges equally or a 2-tuple that provides separate padding for the left (first element) and right (2nd element) edges separately. If not None, this overrides the value inpar
. The value can be negative, which means that the widths are trimmed instead of padded.slitidx (
int
, array_like, optional) – List of indexes (zero-based) to include in the image. If None, all slits not flagged are included.initial (
bool
, optional) – By default, the method will use the tweaked slit edges if they have been defined. If they haven’t been, the initial edges (left_init
andright_init
) are used. To use the nominal edges regardless of the presence of the tweaked edges, set this to True. Seeselect_edges()
.exclude_flag (
str
,list
, optional) – One or more bitmask flag names to ignore when masking Warning – This could conflict with input slitids, i.e. avoid using bothuse_spatial (bool, optional) – If True, use self.spat_id value instead of 0-based indices
flexure (
float
, optional) – If provided, offset each slit by this amount Done in select_edges()
- Returns:
The image with the slit index identified for each pixel.
- Return type:
- property slit_info
THIS NEEDS A DESCRIPTION
- Return type:
- static slit_spat_pos(left, right, nspat)[source]
Return a fidicial, normalized spatial coordinate for each slit.
The fiducial coordinates are given by:
nspec = left.shape[0] (left[nspec//2,:] + right[nspec//2,:])/2/nspat
- Parameters:
left (numpy.ndarray) – Array with left slit edges. Shape is \((N_{\rm spec},N_{\rm slits})\).
right (numpy.ndarray) – Array with right slit edges. Shape is \((N_{\rm spec},N_{\rm slits})\).
nspat (
int
) – Number of pixels in the spatial direction in the image used to trace the slit edges.
- Returns:
Vector with the list of floating point spatial coordinates.
- Return type:
- property slitord_id
Return array of slit_spatId (MultiSlit, SlicerIFU) or ech_order (Echelle) values
- Return type:
- property slitord_txt
Return string indicating if the logs/QA should use “slit” (MultiSlit, SlicerIFU) or “order” (Echelle).
- Returns:
Either ‘slit’ or ‘order’
- Return type:
- spatial_coordinate_image(slitidx=None, full=False, slitid_img=None, pad=None, initial=False, flexure_shift=None)[source]
Generate an image with the normalized spatial coordinate within each slit.
- Parameters:
slitidx (
int
, array_like, optional) – List of indices to include in the image. If None, all slits are included.full (
bool
, optional) – If True, return a full image with the coordinates based on the single provided slit. In this case, slitids must be provided and it can be only one slit! If False, only those coordinates falling in the slit (as determined byslitid_image
or by a call toslit_img()
) will be returned.slitid_img (numpy.ndarray, optional) – Image identifying the slit associated with each pixel. Should have the same shape as expected by the object (
nspec
bynspat
). Pixels not associated with any pixel have a value of -1. This is a convenience parameter that allows previous executions ofslit_img()
to be passed in directly instead of having to repeat the operation. A slit ID image is only required iffull
is False.pad (
float
,int
, optional) – When constructing the slit ID image, use this value to override the value of pad inpar
. Only used ifslitid_img
is not provided directly andfull
is False.initial (
bool
, optional) – By default, the method will use the tweaked slit edges if they have been defined. If they haven’t been, the nominal edges (left
andright
) are used. To use the nominal edges regardless of the presence of the tweaked edges, set this to True. Seeselect_edges()
.
- Returns:
Array specifying the spatial coordinate of pixel in its own slit, scaled to go from 0 to 1. If
full
is True, the image provides the coordinates relative to the left edge for the provided slit over the full image.- Return type:
- spatial_coordinates(initial=False, flexure=None)[source]
Return a fiducial coordinate for each slit.
This is a simple wrapper for
select_edges()
andslit_spat_pos()
.- Parameters:
original (
bool
, optional) – By default, the method will use the tweaked slit edges if they have been defined. If they haven’t been, the nominal edges (left
andright
) are used. To use the nominal edges regardless of the presence of the tweaked edges, set this to True. Seeselect_edges()
.- Returns:
Vector with the list of floating point spatial coordinates.
- Return type:
- version = '1.1.5'
SlitTraceSet data model version.
- pypeit.slittrace.assign_addobjs_alldets(sobjs, calib_slits, spat_flexure, platescale, slitmask_par, find_fwhm)[source]
Loop around all the calibrated detectors to assign RA, DEC and OBJNAME to extracted object and to force extraction of undetected objects.
- Parameters:
sobjs (
SpecObjs
) – List of SpecObj that have been found and traced.calib_slits (numpy.ndarray) – Array of SlitTraceSet with information on the traced slit edges.
spat_flexure (
list
) – List of shifts, in spatial pixels, between this image and SlitTrace.platescale (
list
) – List of platescale for every detector.slitmask_par (
PypeItPar
) – Slitmask PypeIt parameters.find_fwhm (
float
) – Initial guess of the objects fwhm in pixels (used in object finding)
- Returns:
Updated list of spectra that have been found and traced.
- Return type:
- pypeit.slittrace.average_maskdef_offset(calib_slits, platescale, list_detectors)[source]
Loop around all the calibrated detectors to compute the median offset between the expected and measure slitmask position. This info is recorded in the SlitTraceSet datamodel.
- Parameters:
calib_slits (
list
) – List ofSlitTraceSet
objects with information on the traced slit edges.platescale (
float
) – Platescale, must be the same for every detector.list_detectors (numpy.ndarray) – An array that lists the detector numbers of the current spectrograph; see
list_detectors()
. If there are multiple detectors along the dispersion direction, there are ordered along the first axis. For example, all the “bluest” detectors would be inlist_detectors[0]
.
- Returns:
Array of
SlitTraceSet
objects with updated information on the traced slit edges.- Return type:
- pypeit.slittrace.get_maskdef_objpos_offset_alldets(sobjs, calib_slits, spat_flexure, platescale, det_buffer, slitmask_par, dither_off=None)[source]
Loop around all the calibrated detectors to extract information on the object positions expected by the slitmask design and the offsets between the expected and measure slitmask position. This info is recorded in the SlitTraceSet datamodel.
- Parameters:
sobjs (
pypeit.specobjs.SpecObjs
) – List of SpecObj that have been found and tracedcalib_slits (
list
) – List of SlitTraceSet with information on the traced slit edgesspat_flexure (
list
) – List of shifts, in spatial pixels, between this image and SlitTraceplatescale (
list
) – List of platescale for every detectordet_buffer (
int
) – Minimum separation between detector edges and a slit edgeslitmask_par (
pypeit.par.pypeitpar.PypeItPar
) – slitmask PypeIt parametersdither_off (
float
, optional) – dither offset recorded in the header of the observations
- Returns:
List of SlitTraceSet with updated information on the traced slit edges