pypeit.core.scattlight module
Module for sky subtraction
- pypeit.core.scattlight.fine_correction(frame, bpm, offslitmask, method='median', polyord=2, debug=False)[source]
Calculate a fine correction to the residual scattered light of the input frame.
- Parameters:
frame (numpy.ndarray) – Raw 2D data frame (nspec, nspat) to be used to compute the fine correction of the scattered light. This frame should be the raw frame, minus the first estimate of the scattered light that has been derived from the
scattered_light_model()
function.bpm (numpy.ndarray) – 2D boolean array indicating the bad pixels (True=bad), same shape as frame
offslitmask (numpy.ndarray) – A boolean mask indicating the pixels that are on/off the slit (True = off the slit), same shape as frame
method (
str
, optional) –- Method to use to determine the fine correction to the scattered light. Options are:
’median’: Use the median of the off-slit pixels to determine the scattered light
- ’poly’: Use a polynomial fit to the off-slit pixels to determine the scattered light. If this
option is chosen, the polynomial order must be specified. See polyord below.
polyord (
int
, optional) – Polynomial order to use for fitting the residual scattered light in the spatial direction.debug (
bool
, optional) – If True, debug the final model fit that’s been output
- Returns:
scatt_img – A 2D image (nspec, nspat) of the fine correction to the scattered light determined from the input frame.
- Return type:
- pypeit.core.scattlight.mask_slit_regions(offslitmask, centrace, mask_regions=None)[source]
Exclude some user-specified inter-slit regions for the fine correction to the scattered light determination
- Parameters:
offslitmask (numpy.ndarray) – A boolean mask indicating the pixels that are on/off the slit (True = off the slit)
centrace (numpy.ndarray) – A 2D array, shape is (nspec, nslit), containing the central trace of each slit.
mask_regions (
int
,list
, optional) – A list of regions that should be excluded in the fine correction to the scattered light. A zero value indicates that all pixels left of the first slit will be masked. A value of one indicates that all pixels between the first and second slit will be masked, and so forth. A list of these integers can be supplied. If None, no inter-slit regions will be excluded. If an integer is specified, just one region will be excluded.
- Returns:
good_mask – A 2D boolean array indicating the pixels that are suitable to use for the fine scattered light correction.
- Return type:
- pypeit.core.scattlight.pad_frame(_frame, detpad=300)[source]
Clean the edges of the input frame and then pad the frame to avoid edge effects.
- Parameters:
_frame (numpy.ndarray) – Frame to be padded
detpad (int) – Number of pixels to pad the frame on each side
- Returns:
_frame_pad – Padded frame
- Return type:
- pypeit.core.scattlight.scattered_light(frame, bpm, offslitmask, x0, bounds, detpad=300, debug=False)[source]
Calculate a model of the scattered light of the input frame.
- Parameters:
frame (numpy.ndarray) – Raw 2D data frame (nspec, nspat) to be used to compute the scattered light.
bpm (numpy.ndarray) – 2D boolean array indicating the bad pixels (True=bad), same shape as frame
offslitmask (numpy.ndarray) – A boolean mask indicating the pixels that are on/off the slit (True = off the slit), same shape as frame
x0 (numpy.ndarray) – A 1D array containing the best-fitting model parameters
bounds (
tuple
) – A tuple of two elements, containing two numpy.ndarray of the same length as x0. These two arrays contain the lower (first element of the tuple) and upper (second element of the tuple) bounds to consider on the scattered light model parameters.debug (
bool
, optional) – If True, debug the final model fit that’s been output
- Returns:
scatt_img (numpy.ndarray) – A 2D image of the scattered light determined from the input frame. Alternatively, if a constant value is used, a constant floating point value can be returned as well.
modelpar (numpy.ndarray) – A 1D array containing the best-fitting model parameters
success (
bool
) – True if the fit was successful, False otherwise
- pypeit.core.scattlight.scattered_light_model(param, img)[source]
Model used to calculate the scattered light.
The current model to generate the scattered light is a shifted, scaled, and blurred version of the data recorded on the detector. We allow for a shift in the spatial and spectral direction, and the blurring kernel can have different widths in the spectral and spatial directions. Finally, the model includes a polynomial correction in the spectral direction to account for a change in the spectral efficiency of the scattered light relative to the detector image.
- Parameters:
param (numpy.ndarray) –
Model parameters that determine the scattered light based on the input img. Here is a list of the individual parameter meanings:
param[0] = Gaussian kernel width in the spectral direction
param[1] = Gaussian kernel width in the spatial direction
param[2] = Lorentzian kernel width in the spectral direction
param[3] = Lorentzian kernel width in the spatial direction
param[4] = Pixel shift of the scattered light in the spectral direction
param[5] = Pixel shift of the scattered light in the spatial direction
param[6] = Zoom factor of the scattered light (~1)
param[7] = constant offset for scattered light (independent of img)
param[8] = Kernel angle
- param[9] = Relative importance of Gaussian vs Lorentzian.
0 < value < 1 means Lorentzian is weighted more value > 1 means Gaussian is weighted more.
param[10:] = Polynomial scaling coefficients
img (numpy.ndarray) – Raw image that you want to compute the scattered light model. shape is (nspec, nspat) Model used to calculate the scattered light. This function is used to generate a model of the scattered light, based on a set of model parameters that have been optimized using self.scattered_light().
- Returns:
model – Model of the scattered light for the input
- Return type:
- pypeit.core.scattlight.scattered_light_model_pad(param, img, detpad=300)[source]
Construct a scattered light model for the input image, with the model parameters defined by param. This function is used to generate a model of the scattered light, based on a set of model parameters that have first been optimized using scattered_light(). The model is generated on a padded version of the input image, and then trimmed to match the input image size.
- Parameters:
param (numpy.ndarray) –
Model parameters that determine the scattered light based on the input img. Here is a list of the individual parameter meanings:
param[0] = Gaussian kernel width in the spectral direction
param[1] = Gaussian kernel width in the spatial direction
param[2] = Lorentzian kernel width in the spectral direction
param[3] = Lorentzian kernel width in the spatial direction
param[4] = Pixel shift of the scattered light in the spectral direction
param[5] = Pixel shift of the scattered light in the spatial direction
param[6] = Zoom factor of the scattered light (~1)
param[7] = constant offset for scattered light (independent of img)
param[8] = Kernel angle
- param[9] = Relative importance of Gaussian vs Lorentzian.
0 < value < 1 means Lorentzian is weighted more value > 1 means Gaussian is weighted more.
param[10:] = Polynomial scaling coefficients
img (numpy.ndarray) – Image used to generate the scattered light model
detpad (int) – Number of pixels to pad the frame on each side
- Returns:
model – Model of the scattered light for the input
- Return type:
- pypeit.core.scattlight.scattlight_resid(param, wpix, img)[source]
Residual function used to optimize the model parameters
- Parameters:
param (numpy.ndarray) – 1D array of model parameters to use for the fitting function.
wpix (tuple) – A tuple containing the x,y coordinates of the pixels in img to be used for computing the residual.
img (numpy.ndarray) – Data image to be used to compute the residual. Shape is (nspec, nspat)
- Returns:
resid – A 1D vector of the residuals
- Return type: