pypeit.edgetrace module

The primary purpose of this module is to provide the classes/methods used to trace slit edges.

For a command-line script that executes the automatic tracing, use pypeit_trace_edges. As always, for a list of the script options, run:

$ pypeit_trace_edges -h

With a PypeIt Reduction File, a typical execution of the script would be:

$ pypeit_trace_edges -f my_pypeit_file.pypeit

To show the trace results after completing each stage and/or to run in debugging mode, use the –show and/or –debug options:

$ pypeit_trace_edges -f my_pypeit_file.pypeit --debug --show

Programmatically, if you have a PypeIt Reduction File and a path for the reductions (redux_path), an example of how to trace the slits in a single detector is as follows:

# Imports
from pypeit.pypeit import PypeIt
from pypeit import edgetrace
from pypeit.images import buildimage

# Instantiate the PypeIt class to perform the necessary setup
rdx = PypeIt(pypeit_file, redux_path=redux_path)

# Find the trace frames files for a specific calibration group
group = 0
tbl_rows = rdx.fitstbl.find_frames('trace', calib_ID=group, index=True)
files = rdx.fitstbl.frame_paths(tbl_rows)

# Select a detector to trace
det = 1

# Setup the output paths for the trace file; these can be anything but
# the defaults are below
calib_dir = rdx.par['calibrations']['caldir']
setup = rdx.fitstbl['setup'][tbl_rows[0]]
calib_id = rdx.fitstbl['calib'][tbl_rows[0]]

# Skip the bias subtraction, if reasonable; see
# pypeit.biasframe.BiasFrame to construct a bias to subtract from
# the TraceImage
rdx.par['calibrations']['traceframe']['process']['bias'] = 'skip'

# Construct the TraceImage
traceImage = buildimage.buildimage_fromlist(rdx.spectrograph, det,
                                            rdx.par['calibrations']['traceframe'],
                                            files, calib_dir=self.calib_dir,
                                            setup=setup, calib_id=calib_id)

# Then run the edge tracing.  This performs the automatic tracing.
edges = edgetrace.EdgeTraceSet(traceImage, rdx.spectrograph,
                               rdx.par['calibrations']['slitedges'], auto=True)
# You can look at the results using the show method:
edges.show()
# Or in ginga viewer
edges.show(in_ginga=True)
# And you can save the results to a file
edges.to_file()

If you want to instead start without a pypeit file, you could do the following for, e.g., a single unbinned Keck DEIMOS flat-field exposure in a fits file called trace_file:

import os
from pypeit import edgetrace
from pypeit.images import buildimage
from pypeit.spectrographs.util import load_spectrograph

spec = load_spectrograph('keck_deimos')
par = spec.default_pypeit_par()
par['calibrations']['traceframe']['process']['bias'] = 'skip'
# Make any desired changes to the parameters here
det = 3
calib_dir = par['calibrations']['caldir']

# Construct the TraceImage
traceImage = buildimage.buildimage_fromlist(spec, det, par['calibrations']['traceframe'],
                                            [trace_file], calib_dir=self.calib_dir,
                                            setup='A', calib_id=1)

edges = edgetrace.EdgeTraceSet(traceImage, spec, par['calibrations']['slitedges'], auto=True)
edges.to_file()
class pypeit.edgetrace.EdgeTraceBitMask[source]

Bases: BitMask

Mask bits used during slit tracing.

property bad_flags

List the flags that mean the trace is bad.

property exclude_flags

List of flags to exclude when finding bad trace data.

property insert_flags

List of flags used to mark traces inserted for various reasons.

property order_flags

List of flags related to the echelle order.

class pypeit.edgetrace.EdgeTraceSet(traceimg, spectrograph, par, qa_path=None, auto=False, debug=False, show_stages=False)[source]

Bases: CalibFrame

Core class that identifies, traces, and pairs edges in an image to define the slit apertures.

The instantiation of the object can either be used to produce an empty placeholder that you then use multiple times to trace different images, or you can have the tracing begin immediately upon instantiation. For the latter, you must provide (at minimum) the img to trace, and the initialization then run initial_trace() or auto_trace(), depending on the value of auto. To automatically have the instantiation save the results, set save=True on instantiation.

To load an existing calibration file with the result of a trace, use the from_file method:

edges = EdgeTraceSet.from_file(file)

Most commonly, one will use the automatic tracing routine to trace the slit edges; see the description of the steps used during auto-tracing in the docs for auto_trace().

The success of the tracing critically depends on the parameters used. The defaults are tuned for each spectrograph based on testing using data in the PypeIt development suite. See PypeItPar Keywords for the full documentation of the EdgeTracePar parameters.

Finally, note that the design and object data are currently empty, as part of a development path for matching slits traced on the detector to slits expected from provided metadata. Once finished these objects will only contain data for spectrograph output files that provide the relevant metadata.

The datamodel attributes are:

Version: 1.0.1

Attribute

Type

Array Type

Description

PYP_SPEC

str

PypeIt spectrograph name

dispname

str

Spectrograph disperser name.

edge_cen

numpy.ndarray

float, numpy.floating

(Floating-point) Measured spatial coordinate of the edge traces for each spectral pixel. Shape is (Nspec,Ntrace).

edge_err

numpy.ndarray

float, numpy.floating

Error in the measured spatial coordinate edge traces.

edge_fit

numpy.ndarray

float, numpy.floating

The best-fit model result for the trace edge.

edge_msk

numpy.ndarray

numpy.int32

Bitmask for the edge trace positions.

fittype

str

An informational string identifying the type of model used to fit the trace data. Either pca for a PCA decomposition or the polynomial function type and order

left_pca

TracePCA

The PCA decomposition of the left-edge traces. Not defined if PCA performed on all traces, regardless of edge side (i.e., the left_right_pca parameter is False). See TracePCA.

maskdef_id

numpy.ndarray

int, numpy.integer

slitmask ID number for the edge traces. IDs are for, respectively, left and right edges. Only defined if mask-design metadata is available.

nspat

int

Image pixels in the spatial direction.

nspec

int

Image pixels in the spectral direction.

orderid

numpy.ndarray

int, numpy.integer

For echelle spectrographs, this is the order ID number for the edge traces. Negative and positive IDs are for, respectively, left and right edges.

pca

TracePCA

The PCA decomposition of all edge traces. Not defined if PCA separated between left and right traces (i.e., the left_right_pca parameter is True). See TracePCA.

pcatype

str

String identifier for the measurements used to construct the PCA (center or fit)

right_pca

TracePCA

The PCA decomposition of the right-edge traces. Not defined if PCA performed on all traces, regardless of edge side (i.e., the left_right_pca parameter is False). See TracePCA.

sobelsig

numpy.ndarray

float, numpy.floating

Sobel-filtered image used to detect edges

tracebpm

numpy.ndarray

bool, numpy.bool

Bad-pixel mask for trace image

traceid

numpy.ndarray

int, numpy.integer

ID number for the edge traces. Negative and positive IDs are for, respectively, left and right edges.

traceimg

TraceImage

Image used to construct the edge traces; see TraceImage and PypeItImage.

Parameters:
  • traceimg (TraceImage) – Two-dimensional image used to trace slit edges. The object provides the image, the bad-pixel mask, the detector information, and (one of) the original raw file name when matching slits to a design file.

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

  • par (EdgeTracePar) – The parameters used to guide slit tracing. Used to set par.

  • qa_path (str, Path, optional) – Directory for QA output. If None, no QA plots are provided.

  • auto (bool, optional) – Find the edge traces using auto_trace(). If False, the trace data will only be the result of running initial_trace().

  • debug (bool, optional) – Run in debug mode.

  • show_stages (bool, optional) – After ever stage of the auto trace prescription (auto=True), show the traces against the image using show().

