pypeit.spectrographs.slitmask module

Module to define the SlitMask class

class pypeit.spectrographs.slitmask.SlitMask(corners, slitid=None, align=None, science=None, onsky=None, objects=None, posx_pa=None, object_names=None, mask_radec=None)[source]

Bases: object

Generic class for a slit mask that holds the slit positions and IDs.

By default no mask bits are set. Only altered if align or science arguments are provided.

Parameters:
  • corners (array-like) – A list or numpy.ndarray with the list of coordinates. The object must be 2 or 3 dimensional: 1st axis is the slit, 2nd axis are the 4 corners, 3rd axis is the x and y coordinate for the corner. If two dimensional, class assumes there is only one slit. The size of the last two dimensions must always be 4 and 2, respectively. The input order is expected to be clockwise, starting from the top-right corner; i.e., the order is top-right (high x, low y), bottom-right (low x, low y), bottom-left (low x, high y), top-left (high x, high y). The x coordinates are along the spatial direction and the y coordinates are long the spectral direction. The computed length (difference in x), width (difference in y), and position angle (Cartesian angle) is relative to this assumption.

  • slitid (int, array-like, optional) – A list of unique integer IDs for each slit. If None, just a running 0-indexed list is used. Can be a single integer or have shape \((N_{\rm slit},)\).

  • align (bool, numpy.ndarray, optional) – Indicates slit(s) used for on-sky alignment. Can be a single boolean or have shape \((N_{\rm slit},)\).

  • science (bool, numpy.ndarray, optional) – Indicates slit(s) include a target of scientific interest. Can be a single boolean or have shape \((N_{\rm slit},)\).

  • onsky (numpy.ndarray, optional) – 1D or 2D array with on-sky metrics for each slit. The shape of the array must be \((5,)\) or \((N_{\rm slit},5)\), and the order along rows must match slit ID order (see slitid). The five metrics per slit are (1-2) right ascension and declination of the slit center, (3-4) the slit length and width in arcseconds, and (5) the position angle of the slit from N through E in degrees.

  • objects (numpy.ndarray, optional) – List of objects observed as a 1D or 2D array with shape \((9,)\) or \((N_{\rm obj},9)\). The nine elements for each object is the slit id, the object ID, the right ascension and declination of the target, the object name, the object magnitude and band, and the object top and bottom distances from the slit edges. The order of the objects does not have to match that of the slit IDs. Also, there can be slits without objects and slits with multiple objects; however, objects cannot be provided that are not in any slit (i.e., the slit IDs in the first column of this array have to be valid).

corners

See above.

Type:

numpy.ndarray

slitid

See slitid above.

Type:

numpy.ndarray

mask

Mask bits selecting the type of slit.

Type:

numpy.ndarray

onsky

See above.

Type:

numpy.ndarray

objects

See above.

Type:

numpy.ndarray

slitindx

The index that maps from the slit data to the object data. For example:

objslitdb = self.onsky[self.slitindx]

provides the onsky slit data for each object with a shape matched to the relevant entry in self.objects.

Type:

numpy.ndarray

center

The geometric center of each slit.

Type:

numpy.ndarray

top

The top coordinate of each slit.

Type:

numpy.ndarray

bottom

The bottom coordinate of each slit.

Type:

numpy.ndarray

length

The slit length.

Type:

numpy.ndarray

width

The slit width.

Type:

numpy.ndarray

pa

The cartesian rotation angle of the slit in degrees.

Type:

numpy.ndarray

mask_radec

RA, Dec (deg) of the pointing of the mask (approximate center)

Type:

tuple, optional

posx_pa

Sky PA that points to positive x (spatial) on the detector

Type:

float

negx_pa

Sky PA that points to negative x (spatial) on the detector

Type:

float

object_names

Object names

Type:

numpy.ndarray

Raises:

ValueError – Raised if the shape of the input corners array is incorrect or if the number of slit IDs does not match the number of slits provided.

property alignment_slit

Boolean array selecting the alignment slits.

bitmask = <pypeit.spectrographs.slitmask.SlitMaskBitMask object>
is_alignment(i)[source]

Check if specific slit is an alignment slit.

is_science(i)[source]

Check if specific slit should have a science target.

property nslits

The number of slits.

property science_slit

Boolean array selecting the slits with science targets.

class pypeit.spectrographs.slitmask.SlitMaskBitMask[source]

Bases: BitMask

Mask bits used for slit mask design data.

class pypeit.spectrographs.slitmask.SlitRegister(trace_spat, mask_spat, trace_mask=None, guess=[0.0, 1.0], fix=[False, False], bounds=None, penalty=False, maxiter=1, maxsep=None, sigma=5, debug=False, fit=False)[source]

Bases: object

Match trace and slit mask positions using a linear relation.

The class performs a least-squares fit for the following:

\[X_{\rm trace} = o + s\ X_{\rm slit},\]

where \(s\) is a scale factor that nominally converts the units of the input slit-mask positions to pixels and \(o\) is the offset in pixels.

Assuming the trace coordinates are in pixels and the mask coordinates are in mm at the focal plane, the nominal scale should be the plate-scale of the telescope (arcsec/mm) divided by the plate-scale of the detector (arcsec/pixel). The nominal offset is then the negative of the coordinate of the relevant detector edge relative to the detector center in pixels.

Parameters:
  • trace_spat (array-like) – Trace coordinates in the spatial direction, typically provided in pixel index units (but not necessarily index integers).

  • mask_spat (array-like) – Slit-mask coordinates in the spatial direction, typically provided in mm in the focal plane.

  • guess (array-like, optional) – The initial guess for the fit parameters. Parameters are currently an offset (guess[0]) and a scale (guess[1]).

  • fix (array-like, optional) – An array of booleans indicating if the provided guess parameters should be fixed during the fit.

  • bounds (array-like, optional) – An array of the lower and upper bounds for the parameters. Must have shape \((N_{\rm par},2)\), even if some parameters are fixed; currently \(N_{\rm par}=2\). If None, no bounds are imposed (specifically, bounds are set to \(\pm``numpy.inf\).)

  • penalty (bool, optional) – Include a logarithmic penalty for slits that are matched to multiple slits.

  • maxiter (int, optional) – Maximum number of fit iterations to perform. If None, rejection iterations are performed until no points are rejected. If 1, only a single fit is performed without any rejection; i.e., the number of rejection iterations is maxiter-1.

  • maxsep (float, optional) – The maximum allowed separation between the calibrated coordinates of the designed slit position in pixels and the matched trace. If None, rejection is done iteratively using 5-sigma clipping.

  • debug (bool, optional) – Show a plot of the fit residuals after each iteration.

  • fit (bool, optional) – Perform the fit based on the input. If False, all optional entries are ignored and the user has to call the find_best_match() function with those same entries to perform the fit.

trace_spat

Trace coordinates in the spatial direction.

Type:

numpy.ndarray

mask_spat

Slit-mask coordinates in the spatial direction.

Type:

numpy.ndarray

guess_par

The guess model parameters.

Type:

numpy.ndarray

fit_par

Flag that the parameters should be fit.

Type:

numpy.ndarray

bounds

The boundaries imposed on the fit parameters.

Type:

tuple

par

The full parameter set, including any fixed parameters.

Type:

numpy.ndarray

penalty

Flag to apply the penalty function during the fit.

Type:

bool

match_coo

The best matching coordinates of the slit mask design coordinates to the slit trace pixel positions.

Type:

numpy.ndarray

match_index

Indices in the mask_spat that are best matched to the trace_spat coordinates.

Type:

numpy.ndarray

match_separation

Signed difference between the best-matching trace and mask positions in trace units; negative values mean the best-matching mask position is larger than the trace position.

Type:

numpy.ndarray

Raises:

NotImplementedError – Raised if the spatial positions are not 1D arrays.

_setup_to_fit(guess, fix, bounds, penalty, maxsep, sigma)[source]

Setup the necessary attributes for the fit.

find_best_match(guess=[0.0, 1.0], fix=[False, False], bounds=None, penalty=False, maxiter=1, maxsep=None, sigma=5, debug=False)[source]

Find the best match between the trace and slit-mask positions.

Populates match_coo, match_separation, and match_index; the latter is also returned.

Parameters:
  • guess (array-like, optional) – The initial guess for the fit parameters. Parameters are currently an offset (guess[0]) and a scale (guess[1]).

  • fix (array-like, optional) – An array of booleans indicating if the provided guess parameters should be fixed during the fit.

  • bounds (array-like, optional) – An array of the lower and upper bounds for the parameters. Must have shape \((N_{\rm par},2)\), even if some parameters are fixed; currently \(N_{\rm par}=2\). If None, no bounds are imposed (specifically, bounds are set to \(\pm``numpy.inf\).)

  • penalty (bool, optional) – Include a logarithmic penalty for slits that are matched to multiple slits.

  • maxiter (int, optional) – Maximum number of fit iterations to perform. If None, rejection iterations are performed until no points are rejected. If 1, only a single fit is performed without any rejection; i.e., the number of rejection iterations is maxiter-1.

  • maxsep (float, optional) – The maximum allowed separation between the calibrated coordinates of the designed slit position in pixels and the matched trace. If None, rejection is done iteratively using sigma clipping.

  • sigma (float, optional) – The sigma value to use for rejection. If None, it will use the default set by astropy.stats.sigma_clipped_stats.

  • debug (bool, optional) – Show a plot of the fit residuals after each iteration.

Returns:

The index of the slit mask position matched to each trace position.

Return type:

numpy.ndarray

mask_to_trace_coo(par=None)[source]

Compute the mask positions in the trace coordinate system.

This is the core method that converts the mask coordinates to the trace coordinates; any other method that needs these data call this method.

Parameters:

par (array-like, optional) –

The list of (free) parameters. If any parameters are fixed, they should not be included, such that the full parameter set is:

self.par[self.fit_par] = par

If None, use par.

Returns:

The expected pixel positions of the trace given the mask positions and the model parameters.

Return type:

numpy.ndarray

match(par=None, unique=False)[source]

Match each trace to the nearest slit position based on the provided or internal fit parameters (using mask_to_trace_coo()).

Note

Even though this method returns the information identically meant for match_coo, match_separation, and match_index, these internals are not altered by this function.

Parameters:
  • par (numpy.ndarray, optional) – The parameter vector. See mask_to_trace_coo().

  • unique (bool, optional) – Force the set of matched indices to be unique. This can lead to large separations, which can be used to find traced slits that are inconsistent with the mask design. If the function runs out of available mask slits to match to, it will set the remaining indices to -1 and the separation to 9999.

Returns:

(1) the mask coordinates in the detector frame for all mask coordinates; (2) the signed minimum difference between each trace and any mask position; and (3) the index in the array of mask coordinates that provide the best match to each trace coordinate. The indices the latter are only unique if requested.

Return type:

Returns three numpy.ndarray objects

minimum_separation(par=None)[source]

Return the minimum trace and mask separation for each trace.

This is the method that returns the residuals that are minimized by scipy.optimize.least_squares. Those residuals are provided by a call to match() without forcing the matching pairs to be unique.

The minimum separation can be penalized (if penalty is True), which multiplies each separation by \(2^{dN}\) if there are non-unique matches slit-to-trace matches. Here, \(dN\) is the difference between the number of traces and the number of uniquely matched mask positions; i.e., if two traces are matched to the same mask position, the separation is increase by a factor of 2.

Residuals are only returned for the unmasked trace positions; see trace_mask.

Parameters:

par (numpy.ndarray, optional) – The parameter vector. See mask_to_trace_coo().

Returns:

The signed difference between the trace coordinates and its most closely associated slit position (trace-mask).

Return type:

numpy.ndarray