traceimg

(TraceImage): See argument list.

spectrograph

(Spectrograph): See argument list.

par

See argument list.

Type:

EdgeTracePar

files

The list of raw files used to construct the trace image (img). Only defined if argument img in initial_trace() or auto_trace() is a TraceImage object.

Type:

list

img

Convenience for now.

Type:

numpy.ndarray

det

See argument list.

Type:

int

sobelsig

Sobel-filtered image used to detect left and right edges of slits.

Type:

numpy.ndarray)

sobelsig_left

Lazy-loaded version of sobelsig that clips the features related to right edges. Only kept for convenience.

Type:

numpy.ndarray

sobelsig_right

Lazy-loaded version of sobelsig that clips the features related to left edges. Only kept for convenience.

Type:

numpy.ndarray

nspec

Number of spectral pixels (rows) in the trace image (axis=0).

Type:

int

nspat

Number of spatial pixels (columns) in the trace image (axis=1).

Type:

int

traceid

The list of unique trace IDs.

Type:

numpy.ndarray

edge_img

An integer array with the spatial pixel nearest to each trace edge. This is identically:

self.edge_img = np.round(self.edge_cen
                         if self.edge_fit is None
                         else self.edge_fit).astype(int)
Type:

numpy.ndarray

edge_cen

A floating-point array with the location of the slit edge for each spectral pixel as measured from the trace image. Shape is \((N_{\rm spec},N_{\rm trace})\).

Type:

numpy.ndarray

edge_err

Error in slit edge locations; measurements without errors have their errors set to -1.

Type:

numpy.ndarray

edge_msk

An integer array with the mask bits assigned to each trace centroid; see EdgeTraceBitMask.

Type:

numpy.ndarray

edge_fit

A model fit to the edge_cen data.

Type:

numpy.ndarray

fittype

An informational string identifier for the type of model used to fit the trace data.

Type:

str

pca

Result of a PCA decomposition of the edge traces, used to predict new traces. This can either be a single pypeit.tracepca.TracePCA object or a list of two pypeit.tracepca.TracePCA objects if the PCA decomposition is peformed for the left (pca[0]) and right (pca[1]) traces separately.

Type:

list, TracePCA

pcatype

An informational string indicating which data were used in the PCA decomposition, ‘center’ for edge_cen or ‘fit’ for edge_fit.

Type:

str

design

Collated slit-mask design data matched to the edge traces.

Type:

astropy.table.Table

objects

Collated object ID and coordinate information matched to the design table.

Type:

astropy.table.Table

qa_path

Directory for QA output. If None, no QA plots are provided.

Type:

Path

log

A list of strings indicating the main methods applied when tracing.

Type:

list

maskdef_id

An integer array with the slitmask IDs assigned to each trace.

Type:

numpy.ndarray

omodel_bspat

A floating-point array with the location of the slit LEFT edge, averaged along the spectral direction, predicted by the optical model (before x-correlation with traced edges)

Type:

numpy.ndarray

omodel_tspat

A floating-point array with the location of the slit RIGHT edge, averaged along the spectral direction, predicted by the optical model (before x-correlation with traced edges)

Type:

numpy.ndarray

coeff_b

A floating-point array with the coefficients (offset, scale) of the x-correlation between LEFT edges predicted by the optical model and the ones traced on the image.

Type:

numpy.ndarray

coeff_t

A floating-point array with the coefficients (offset, scale) of the x-correlation between RIGHT edges predicted by the optical model and the ones traced on the image.

Type:

numpy.ndarray

maskfile

Full path to the file used to extract slit-mask information

Type:

str

_base_header(hdr=None)[source]

Construct the baseline header for all HDU extensions.

This appends the EdgeTracePar and EdgeTraceBitMask data to the headers of all HDU extensions. This is overkill, but avoids overriding pypeit.datamodel.DataContainer.to_hdu() by just including the data in all HDU extensions.

Parameters:

hdr (astropy.io.fits.Header, optional) – Baseline header to add to all returned HDUs. If None, set by pypeit.io.initialize_header().

Returns:

Header object to include in all HDU extensions.

Return type:

astropy.io.fits.Header

_bundle()[source]

Bundle the object for writing using to_hdu().

_fill_design_table(maskdef_id, cc_params_b, cc_params_t, omodel_bspat, omodel_tspat, spat_id)[source]

Fill design based on the results of the design registration.

The design is an astropy.table.Table with 13 columns:
  • ‘TRACEID’: Trace ID Number

  • ‘TRACESROW’: Spectral row for provided left and right edges

  • ‘TRACELPIX’: Spatial pixel coordinate for left edge

  • ‘TRACERPIX’: Spatial pixel coordinate for right edge

  • ‘MASKDEF_ID’: Slit ID Number from slit-mask design

  • ‘SLITLMASKDEF’: Left edge of the slit in pixel from slitmask design before x-correlation

  • ‘SLITRMASKDEF’: Right edge of the slit in pixel from slitmask design before x-correlation

  • ‘SLITRA’: Right ascension of the slit center (deg)

  • ‘SLITDEC’: Declination of the slit center (deg)

  • ‘SLITLEN’: Slit length (arcsec)

  • ‘SLITWID’: Slit width (arcsec)

  • ‘SLITPA’: Slit position angle on sky (deg from N through E)

  • ‘ALIGN’: Slit used for alignment (1-yes; 0-no), not target observations.

And three .meta info:
  • ‘MASKFILE’: name of file with the slitmask info

  • ‘MASKOFFL’, ‘MASKOFFR’: The coefficient ‘offset’ of the x-correlation between edges predicted by

    the slitmask design and the one traced on the image. One value per each edge side.

  • ‘MASKSCLL’, ‘MASKSCLR’: The coefficient ‘scale’ of the x-correlation between edges predicted by

    the slitmask design and the one traced on the image. One value per each edge side

  • ‘MASKRMSL’, ‘MASKRMSR’: The RMS of the x-correlation between edges predicted by the slitmask design

    and the one traced on the image. One value per each edge side

Parameters:
  • maskdef_id (numpy.ndarray) – Slit ID number from slit-mask design matched to traced slits.

  • cc_params_b (tuple) – Three parameters of the cross-correlation (2 coefficients and RMS) between slit-mask design and traced edges for the left and right edges.

  • cc_params_t (tuple) – Three parameters of the cross-correlation (2 coefficients and RMS) between slit-mask design and traced edges for the left and right edges.

  • omodel_bspat (numpy.ndarray) – Left and right spatial position of the slit edges from optical model

  • omodel_tspat (numpy.ndarray) – Left and right spatial position of the slit edges from optical model

  • spat_id (numpy.ndarray) – ID assigned by PypeIt to each slit. same as in SlitTraceSet.

_fill_objects_table(maskdef_id)[source]

Fill objects based on the result of the design registration.

The objects is an astropy.table.Table with 5 columns:
  • ‘OBJID’: Object ID Number

  • ‘OBJRA’: Right ascension of the object (deg)

  • ‘OBJDEC’: Declination of the object (deg)

  • ‘OBJNAME’: Object name assigned by the observer

  • ‘OBJMAG’: Object magnitude provided by the observer

  • ‘OBJMAG_BAND’: Band of the magnitude provided by the observer

  • ‘MASKDEF_ID’: Slit ID Number from slit-mask design

  • ‘OBJ_TOPDIST’: Projected distance (in arcsec) of the object from the left

    edge of the slit (in PypeIt orientation)

  • ‘OBJ_BOTDIST’: Projected distance (in arcsec) of the object from the right

    edge of the slit (in PypeIt orientation)

  • ‘TRACEID’: Row index that matches ‘TRACEID’ in the design table

Parameters:

maskdef_id (numpy.ndarray) – Slit ID number from slit-mask design matched to traced slits.

_flag_edges(trace_cen, indx, flg)[source]

Convenience function for slit checking, which performs the same operations to the set of edges flagged for the many reasons iterated through in check_synced().

This modifies edge_msk directly.

Parameters:
  • trace_cen (numpy.ndarray) – The 2D array with the edge traces.

  • indx (numpy.ndarray) – The boolean array selecting the edges to be flagged.

  • flg (str) – The bit flag to be assigned.

_get_insert_locations()[source]

Find where edges need to be inserted.

This only determines where the left-right ordering of the traces implies that a trace needs to be inserted. Where the trace is inserted and with what shape is determined by _get_reference_locations() and sync(), respectively.

Returns:

Three numpy.ndarray objects are returned:

  • An integer vector identifying the type of side for the fully synchronized edge set. Elements of the vector should alternate left (-1) and right (1).

  • A boolean array selecting the edges in the returned list of sides that should be added to the existing trace set.

  • An array with the indices in the existing trace arrays where the new traces should be inserted.

Return type:

tuple

_get_reference_locations(trace_cen, add_edge)[source]

Determine the reference locations for traces to add during the left-right synchronization.

The positions of new traces are determined either by the median slit length of the existing left-right pairs, or based on the slit length of the nearest left-right pair. This function only determines the positions for the new traces at the reference spetral location (row). The shape is determined by sync().

Used parameters from par (EdgeTracePar) are sync_center, sync_to_edge, and gap_offset.

Parameters:
  • trace_cen (numpy.ndarray) – Trace data to use for determining new edge locations.

  • add_edge (numpy.ndarray) – Boolean array indicating that a trace in the new array is an added trace. The number of False entries in add_edge should match the length of the 2nd axis of trace_cen.

Returns:

Reference positions for all edge traces, both for the existing and new traces.

Return type:

numpy.ndarray

_masked_single_slit(trace_cen)[source]

Handle masking a single slit as too short.

This is largely a kludge that was necessary to address the features seen in the Keck_LRIS_blue long-slit data in the dev suite.

Todo

We need to understand how often a single slit is masked and if this is always the right way to deal with it.

Parameters:

trace_cen (numpy.ndarray) – Trace data to use for determining new edge locations.

_parse_exclude_regions()[source]

Parse the exclude_regions parset.

Returns:

Returns two arrays with the starting and ending pixels of the regions to exclude in this detector

Return type:

tuple

_reinit_trace_data()[source]

Convenience method to set all attributes related to trace data to None.

_reset_pca(rebuild)[source]

” Reset the PCA decomposition.

The PCA is reset by either rebuilding it with the previous set of parameters (rebuild is True) or removing it (setting the relevant attributes to None when rebuild is False).

_side_dependent_sobel(side)[source]

Return the Sobel-filtered image relevant to tracing the given side.

This is primarily a wrapper for prepare_sobel_for_trace(); the boxcar smoothing is always 5 pixels. Barring a instantiation of the object, this calculation is only done once per side by “lazy loading” sobelsig_left or sobelsig_right.

Parameters:

side (str) – The side to return; must be 'left' or 'right'.

Returns:

The manipulated Sobel image relevant to tracing the specified edge side.

Return type:

numpy.ndarray

Raises:

PypeItError – Raised if side is not 'left' or 'right'.

add_user_traces(user_traces, method='straight')[source]

Add traces for user-defined slits.

Parameters:
  • user_slits (list) – A list of lists with the coordinates for the new traces. For each new slit, the list provides the spectral coordinate at which the slit edges are defined and the left and right spatial pixels that the traces should pass through. I.e., [664, 323, 348] mean construct a left edge that passes through pixel (664,323) (ordered spectral then spatial) and a right edge that passes through pixel (664,348).

  • method (str, optional) –

    The method used to construct the traces. Options are:

    • 'straight': Simply insert traces that have a constant spatial position as a function of spectral pixel.

    • 'nearest': Constrain the trace to follow the same form as an existing trace that is nearest the provided new trace coordinates.

    • 'pca': Use the PCA decomposition of the traces to predict the added trace. If the PCA does not currently exist, the function will try to (re)build it.

auto_trace(bpm=None, debug=False, show_stages=False)[source]

Execute a fixed series of methods to automatically identify and trace slit edges.

The current algorithm is:

  • Detect and follow slit edges using initial_trace().

  • Refine the measured centroids of the edge locations using centroid_refine().

  • Fit the measured centroids with a polynomial using fit_refine(), which is basically a wrapper for pypeit.core.trace.fit_trace(). If a PCA decomposition of the traces is not possible because there are too few traces, the next two steps are skipped (skipping down to edge synchronization), and the final measured slit edge traces are the result of this step.

  • Construct a PCA decomposition of the fitted forms for the traces using pca_refine(). Traces are either decomposed as a group or split into separate left and right decompositions.

  • Use the PCA decomposition to rectify the trace image such that the slit edges should be aligned, collapse the image to get a high-signal-to-noise detection of each slit edge, and use these locations to predict, remeasure, and refit the slit edges; see peak_refine(). The final measured slit edge traces are based on the result of peak_refine().

  • Synchronize the left and right traces into pairs that define slit apertures using sync().

  • Use maskdesign_matching() to match the slit edge traces found with the ones predicted by the slit-mask design.

  • Use add_user_traces() and rm_user_traces() to add and remove traces as defined by the user-provided lists in the par.

Used parameters from par (EdgeTracePar) are use_maskdesign.

Parameters:
  • bpm (numpy.ndarray, optional) – Bad-pixel mask for the trace image. Must have the same shape as img. If None, all pixels are assumed to be valid.

  • debug (bool, optional) – Run in debug mode.

  • show_stages (bool, optional) –

    After ever stage of the auto trace, execute show() to show the results. These calls to show are always:

    self.show()
    

bitmask = <pypeit.edgetrace.EdgeTraceBitMask object>

Bit interpreter for the edge tracing mask.

bound_detector()[source]

Insert traces at both detector boundaries.

Accounting for the requested detector buffer, a left and right trace are placed at the detector edges. Traces are masked as user-inserted.

Only used parameter from par (EdgeTracePar) is det_buffer.

build_pca(use_center=False, debug=False)[source]

Build a PCA model of the current edge data.

Primarily a wrapper that instantiates pca, or left_pca and right_pca if left and right traces are analyzed separately. All of these objects will be a pypeit.tracepca.TracePCA, if instantiated. After executing this method, traces can be predicted by pypeit.tracepca.TracePCA.predict() for the relevant PCA decomposition; see predict_traces().

If no parametrized function has been fit to the trace data or if specifically requested (see use_center), the PCA is based on the measured trace centroids (edge_cen); othwerwise, the PCA uses the parametrized trace fits (edge_fit).

The reference spectral row used for the decomposition (see pypeit.tracepca.TracePCA) is set by pypeit.core.trace.most_common_trace_row() using the existing mask. If treating left and right traces separately, the reference spectral row is the same for both PCA decompositions.

Used parameters from par (EdgeTracePar) are fit_min_spec_length, left_right_pca, pca_n, pca_var_percent, pca_function, pca_order, pca_sigrej, pca_maxrej, and pca_maxiter.

Parameters:
  • use_center (bool, optional) – Use the center measurements for the PCA decomposition instead of the functional fit to those data. This is only relevant if both are available. If no fits have been performed, the function will automatically use the center measurements.

  • debug (bool, optional) – Run in debug mode.

calib_file_format = 'fits.gz'

Calibration frame file format.

calib_type = 'Edges'

Root name for the calibration frame file.

can_pca()[source]

Determine if traces are suitable for PCA decomposition.