show(par=None, maxsep=None, sigma=None, unique=True, minmax=None, synced=False)[source]

Plot the fit residuals.

Parameters:
  • par (numpy.ndarray, optional) – The parameter vector. See mask_to_trace_coo().

  • maxsep (float, optional) – The maximum allowed separation between the calibrated coordinates of the designed slit position in pixels and the matched trace. If None, use maxsep; see find_best_match().

  • sigma (float, optional) – The sigma value to use for rejection. If None, use sigma; see find_best_match().

  • unique (bool, optional) – Force the set of matched indices to be unique; see match().

  • minmax (array-like, optional) – A two-element array with the minimum and maximum coordinate value to match to the trace data; see trace_mismatch().

  • synced (bool, optional) – The mask coordinates being matched to are synced left-to-right in adjacent pairs. I.e., the indices of left edges are all even and the indices of all right edges are odd.

trace_mismatch(maxsep=None, sigma=None, minmax=None, synced=False)[source]

Return the mismatches between the mask and trace positions.

Based on the best-fitting (or fixed) offset and scale parameters, match() is executed, forcing the slit-mask and trace positions pairs to be uniquely matched.

The set of slit-mask positions without a matching trace are identified by finding those slits in the range relevant to the list of trace coordinates (see minmax), but without a matching trace index.

Todo

explain synced adjustment

The set of mask-to-trace matches are identified as “bad” if they meet any of the following criteria:

  • The trace has not been masked (see trace_mask)

  • A unique match could not be found (see match())

  • The absolute value of the separation is larger than the provided maxsep (when maxsep is not None).

  • The separation is rejected by a sigma-clipping (see sigma)

Note that there is currently no argument that disables the determination of bad traces. However, bad traces are simply returned by the method; this function changes none of the class attributes.

Parameters:
  • maxsep (float, optional) – The maximum allowed separation between the calibrated coordinates of the designed slit position in pixels and the matched trace. If None, use maxsep; see find_best_match().

  • sigma (float, optional) – The sigma value to use for rejection. If None, use sigma; see find_best_match().

  • minmax (array-like, optional) – A two-element array with the minimum and maximum coordinate value to match to the trace data. If None, this is determined from trace_spat and the standard deviation of the fit residuals.

  • synced (bool, optional) – The mask coordinates being matched to are synced left-to-right in adjacent pairs. I.e., the indices of left edges are all even and the indices of all right edges are odd.

Returns:

(1) the indices of mask positions without a matching trace position and (2) the list of trace positions identified as “bad.”

Return type:

Two numpy.ndarray objects are returned

pypeit.spectrographs.slitmask.build_slit_function(edges, size=None, oversample=1, sigma=None)[source]

Construct a unit normalized slit function

pypeit.spectrographs.slitmask.correct_slitpa(slitpa, maskpa)[source]

Flip 180 degree the slit PA if the value recorded in the slitmask design is more than +/-90 degree from the slitmask PA.

Parameters:
  • slitpa (float or numpy.ndarray) – position angle of the slits.

  • maskpa – (float): position angle of the slitmask.

Returns:

flipped slitpa, if it is more than +/-90 from the maskpa, otherwise unchanged slitpa.

Return type:

float or numpy.ndarray

pypeit.spectrographs.slitmask.load_keck_deimoslris(filename: str, instr: str)[source]

Load up the mask design info from the header of the file provided

Parameters:
  • filename (str) –

  • instr (str) – Name of spectrograph Allowed are keck_lris_xxx, keck_deimos

Returns:

[description]

Return type:

[type]

pypeit.spectrographs.slitmask.positive_pa(pa: float)[source]

Modify input pa to be positive (0-360)

Parameters:

pa (float) – [description]

Returns:

[description]

Return type:

[type]

pypeit.spectrographs.slitmask.slit_function_length(edges, oversample=1)[source]
pypeit.spectrographs.slitmask.xc_trace(x_trace, x_design, pix_per_mm)[source]

Use a cross-correlation to find the offset