The criterion is that a minimum number of traces (pca_min_edges) must cover more than the fraction of the full spectral range specified by fit_min_spec_length in par. Traces that are inserted are ignored.

If the PCA decomposition will be performed on the left and right traces separately, the function will return False if there are fewer than the minimum left or right edge traces.

Used parameters from par (EdgeTracePar) are fit_min_spec_length, left_right_pca, and pca_min_edges.

Warning

This function calls check_traces() using fit_min_spec_length to flag short traces, meaning that edge_msk will can be altered by this call.

Returns:

Flag that traces meet criterion for PCA decomposition.

Return type:

bool

centroid_refine(follow=True, start_indx=None, continuous=False, use_fit=False)[source]

Refine the edge positions using a moment analysis and assess the results.

The starting point for all trace positions is edge_cen, unless use_fit is True.

The method runs in two primary modes, depending on the value of follow:

When follow=False, the function simply executes masked_centroid() to recenter all the locations for the left and right edges at each spectral position (row) independently. The maximum shift between any input and output centroid is max_shift_abs from par.

When follow=True, the method uses follow_centroid() to recenter each edge starting from the start_indx spectral index position (row) and then follows the centroid to higher and lower spectral rows. In this case, the center of the aperture used for each spectral row is the centroid of the trace measured for the preceding spectral row. The maximum shift between the input and output center for the first position analyzed is set by max_shift_abs and the maximum shift for each subsequent spectral row is set by max_shift_adj; both parameters are provided in par. In this approach, it’s typically best to let the method determine the starting spectral row instead of providing it directly. If left to its own devices, it will iterate through all the traces collecting those that all cross a specific spectral row into groups that it can follow simultaneously. If a starting specral row is provided directly, all traces must cross that row.

Regardless of the value of follow, masked_centroid() is run with uniform weighting and an aperture width set to be twice fwhm_uniform from par.

During the refinement, check_traces() after each new trace is refined, indentifying repeat and errorneous traces. (The minimum length of the trace passed to check_traces() is set by det_min_spec_length in par.) If clip in par is True, all bad traces are removed (see remove_traces()).

Other used parameters from par (EdgeTracePar) are max_shift_abs, max_shift_adj, and max_spat_error.

Warning

  • This function modifies the internal trace arrays in place.

  • Because this changes edge_cen and edge_err, any model fitting of these data are erased by this function! I.e., edge_fit and fittype are set to None.

  • This always removes the PCA if it exists.

Parameters:
  • follow (bool, optional) – Perform the centroiding first at a single row (see start_indx) and then move to higher and lower spectral rows in series. See follow_centroid().

  • start_indx (int, optional) – The index of the starting spectral row when following the trace between adjacent rows; see follow. If None, the starting spectral row is set by finding the spectral row that crosses the most unmasked trace positions; see most_common_trace_row(). Value is ignored if follow=False.

  • continuous (bool, optional) – Keep only the continuous part of the traces from the starting spectral row. See follow_centroid().

  • use_fit (bool, optional) – Use the fitted traces as the starting point for the refinement. Otherwise, uses :edge_cen. If True and edge_fit is None, the method will raise an exception.

check_synced(rebuild_pca=False)[source]

Quality check and masking of the synchronized edges.

Before executing this method, the slit edges must be synchronized (see sync()) and ordered spatially in left-right pairs (see spatial_sort()); only the former is checked explicitly. Any traces fully masked as bad (see clean_traces()) are removed, along with its synchronized partner.

Used parameters from par (EdgeTracePar) are minimum_slit_gap, minimum_slit_length, minimum_slit_length_sci, and length_range.

Checks are:

  • Any trace falling off the edge of the detector is masked (see EdgeTraceBitMask). This is the only check performed by default (i.e., when no keyword arguments are provided).

  • Traces that form slit gaps (the median difference between the right and left traces of adjacent slits) that are below an absolute tolerance are removed and the two relevant slits are merged. This is done before the checks of the slit length below such that the merged slit is assessed in any expected slit length constraints.

  • Traces that form a slit with a length (the median difference between the left and right edges) below an absolute tolerance (i.e., right-left < atol) are masked as SHORTSLIT (see EdgeTraceBitMask). The absolute tolerance is set using the platescale provided by the spectrograph class, the spatial binning (from binning), and the minimum slit length in arcsec (minimum_slit_length in par).

  • Traces that form a slit with a length (the median difference between the left and right edges) below an absolute tolerance (i.e., right-left < atol) for a science slit are masked as either BOXSLIT or SHORTSLIT (see EdgeTraceBitMask), depending on the value of minimum_slit_length: if minimum_slit_length is None, all are flagged as BOXSLIT; otherwise, BOXSLITs are those that are larger than minimum_slit_length but smaller than minimum_slit_length_sci. The absolute tolerance is set using the platescale provided by the spectrograph class, the spatial binning (from binning), and the minimum slit length in arcsec for the science slits (minimum_slit_length_sci in par).

  • Traces that form a slit with a length that is abnormal relative to the median width are masked; i.e., abs(log((right[0]-left[0])/median(right-left))) > log(1+rtol), where rtol is identically length_range in par.

Parameters:

rebuild_pca (bool, optional) – If the pca exists and the masking or number of traces changes because of the performed checks, redo the PCA decomposition using the new traces and trace masks and the previous parameter set.

Returns:

Flag that checking and cleaning traces maintained left-right syncronization. If False, traces will need to be re-syncronized.

Return type:

bool

check_traces(cen=None, msk=None, subset=None, min_spatial=None, max_spatial=None, minimum_spec_length=None)[source]

Validate new or existing trace data.

This function only changes the mask status of the traces and returns flags for traces that are good or bad. Traces are flagged as bad if they meet any of the following conditions:

  • If and only if subset is provided, the selected subset of traces are checked against all other traces to find repeats. Repeats are identified as those tracing the same side of a slit and having a minimum spatial separation of less than the provided match_tol parameter in par. Beware that this should be done by providing cen directly.

  • Traces that do not extend for at least some fraction of the detector (see minimum_spec_length).

  • Traces with spatial positions at the central spectral channel that are above or below the provided threshold (see min_spatial and max_spatial).

  • Traces that fully land off the edge of the detector; see trace_pixels_off_detector().

The only used parameter from par (EdgeTracePar) is match_tol.

Warning

msk is edited in-place

Parameters:
  • cen (numpy.ndarray, optional) – The adjusted center of the refined traces. Shape is \((N_{\rm spec}, N_{\rm refine},)\). If None, use edge_cen.

  • msk (numpy.ndarray, optional) – The mask bits for the adjusted center of the refined traces. Shape is \((N_{\rm spec}, N_{\rm refine},)\). If None, use edge_msk. This is edited in-place!

  • subset (numpy.ndarray, optional) – Boolean array selecting the traces to compare. Shape is \((N_{\rm trace},)\), with \(N_{\rm refine}\) True values. It is expected that all the traces selected by a subset must be from the same slit side (left or right). If None, no repeat traces can be identified.

  • min_spatial (int, optional) – Clip traces that hit this minimum spatial index (column) at the center spectral row (self.nspec//2). If None, no traces clipped.

  • max_spatial (int, optional) – Clip traces that hit this maximum spatial index (column) at the center spectral row (self.nspec//2). If None, no traces clipped.

  • minimum_spec_length (float, optional) – The minimum number of spectral rows in an edge trace.

Returns:

Returns two boolean arrays selecting the good and bad traces. Both are returned because neither will be true for traces not selected by subset if it is provided. The shape of both arrays is \((N_{\rm trace},)\).

Return type:

tuple

clean_traces(force_flag=None, rebuild_pca=True, sync_mode='ignore', assume_synced=False)[source]

Remove any traces that are fully masked as bad.

Traces selected for removal must be fully masked; see fully_masked_traces().

By default, flags used to select bad trace measurements are those provided by EdgeTraceBitMask.bad_flags(), and those flags excluded from designating a trace as bad are provided by EdgeTraceBitMask.exclude_flags(). To force removal of traces with certain flags, regardless of these two groups, use force_flag. See EdgeTraceBitMask for list of flags.

Parameters:
  • force_flag (str, list, optional) – Force inclusion of these flags in assessing whether or not a trace is fully masked.

  • rebuild_pca (bool, optional) – Rebuild the PCA decomposition of the traces based on the loaded trace data. Passed directly to remove_traces(), which is only called if traces are in fact removed.

  • sync_mode (str, optional) – If the traces are left-right synchronized (see sync() and is_synced()), use this method to deal with edges paired with those to be removed. See synced_selection().

  • assume_synced (bool, optional) – Assume the set of traces is synced. See synced_selection().

current_trace_locations()[source]

Return an image with the trace IDs at the locations of each edge in the original image.

datamodel = {'PYP_SPEC': {'descr': 'PypeIt spectrograph name', 'otype': <class 'str'>}, 'dispname': {'descr': 'Spectrograph disperser name.', 'otype': <class 'str'>}, 'edge_cen': {'atype': (<class 'float'>, <class 'numpy.floating'>), 'descr': '(Floating-point) Measured spatial coordinate of the edge traces for each spectral pixel.  Shape is (Nspec,Ntrace).', 'otype': <class 'numpy.ndarray'>}, 'edge_err': {'atype': (<class 'float'>, <class 'numpy.floating'>), 'descr': 'Error in the measured spatial coordinate edge traces.', 'otype': <class 'numpy.ndarray'>}, 'edge_fit': {'atype': (<class 'float'>, <class 'numpy.floating'>), 'descr': 'The best-fit model result for the trace edge.', 'otype': <class 'numpy.ndarray'>}, 'edge_msk': {'atype': <class 'numpy.int32'>, 'descr': 'Bitmask for the edge trace positions.', 'otype': <class 'numpy.ndarray'>}, 'fittype': {'descr': 'An informational string identifying the type of model used to fit the trace data.  Either ``pca`` for a PCA decomposition or the polynomial function type and order', 'otype': <class 'str'>}, 'left_pca': {'descr': 'The PCA decomposition of the left-edge traces.  Not defined if PCA performed on all traces, regardless of edge side (i.e., the left_right_pca parameter is False).  See :class:`~pypeit.tracepca.TracePCA`.', 'otype': <class 'pypeit.tracepca.TracePCA'>}, 'maskdef_id': {'atype': (<class 'int'>, <class 'numpy.integer'>), 'descr': 'slitmask ID number for the edge traces. IDs are for, respectively, left and right edges.  Only defined if mask-design metadata is available.', 'otype': <class 'numpy.ndarray'>}, 'nspat': {'descr': 'Image pixels in the spatial direction.', 'otype': <class 'int'>}, 'nspec': {'descr': 'Image pixels in the spectral direction.', 'otype': <class 'int'>}, 'orderid': {'atype': (<class 'int'>, <class 'numpy.integer'>), 'descr': 'For echelle spectrographs, this is the order ID number for the edge traces.  Negative and positive IDs are for, respectively, left and right edges.', 'otype': <class 'numpy.ndarray'>}, 'pca': {'descr': 'The PCA decomposition of all edge traces.  Not defined if PCA separated between left and right traces (i.e., the left_right_pca parameter is True).  See :class:`~pypeit.tracepca.TracePCA`.', 'otype': <class 'pypeit.tracepca.TracePCA'>}, 'pcatype': {'descr': 'String identifier for the measurements used to construct the PCA (center or fit)', 'otype': <class 'str'>}, 'right_pca': {'descr': 'The PCA decomposition of the right-edge traces.  Not defined if PCA performed on all traces, regardless of edge side (i.e., the left_right_pca parameter is False).  See :class:`~pypeit.tracepca.TracePCA`.', 'otype': <class 'pypeit.tracepca.TracePCA'>}, 'sobelsig': {'atype': (<class 'float'>, <class 'numpy.floating'>), 'descr': 'Sobel-filtered image used to detect edges', 'otype': <class 'numpy.ndarray'>}, 'tracebpm': {'atype': (<class 'bool'>, <class 'numpy.bool_'>), 'descr': 'Bad-pixel mask for trace image', 'otype': <class 'numpy.ndarray'>}, 'traceid': {'atype': (<class 'int'>, <class 'numpy.integer'>), 'descr': 'ID number for the edge traces.  Negative and positive IDs are for, respectively, left and right edges.', 'otype': <class 'numpy.ndarray'>}, 'traceimg': {'descr': 'Image used to construct the edge traces; see :class:`~pypeit.images.buildimage.TraceImage` and :class:`~pypeit.images.pypeitimage.PypeItImage`.', 'otype': <class 'pypeit.images.buildimage.TraceImage'>}}

DataContainer datamodel.

static empty_design_table(rows=None)[source]

Construct an empty design table.

Parameters:

rows (int, optional) – Number of table rows for each table column, expected to be the number of matched slits. If None, the table has empty columns.

Returns:

Instance of the empty design table.

Return type:

astropy.table.Table

static empty_objects_table(rows=None)[source]

Construct an empty objects table.

Parameters:

rows (int, optional) – Number of table rows for each table column, expected to be the number of objects. If None, the table has empty columns.

Returns:

Instance of the empty object table.

Return type:

astropy.table.Table

fit_refine(weighting='uniform', debug=False, idx=None)[source]

Iteratively re-measure and fit a functional form to the edge locations.

Before fitting the traces, check_traces() is used to flag traces that do not meet a minimum length required for fitting (set by fit_min_spec_length in par).

After this, the function is primarily a wrapper for fit_trace(), run once per edge side (left and right). Both the measured centers and the fitted model are modified by fit_trace(), such that edge_cen, edge_err, edge_msk, edge_fit, edge_fit_type, and edge_img are all modified by this method. Note that only traces that are not fully flagged are fit.

Used parameters from par (EdgeTracePar) are max_shift_abs, max_spat_error, fit_function, fit_order, fwhm_uniform, fwhm_gaussian, fit_maxdev, fit_maxiter, fit_niter, and fit_min_spec_length.

Parameters:
  • weighting (str, optional) – The weighting to apply to the position within each integration window (see fit_trace()).

  • debug (bool, optional) – Run fit_trace() in debug mode.

  • idx (numpy.ndarray, optional) – Array of strings with the IDs for each object. Used only if debug is true for the plotting (see fit_trace()).

property flagged_bits

List the flags that have been tripped.

Returns:

List of the unique flags tripped by any trace datum.

Return type:

list

classmethod from_hdu(hdu, hdu_prefix=None, chk_version=True)[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 muliple DataContainer objects contained by EdgeTraceSet.

Parameters:
fully_masked_traces(flag=None, exclude=None)[source]

Find fully masked edge traces.

Traces are identified as masked by one or more of the flags in flag and explicitly not flagged by any flag in exclude over the fully spectral range of the detector.

Parameters:
  • flag (str, list, optional) – The bit mask flags to select. If None, any flags are used. See pypeit.bitmask.BitMask.flagged().

  • exclude (str, list, optional) – A set of flags to explicitly exclude from consideration as a masked trace. I.e., if any spectral pixel in the trace is flagged with one of these flags, it will not be considered a fully masked trace. This is typically used to exclude inserted traces from being considered as a bad trace.

Returns:

Boolean array selecting traces that are flagged at all spectral pixels.

Return type:

numpy.ndarray

get_slits()[source]

Use the data to instatiate the relevant SlitTraceSet object.

The traced edges must have first been organized into slits; see sync().

For fixed-format echelle spectrographs, the method automatically calls match_order() and will only include the “slits” that have been correctly matched to known orders.

The SlitTraceSet object will use calibration naming keys used by this parent EdgeTraceSet object.

The SlitTraceSet object will have an astropy.table.Table resulted from merging design and objects together.

Returns:

Object holding the slit traces.

Return type:

SlitTraceSet

good_traces(include_box=False, good_orders=False)[source]

Select traces that are not fully masked.

This is just a wrapper for fully_masked_traces().

Parameters:
  • include_box (bool, optional) – Include any identified box slits in the list of good traces.

  • good_orders (bool, optional) – Only include those traces that are well-matched to echelle orders. This is effectively ignored by multi-slit spectrographs because none of the order flags should have been tripped for them.

Returns:

Boolean array with shape \((N_{\rm trace},)\) flagging good traces.

Return type:

numpy.ndarray

initial_trace(bpm=None)[source]

Perform the initial trace of the image.

This effectively reinstantiates the object and must be the first method called for tracing an image. The algorithm:

  • Lightly boxcar smooths the trace image spectrally.

  • Replaces pixel columns and rows that are substantially masked, if a bad-pixel mask is provided (see bpm).

  • Applies a Sobel filter to the trace image along the spatial axis (columns) to detect slit edges using steep positive gradients (left edges) and steep negative gradients (right edges). See pypeit.core.trace.detect_slit_edges().

  • Follows the detected left and right edges along spectrally adjacent pixels to identify coherent traces. See pypeit.core.trace.identify_traces().

  • Initializes the attributes that provide the trace position for each spectral pixel based on these results.

Used parameters from par (EdgeTracePar) are filt_iter, sobel_mode, edge_thresh, and follow_span.

Parameters:

bpm (numpy.ndarray, optional) – Bad-pixel mask for the trace image. Must have the same shape as img. If None, all pixels are assumed to be valid.

insert_traces(side, trace_cen, loc=None, mode='user', resort=True, nudge=True)[source]

Insert/append a set of edge traces.

New traces to add are first nudged away from the detector edge (see nudge_traces()) according to parameters max_nudge and det_buffer from par (EdgeTracePar). They are then inserted or appended to the existing traces and masked according to the provided mode. The traces are added to both the measured centroid list and the fitted model data. Then the full list of traces can be resorted spatially, according to the provided resort.

Typically, the inserted traces will be masked, which means that any existing PCA decomposition will be unchanged. However, if mode is None, these inserted traces would be used in the construction of the PCA.

Warning

If the traces can be nudged away from the detector edge, the offset can, e.g., place an inserted left edge to the right of its associated right edge. This possibility is currently not handled by this function.

Parameters:
  • side (int, numpy.ndarray) – Side for each trace to be added: -1 for left, 1 for right. Shape is \((N_{\rm new},)\).

  • trace_cen (numpy.ndarray) – Array with one or more vectors of trace locations. Can be 1D or 2D with shape \((N_{\rm spec},)\) or \((N_{\rm spec}, N_{\rm new})\), respectively.

  • loc (int, numpy.ndarray, optional) – Indices in the current trace arrays at which to insert the new traces; see numpy.insert. If None, traces are appended.

  • mode (str, optional) –

    Mode used for generating the traces to insert used to flag the traces. Options are:

    • None: Traces are simply inserted without flagging.

    • 'user': Traces are the result of a user request.

    • 'sync': Traces were generated by synchronizing left and right traces.

    • 'mask': Traces were generated based on the expected slit positions from mask design data.

    • 'order': Traces are the expected location of an echelle order.

  • resort (bool, optional) – Resort the traces in the spatial dimension; see spatial_sort().

  • nudge (bool, optional) – Allow the traces to be nudged away from the detector edge according to par and nudge_traces().

internals = ['calib_id', 'calib_key', 'calib_dir', 'spectrograph', 'par', 'qa_path', 'edge_img', 'sobelsig_left', 'sobelsig_right', 'design', 'objects', 'log', 'omodel_bspat', 'omodel_tspat', 'cc_params_b', 'cc_params_t', 'maskfile', 'slitmask', 'success']

Attributes kept separate from the datamodel.

property is_empty

Flag that object has no trace data.

property is_left

Boolean array selecting the left traces.

property is_right

Boolean array selecting the right traces.

property is_synced

Boolean that left-right traces are synchronized into slits.

For this to be true:

  • all traces must be valid; i.e., np.all(self.good_traces(include_box=True)) must be True.

  • The sorted traces must be ordered left-right-left-right-…

maskdesign_matching(debug=False)[source]

# TODO - break up this method (too long) Match slit info from the mask design data to the traced slits.

Use of this method requires:
  • a PCA decomposition is available,

  • spectrograph has a viable get_maskdef_slitedges method to read edge trace locations predicted by the slitmask design. This data can be pulled from one of the files used to construct the trace image.

The method uses a collection of scripts in pypeit.core.slitdesign_matching which are taken from DEEP2 IDL-based pipeline for DEIMOS data.

NOTE this method was first written to deal with edge traces predicted by DEIMOS optical model, but subsequently adapted to be used with other spectrographs that don’t use optical models.

Parameters:

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

match_order(reference_row=None)[source]

Match synchronized slits to the expected echelle orders.

This function will fault if called for non-Echelle spectrographs!

For Echelle spectrographs, this finds the best matching order for each left-right trace pair; the traces must have already been synchronized into left-right pairs. Currently, this is a very simple, non-optimized match:

  • The closest order from self.spectrograph.order_spat_pos is matched to each slit.

  • Any slit that cannot be matched to an order – either because there are more “slits” than orders or because the separation is larger than the provided tolerance – is flagged as ORDERMISMATCH.

  • A warning is issued if the number of valid matches is not identical to the number of expected orders (self.spectrograph.norders). The warning includes the list of orders that were not identified.

The match tolerance is set by the parameter order_match. An offset can be applied to improve the match via the parameter order_offset; i.e., this should minimize the difference between the expected order positions and self.slit_spatial_center() + self.par['order_offset']. Both order_match and order_offset are given in fractions of the detector size along the spatial axis.

The result of this method is to instantiate orderid.

Parameters:

reference_row (int, optional) – The spectral pixel (row) used to generate spatial positions of the orders to match against the expected positions. If None, use the PCA reference row if a PCA exists or the central row (i.e., self.nspec//2), otherwise.

Returns:

The median offset in the relative of the detector size between the archived order positions and those measured via the edge tracing.

Return type:

float

Raises:

PypeItError – Raised if the number of orders, the order number, or their spatial locations are not defined for an Echelle spectrograph. Also raised if the edge traces are not synced.

merge_traces(merge_frac=0.5, refit=True, debug=False)[source]

Merge traces based on their spatial separation.

Traces are merged if:

  • they are from the same slit side (left or right)

  • the fitted trace locations of two different traces are separated by less than match_tol (from par) for more than merge_frac of the spectral range of the detector.

Since the merging is based on the fitted trace location, the traces must have been fit beforehand; see fit_refine().

When there are traces found that should be merged, the unmasked centroid measurements from the shorter trace(s) are added to the longest trace (meaning that the shorter trace is effectively removed). The traces are then resorted and given new trace IDs.

Warning

  • This automatically removes any existing PCA decomposition

  • If functional fits exist and refit is False, all fits are removed.

  • If traces are merged, a warning is issued that the traces may no longer be left-right synchronized.

Parameters:
  • merge_frac (float, optional) – The fraction of the spectral length of the detector that should be less than match_tol in par used to find traces to merge.

  • refit (bool, optional) – If fitted models exist and traces are merged, run fit_refine() with its default arguments after the traces have been merged to refit the trace data.

  • debug (bool, optional) – Only passed to fit_refine() if refit is True.

property nslits
property ntrace

The number of edges (left and right) traced.

nudge_traces(trace_cen)[source]

Nudge traces away from the detector edge.

Traces are shifted spatially, up to a maximum value (max_nudge), to be no closer than a minimum number (det_buffer) pixels from the detector edges. Both parameters are pulled from par (EdgeTracePar). No limit is imposed on the size of the shift if max_nudge is None.

Warning

A failure mode that is not dealt with is when multiple traces fall off the detector and are nudged to nearly the same location.

Parameters:

trace_cen (numpy.ndarray) – Array with trace locations to adjust. Must be 2D with shape \((N_{\rm spec}, N_{\rm trace})\).

Returns:

The nudged traces.

Return type:

numpy.ndarray

order_refine(debug=False)[source]

For echelle spectrographs, attempt to add any orders that are not present in the current set of edges.

order_refine_fixed_format(reference_row, debug=False)[source]

Refine the order locations for fixed-format Echelles.

order_refine_free_format(reference_row, combined_order_tol=1.8, debug=False)[source]

Refine the order locations for “free-format” Echelles.

Traces must be synced before reaching here.

order_refine_free_format_qa(cen, combined_orders, width, gap, width_fit, gap_fit, order_cen, order_missing, ofile=None)[source]

QA plot for order placements

output_float_dtype

Regardless of datamodel, output floating-point data have this fixed bit size.

alias of float32

pca_refine(use_center=False, debug=False, force=False)[source]

Use a PCA decomposition to refine the traces.

If no parametrized function has been fit to the trace data or if specifically requested (see use_center), the PCA is based on the measured trace centroids (edge_cen); othwerwise, the PCA uses the parametrized trace fits (edge_fit).

If needed or forced to, this first executes build_pca() and then uses predict_traces() to use the PCA to reset the trace data.

Only used parameter from par (EdgeTracePar) is left_right_pca.

Parameters:
  • use_center (bool, optional) – Use the center measurements for the PCA decomposition instead of the functional fit to those data. This is only relevant if both are available. If no fits have been performed, the function will automatically use the center measurements.

  • debug (bool, optional) – Run in debug mode.

  • force (bool, optional) – Force the recalculation of the PCA even if it has already been done.

peak_refine(rebuild_pca=False, debug=False)[source]

Refine the trace by isolating peaks and troughs in the Sobel-filtered image.

This function requires that the PCA model exists; see build_pca() or pca_refine(). It is also primarily a wrapper for peak_trace(). See the documentation of that function for the explanation of the algorithm.

If the left and right traces have separate PCA decompositions, this function makes one call to peak_trace() for each side. Otherwise, a single call is made to peak_trace() where both the peak and troughs in sobelsig are detected and traced.

Optionally, the code will match and compare the traces found and fit by peak_trace() to the original traces. If the RMS difference between the matched traces is large, they can be removed (see trace_rms_tol in EdgeTracePar).

Note that this effectively reinstantiates much of the object attributes, including traceid, edge_cen, edge_err, edge_msk, edge_img, edge_fit, and fittype.

Used parameters from par (EdgeTracePar) are left_right_pca, edge_thresh, smash_range, edge_detect_clip, trace_median_frac, trace_thresh, trace_rms_tol, fit_function, fit_order, fwhm_uniform, niter_uniform, fwhm_gaussian, niter_gaussian, fit_maxdev, and fit_maxiter.

Parameters:
  • rebuild_pca (bool, optional) – This method fundamentally resets the trace data, meaning that the PCA is no longer valid. Use this boolean to have the method rebuild the PCA based on the refined traces. Note that the PCA is not then used to reset the fitted trace data; i.e., edge_fit remains based on the output of peak_trace().

  • debug (bool, optional) – Run in debug mode.

Raises:

PypeItError – Raised if pca is not defined.

predict_traces(edge_cen, side=None)[source]

Use the PCA decomposition to predict traces.

The PCA decomposition must be available via pca or left_pca and right_pca; see build_pca(). This is a convenience method to handle the PCA predictions given that left and right traces can be decomposed separately or simultaneously.

Parameters:
  • edge_cen (float, numpy.ndarray) – A single value or 1D array with the spatial location (column) for 1 or more traces to predict. The predicted traces will pass through these spatial positions (columns) and the reference spectral row set for the PCA decomposition; see build_pca().

  • side (float, numpy.ndarray, optional) – A single value or 1D integer array indicating the edge side to be predicted; -1 for left and 1 for right. Must be the same length as edge_cen. This is only used if the PCA is side-dependent, and must be provided in the case that it is (see left_right_pca in par).

Returns:

A 1D or 2D array of size nspec by the length of the position array provided. If a single coordinate is provided, a single trace vector is returned.

Return type:

numpy.ndarray

qa_plot(fileroot=None, min_spat=20)[source]

Build a series of QA plots showing the edge traces.

Warning

This method is slow.

Parameters:
  • fileroot (str, optional) – Root name for the output files. The number of output files depends on the layout and the number of traces found. If None, plots are displayed interactively.

  • min_spat (int, optional) – Minimum number of spectral pixels to plot for each trace. If None, set to twice the difference between the minimum and maximum centroid of the plotted trace.

rectify(flux, bpm=None, extract_width=None, mask_threshold=0.5, side='left')[source]

” Rectify the provided image based on the current edge trace PCA model.

The is primarily a wrapper for pypeit.sampling.rectify_image(); see its documentation for more detail.

Used parameters from par (EdgeTracePar) are left_right_pca.

Parameters:
  • flux (numpy.ndarray) – The 2D image to rectify. Its shape should match the image used to construct the edge traces: \((N_{\rm spec}, N_{\rm spat})\).

  • bpm (numpy.ndarray, optional) – Boolean bad-pixel mask for pixels to ignore in input image. If None, no pixels are masked in the rectification. If provided, shape must match flux.

  • extract_width (float, optional) – The width of the extraction aperture to use for the image rectification. When using extraction to rectify the image, flux conservation is not as accurate. If None, the image rectification is performed using pypeit.sampling.Resample along each row (spectral position).

  • mask_threshold (float, optional) – Either due to mask or the bounds of the provided flux, pixels in the rectified image may not be fully covered by valid pixels in flux. Pixels in the output image with less than this fractional coverage by input pixels are flagged in the output.

  • side (str, optional) – If the PCA decomposition was performed for the left and right traces separately, this selects the PCA to use when constructing the rectification coordinate system. This is ignored if both were used in the PCA.

Returns:

Two numpy.ndarray objects are returned both with shape \((N_{\rm spec}, N_{\rm spat})\), the rectified image and its boolean bad-pixel mask.

remove_traces(indx, resort=True, rebuild_pca=False)[source]

Remove a set of traces.

This method directly alters the attributes containing the trace data.

Warning

If the traces are left-right synchronized, you should decide how you want to treat the traces as pairs. See synced_selection(). For example, to remove both traces in a pair if either are selected, run:

edges.remove_traces(edges.synced_selection(indx, mode='both'))
Parameters:
  • indx (array-like) – The boolean array with the traces to remove. Length must be \((N_{\rm trace},)\).

  • resort (bool, optional) – Re-sort the traces and trace IDs to be sequential in the spatial direction. See spatial_sort().

  • rebuild_pca (bool, optional) – If the pca exists, rebuild it using the new traces and the previous parameter set. If False, any existing PCA model is removed.

rm_user_traces(rm_traces)[source]

Parse the user input traces to remove

Parameters:

rm_user_traces (list) – y_spec, x_spat pairs

Returns:

show(include_error=False, thin=10, in_ginga=False, include_img=True, include_sobel=False, img_buffer=100, flag=None, idlabel=True, slits=None, original=False, title=None)[source]

Show a scatter plot of the current trace data and fit, if it’s available.

Parameters:
  • include_error (bool, optional) – Show the errors on the measurements

  • thin (int, optional) – Thin the data plotted by plotting every thin measurement in the spectral direction. Default is to plot all data; to show every other datum, set thin=2.

  • in_ginga (bool, optional) – Display the trace data in a ginga. If a ginga window is not open, this instantiates one; otherwise, the existing ginga window is cleared. The trace image is shown in one window and the sobel image is shown in a second. The edge traces themselves are shown in both windows. Any masking is ignored except that any traces that are fully masked (see fully_masked_traces()) are not shown. If a SlitTraceSet object is not provided, the data shown is the modeled results for the trace if it exists, and the measured trace locations otherwise.

  • include_img (bool, optional) – Overlay the trace data on the trace image; mutually exclusive with include_sobel.

  • include_sobel (bool, optional) – Overlay the trace data on the Sobel-filtered trace image; mutually exclusive with include_img.

  • img_buffer (int, optional) – Buffer for plot window with respect to trace image size.

  • flag (str, list, optional) – One or more flags to use when selecting masked data to include in the plot; unmasked data are always plotted. If None, no masked data is plotted. If ‘any’, trace locations masked for any reason are plotted. Otherwise, only trace locations flagged by the specified bits are plotted.

  • idlabel (bool, optional) – Label each trace by their ID numbers. If the data is from an echelle spectrograph and the echelle order has been matched to the synchronized slit-trace pairs, the ID labels are based on the order matching.

  • slits (SlitTraceSet, optional) – Slits to plot instead of those kept internally. Note that, if this is provided, the “modeled” and “measured” slit locations are identical.

  • original (bool, optional) – When slits are provided and tweaked slits are available, this selects which traces to show. If True, show the original slit traces; if False, show the tweaked ones.

  • title (str, optional) – Title for the matplotlib plot (ignored if shown in ginga). If None, plot is not given a title.

slit_spatial_center(normalized=True, spec=None, use_center=False, include_box=False)[source]

Return the spatial coordinate of the center of each slit.

The slit edges must be left-right synchronized.

Parameters:
  • normalized (bool, optional) – Return coordinates normalized by the size of the detector.

  • spec (int, optional) – Spectral position (row) at which to return the spatial position. If None, use the PCA reference row if a PCA exists or the central row (i.e., self.nspec//2), otherwise.

  • use_center (bool, optional) – Use the measured centroids to define the slit edges even if the slit edges have been otherwise modeled.

  • include_box (bool, optional) – Include box slits in the calculated coordinates.

Returns:

Spatial coordinates of the slit centers in pixels or in fractions of the detector. Masked locations are for bad/excluded slits.

Return type:

numpy.ma.MaskedArray

spatial_sort(use_mean=False, use_fit=True)[source]

Sort the traces spatially.

The coordinates used for the sorting are either the measured centroids or the fitted parameterization (see use_fit). The fiducial coordinates that are sorted are either the mean of the unmasked coordinates over all spectral rows or the unmasked coordinates at a specified reference spectral row (see use_mean).

The trace IDs are also reassigned to be sorted spatially; i.e., the trace IDs for three synced slits would be [-1, 1, -2, 2, -3, 3].

All attributes are edited in-place.

Parameters:
  • use_mean (bool, optional) – Sort according to the mean of the masked spatial positions. If False, the spatial position at a reference spectral row is used, where the reference spectral row is either the same as used by the PCA (if available) or the result of most_common_trace_row() using the current trace mask.

  • use_fit (bool, optional) – Sort according to the fit positions instead of the measured positions. Otherwise, only use the fit positions if they’re available and the measured location is masked.

sync(rebuild_pca=True, debug=False)[source]

Match left and right edge traces to construct slit edge pairs.

Synchronization of the slits proceeds as follows:

  • Any fully masked traces are removed using clean_traces().

  • If that operation removes all traces (see is_empty()), two traces are added at the edge of the detector.

  • The spatial order of the traces is sorted (see spatial_sort()).

  • At this point, if the traces are already synced (see is_synced()), the synchronization is checked using check_synced() and the method finishes.

  • If the traces need to be synced, the method determines which traces need an inserted trace to make a pair (see _get_insert_locations()).

  • The reference locations for the new traces are determined by _get_reference_locations().

  • The shape of the inserted traces is set according to the sync_predict parameter in par. The shape is either predicted by the PCA decomposition or taken to be exactly the same shape as the nearest left or right edge.

  • These new traces are then inserted using insert_traces() and flagged as being inserted by the synchronization process.

  • Assuming there wasn’t an error in the insertion scheme, the synchronized traces are then checked by check_synced() and the method finishes.

Before the last step above, the synchronization is checked. Certain corner cases can lead to catastrophic errors in where the inserted traces are placed such that the left-right ordering of the synchronized traces is incorrect and exception is raised.

Used parameters from par (EdgeTracePar) are det_buffer, left_right_pca, and sync_predict.

Warning

Synchronizing the left and right edges requires that traces that are fully masked as bad must be removed (this does not include traces that are “masked” as having been deliberately inserted), and these traces are removed regardless of the user-specified clip in par.

Parameters:
  • rebuild_pca (bool, optional) – If the pca exists and traces are removed (see check_synced()), rebuild the PCA using the new traces and the previous parameter set. Note that inserted traces are not included in the PCA decomposition.

  • debug (bool, optional) – Run in debug mode.

Returns:

Returns the status of the syncing. True means success.

Return type:

bool

synced_selection(indx, mode='ignore', assume_synced=False)[source]

Coordinate selection of traces based on their synchronization status.

Parameters:
  • indx (numpy.ndarray) – Boolean vector selecting a set of traces. The shape should be \((N_{\rm trace},)\).

  • mode (str, optional) –

    If the traces are left-right synchronized (see sync() and is_synced()), this method dictates how the selection of traces should be changed to deal with edges that have been synchronized into left-right pairs. Methods are:

    • 'ignore': Ignore the synchronization and just return the input.

    • 'both': If one of the traces in a pair is selected, select both.

    • 'neither': If only one of the traces in a pair is selected, deselect both.

  • assume_synced (bool, optional) – Assume the set of traces is synced when identifying which traces to select/de-select if mode is either 'both' or 'neither'. This essentially ignores the sync status of the traces and will fault if the number of traces is not even. If False, the sync status is checked using is_synced().

Raises:

PypeItError – Raised if the input indx doesn’t have the correct length, if mode is not 'ignore' and the edges have not been left-right synchronized and assume_synced is False, or if the provided mode is not valid.

Returns:

The updated boolean vector selecting traces based on the input mode and the synchronization status.

Return type:

numpy.ndarray

to_hdu(**kwargs)[source]

Construct the HDUList to write.

This overrides pypeit.datamodel.DataContainer.to_hdu() and has the same calling sequence (i.e., kwargs are passed directly to the base class method). This is needed to deal with the multiple TracePCA objects included in the output.

Returns:

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

Return type:

list, astropy.io.fits.HDUList

trace_pixels_off_detector(cen=None)[source]

Determine if the traces land outside of the detector and the defined edge buffer (det_buffer in par).

Parameters:

cen (numpy.ndarray, optional) – Array with the trace locations. Shape is \((N_{\rm spec}, N_{\rm trace})\). If None, uses edge_cen.

Returns:

Boolean array with shape \((N_{\rm spec}, N_{\rm trace})\) selecting pixels that are (True) or are not (False) “off” the detector.

Return type:

numpy.ndarray

Raises:

PypeItError – Raised if cen is not provided and edge_cen is None.

version = '1.0.1'

DataContainer datamodel version.