pypeit.flatfield module
Implements the flat-field class.
- class pypeit.flatfield.FlatField(rawflatimg, spectrograph, flatpar, slits, wavetilts=None, wv_calib=None, slitless=False, spat_illum_only=False, qa_path=None, calib_key=None)[source]
Bases:
object
Builds pixel-level flat-field and the illumination flat-field.
For the primary methods, see
run()
.- Parameters:
rawflatimg (
PypeItImage
) – Processed, combined set of pixelflat imagesspectrograph (
Spectrograph
) – The Spectrograph instance that sets the instrument used to take the observations.flatpar (
FlatFieldPar
) – User-level parameters for constructing the flat-field corrections. If None, the default parameters are used.slits (
SlitTraceSet
) – The current slit traces.wavetilts (
WaveTilts
, optional) – The current fit to the wavelength tilts. I can be None, for example, if slitless is True.wv_calib (
WaveCalib
, optional) – Wavelength calibration object. It can be None, for example, if slitless is True.slitless (bool, optional) – True if the input rawflatimg is a slitless flat. Default is False.
spat_illum_only (bool, optional) – Only perform the spatial illumination calculation, and ignore the 2D bspline fit. This should only be set to true if you want the spatial illumination profile only. If you want to simultaneously generate a pixel flat and a spatial illumination profile from the same input, this should be False (which is the default).
qa_path (str, optional) – Path to QA directory
- rawflatimg
- Type:
- mspixelflat
Normalized flat
- Type:
- msillumflat
Illumination flat
- Type:
- flat_model
Model of the flat
- Type:
- spec_illum
Image of the relative spectral illumination for a multislit spectrograph
- Type:
- extract_structure(rawflat_orig, slit_trim=3)[source]
Generate a relative scaling image for a slicer IFU. All slits are scaled relative to a reference slit, specified in the spectrograph settings file.
- Parameters:
rawflat_orig (numpy.ndarray) – The original raw image of the flatfield
slit_trim (int, optional) – Trim the slit edges by this number of pixels during the fitting. Note that the fit will be evaluated on the pixels indicated by onslit_tweak. A positive number trims the slit edges, a negative number pads the slit edges.
- Returns:
ff_struct – An image containing the detector structure (i.e. the raw flatfield image divided by the spectral and spatial illumination profile fits).
- Return type:
- fit(spat_illum_only=False, doqa=True, debug=False)[source]
Construct a model of the flat-field image.
For this method to work,
rawflatimg
must have been previously constructed.The method loops through all slits provided by the
slits
object, except those that have been masked (i.e., slits withself.slits.mask == True
are skipped). For each slit:Collapse the flat-field data spatially using the wavelength coordinates provided by the fit to the arc-line traces (
WaveTilts
), and fit the result with a bspline. This provides the spatially-averaged spectral response of the instrument. The data used in the fit is trimmed toward the slit spatial center via theslit_trim
parameter inflatpar
.Use the bspline fit to construct and normalize out the spectral response.
Collapse the normalized flat-field data spatially using a coordinate system defined by the left slit edge. The data included in the spatial (illumination) profile calculation is expanded beyond the nominal slit edges using the
slit_illum_pad
parameter inflatpar
. The raw, collapsed data is then median filtered (seespat_samp
inflatpar
) and Gaussian filtered; seepypeit.core.flat.illum_filter()
. This creates an empirical, highly smoothed representation of the illumination profile that is fit with a bspline using thespatial_fit()
method. The construction of the empirical illumination profile (i.e., before the bspline fitting) can be done iteratively, where each iteration sigma-clips outliers; see theillum_iter
andillum_rej
parameters inflatpar
andpypeit.core.flat.construct_illum_profile()
.If requested, the 1D illumination profile is used to “tweak” the slit edges by offsetting them to a threshold of the illumination peak to either side of the slit center (see
tweak_slits_thresh
inflatpar
), up to a maximum allowed shift from the existing slit edge (seetweak_slits_maxfrac
inflatpar
). Seetweak_slit_edges()
. If tweaked, thespatial_fit()
is repeated to place it on the tweaked slits reference frame.Use the bspline fit to construct the 2D illumination image (
msillumflat
) and normalize out the spatial response.Fit the residuals of the flat-field data that has been independently normalized for its spectral and spatial response with a 2D bspline-polynomial fit. The order of the polynomial has been optimized via experimentation; it can be changed but you should use extreme caution when doing so (see
twod_fit_npoly
). The multiplication of the 2D spectral response, 2D spatial response, and joint 2D fit to the high-order residuals define the final flat model (flat_model
).Finally, the pixel-to-pixel response of the instrument is defined as the ratio of the raw flat data to the best-fitting flat-field model (
mspixelflat
)
This method is the primary method that builds the
FlatField
instance, constructingmspixelflat
,msillumflat
, andflat_model
. All of these attributes are altered internally. If the slit edges are to be tweaked using the 1D illumination profile (tweak_slits
inflatpar
), the tweaked slit edge arrays in the internalSlitTraceSet
object,slits
, are also altered.Used parameters from
flatpar
(FlatFieldPar
) arespec_samp_fine
,spec_samp_coarse
,spat_samp
,tweak_slits
,tweak_slits_thresh
,tweak_slits_maxfrac
,rej_sticky
,slit_trim
,slit_illum_pad
,illum_iter
,illum_rej
, andtwod_fit_npoly
,saturated_slits
.Revision History:
11-Mar-2005 First version written by Scott Burles.
2005-2018 Improved by J. F. Hennawi and J. X. Prochaska
3-Sep-2018 Ported to python by J. F. Hennawi and significantly improved
- Parameters:
spat_illum_only (
bool
, optional) – If true, only the spatial illumination profile will be calculated. The 2D bspline fit will not be performed. This is primarily used to build an illumflat.doqa (
bool
, optional) – Save the QA?debug (
bool
, optional) – Show plots useful for debugging. This will block further execution of the code until the plot windows are closed.
- property nslits
Return the number of slits. Pulled directly from
slits
, if it exists.
- run(doqa=False, debug=False, show=False)[source]
Generate normalized pixel and illumination flats.
This is a simple wrapper for the main flat-field methods:
The method is a simple wrapper for
fit()
andshow()
.- Parameters:
- Returns:
Container with the results of the flat-field analysis.
- Return type:
- show(wcs_match=True)[source]
Show all of the flat field products in ginga.
- Parameters:
wcs_match (
bool
, optional) – Match the WCS of the flat-field images
- spatial_fit(norm_spec, spat_coo, median_slit_width, spat_gpm, gpm, debug=False)[source]
Perform the spatial fit
- Parameters:
norm_spec (numpy.ndarray)
spat_coo (numpy.ndarray) – Spatial coordinate array
median_slit_width (
float
)spat_gpm (numpy.ndarray)
gpm (numpy.ndarray)
debug (bool, optional)
- Returns:
- 7 objects
exit_status (int):
spat_coo_data
spat_flat_data
spat_bspl (
bspline
): Bspline model of the spatial fit. Used for illumflatspat_gpm_fit
spat_flat_fit
spat_flat_data_raw
- Return type:
- spatial_fit_finecorr(normed, onslit_tweak, slit_idx, slit_spat, gpm, slit_trim=3, tolerance=0.1, doqa=False)[source]
Generate a relative scaling image for a slicer IFU. All slits are scaled relative to a reference slit, specified in the spectrograph settings file.
- Parameters:
normed (numpy.ndarray) – Raw flat field image, normalized by the spectral and spatial illuminations.
onslit_tweak (numpy.ndarray) – mask indicating which pixels are on the slit (True = on slit)
slit_idx (int) – Slit number (0-indexed)
slit_spat (int) – Spatial ID of the slit
gpm (numpy.ndarray) – Good pixel mask
slit_txt (str) – if pypeline is “Echelle”, then slit_txt should be set to “order”, otherwise, use “slit”
slit_trim (int, optional) – Trim the slit edges by this number of pixels during the fitting. Note that the fit will be evaluated on the pixels indicated by onslit_tweak. A positive number trims the slit edges, a negative number pads the slit edges.
tolerance (float, optional) – Tolerance for the relative scaling of the slits. A value of 0.1 means that the relative scaling of the slits must be within 10% of unity. Any data outside of this tolerance will be masked.
doqa (
bool
, optional:) – Save the QA?
- Returns:
illumflat_finecorr – An image (same shape as normed) containing the fine correction to the spatial illumination profile
- Return type:
- spectral_illumination(gpm=None, debug=False)[source]
Generate a relative scaling image for a slicer IFU. All slits are scaled relative to a reference slit, specified in the spectrograph settings file.
- Parameters:
gpm (numpy.ndarray, None) – Good pixel mask
debug (bool) – Debug the routine
- Returns:
scale_model – An image containing the appropriate scaling
- Return type:
- tweak_slit_edges(left, right, spat_coo, norm_flat, method='threshold', thresh=0.93, maxfrac=0.1, debug=False)[source]
Tweak the slit edges based on the normalized slit illumination profile.
- Parameters:
left (numpy.ndarray) – Array with the left slit edge for a single slit. Shape is \((N_{\rm spec},)\).
right (numpy.ndarray) – Array with the right slit edge for a single slit. Shape is \((N_{\rm spec},)\).
spat_coo (numpy.ndarray) – Spatial pixel coordinates in fractions of the slit width at each spectral row for the provided normalized flat data. Coordinates are relative to the left edge (with the left edge at 0.). Shape is \((N_{\rm flat},)\). Function assumes the coordinate array is sorted.
norm_flat (numpy.ndarray) – Normalized flat data that provide the slit illumination profile. Shape is \((N_{\rm flat},)\).
method (
str
, optional) –Method to use for tweaking the slit edges. Options are:
'threshold'
: Use the threshold to set the slit edge and then shift it to the left or right based on the illumination profile.'gradient'
: Use the gradient of the illumination profile to set the slit edge and then shift it to the left or right based on the illumination profile.
thresh (
float
, optional) – Threshold of the normalized flat profile at which to place the two slit edges.maxfrac (
float
, optional) – The maximum fraction of the slit width that the slit edge can be adjusted by this algorithm. Ifmaxfrac = 0.1
, this means the maximum change in the slit width (either narrowing or broadening) is 20% (i.e., 10% for either edge).debug (
bool
, optional) – Show flow interrupting plots that show illumination profile in the case of a failure and the placement of the tweaked edge for each side of the slit regardless.
- Returns:
Returns six objects:
The threshold used to set the left edge
The fraction of the slit that the left edge is shifted to the right
The adjusted left edge
The threshold used to set the right edge
The fraction of the slit that the right edge is shifted to the left
The adjusted right edge
- Return type:
- class pypeit.flatfield.FlatImages(pixelflat_raw=None, pixelflat_norm=None, pixelflat_bpm=None, pixelflat_model=None, pixelflat_spat_bsplines=None, pixelflat_finecorr=None, pixelflat_spec_illum=None, pixelflat_waveimg=None, illumflat_raw=None, illumflat_spat_bsplines=None, illumflat_bpm=None, illumflat_finecorr=None, PYP_SPEC=None, spat_id=None)[source]
Bases:
CalibFrame
Container for the processed flat-field calibrations.
All of the items in the datamodel are required for instantiation, although they can be None (but shouldn’t be).
The datamodel attributes are:
Version: 1.1.2
Attribute
Type
Array Type
Description
PYP_SPEC
str
PypeIt spectrograph name
illumflat_bpm
Mirrors SlitTraceSet mask for flat-specific flags
illumflat_finecorr
PypeIt 2D polynomial fits to the fine correction of the spatial illumination profile
illumflat_raw
Processed, combined illum flats
illumflat_spat_bsplines
B-spline models for illum flat; see
bspline
pixelflat_bpm
Mirrors SlitTraceSet mask for flat-specific flags
pixelflat_finecorr
PypeIt 2D polynomial fits to the fine correction of the spatial illumination profile
pixelflat_model
Model flat
pixelflat_norm
Normalized pixel flat
pixelflat_raw
Processed, combined pixel flats
pixelflat_spat_bsplines
B-spline models for pixel flat; see
bspline
pixelflat_spec_illum
Relative spectral illumination
pixelflat_waveimg
Waveimage for pixel flat
spat_id
Slit spat_id
- _bundle()[source]
Over-write default _bundle() method to write one HDU per image. Any extras are in the HDU header of the primary image.
- Returns:
A list of dictionaries, each list element is written to its own fits extension. See the description above.
- Return type:
- classmethod _parse(hdu, ext=None, transpose_table_arrays=False, hdu_prefix=None, **kwargs)[source]
Override base-class function to deal with the many idiosyncracies of the datamodel.
See
_parse()
for the argument descriptions, and returned items.
- calib_file_format = 'fits'
The extension and file format of the output file. Should be
'fits'
or'fits.gz'
(for gzipped output).
- calib_type = 'Flat'
The type of the calibration frame, primarily used to set the name of the output file.
- datamodel = {'PYP_SPEC': {'descr': 'PypeIt spectrograph name', 'otype': <class 'str'>}, 'illumflat_bpm': {'atype': <class 'numpy.integer'>, 'descr': 'Mirrors SlitTraceSet mask for flat-specific flags', 'otype': <class 'numpy.ndarray'>}, 'illumflat_finecorr': {'atype': <class 'pypeit.core.fitting.PypeItFit'>, 'descr': 'PypeIt 2D polynomial fits to the fine correction of the spatial illumination profile', 'otype': <class 'numpy.ndarray'>}, 'illumflat_raw': {'atype': <class 'numpy.floating'>, 'descr': 'Processed, combined illum flats', 'otype': <class 'numpy.ndarray'>}, 'illumflat_spat_bsplines': {'atype': <class 'pypeit.bspline.bspline.bspline'>, 'descr': 'B-spline models for illum flat; see :class:`~pypeit.bspline.bspline.bspline`', 'otype': <class 'numpy.ndarray'>}, 'pixelflat_bpm': {'atype': <class 'numpy.integer'>, 'descr': 'Mirrors SlitTraceSet mask for flat-specific flags', 'otype': <class 'numpy.ndarray'>}, 'pixelflat_finecorr': {'atype': <class 'pypeit.core.fitting.PypeItFit'>, 'descr': 'PypeIt 2D polynomial fits to the fine correction of the spatial illumination profile', 'otype': <class 'numpy.ndarray'>}, 'pixelflat_model': {'atype': <class 'numpy.floating'>, 'descr': 'Model flat', 'otype': <class 'numpy.ndarray'>}, 'pixelflat_norm': {'atype': <class 'numpy.floating'>, 'descr': 'Normalized pixel flat', 'otype': <class 'numpy.ndarray'>}, 'pixelflat_raw': {'atype': <class 'numpy.floating'>, 'descr': 'Processed, combined pixel flats', 'otype': <class 'numpy.ndarray'>}, 'pixelflat_spat_bsplines': {'atype': <class 'pypeit.bspline.bspline.bspline'>, 'descr': 'B-spline models for pixel flat; see :class:`~pypeit.bspline.bspline.bspline`', 'otype': <class 'numpy.ndarray'>}, 'pixelflat_spec_illum': {'atype': <class 'numpy.floating'>, 'descr': 'Relative spectral illumination', 'otype': <class 'numpy.ndarray'>}, 'pixelflat_waveimg': {'atype': <class 'numpy.floating'>, 'descr': 'Waveimage for pixel flat', 'otype': <class 'numpy.ndarray'>}, 'spat_id': {'atype': <class 'numpy.integer'>, 'descr': 'Slit spat_id', 'otype': <class 'numpy.ndarray'>}}
Default datamodel for any
CalibFrame
. Derived classes should instantiate their datamodels by first inheriting from the base class. E.g.:class ArcFrame(CalibFrame): datamodel = {**CalibFrame.datamodel, ...}
- fit2illumflat(slits, frametype='illum', finecorr=False, initial=False, spat_flexure=None)[source]
Construct the model flat using the spatial bsplines.
- Parameters:
slits (
SlitTraceSet
) – Definition of the slit edgesframetype (str, optional) – The frame type should be ‘illum’ to return the illumflat version, or ‘pixel’ to return the pixelflat version.
finecorr (bool, optional) – Return the fine correction bsplines (finecorr=True), or the zeroth order correction (finecorr=False)
initial (bool, optional) – If True, the initial slit edges will be used
spat_flexure (float, optional) – Spatial flexure in pixels
- Returns:
An image of the spatial illumination profile for all slits.
- Return type:
- get_bpmflats(frametype='pixel')[source]
Get the processed bad-pixel mask.
- Parameters:
frametype (
str
, optional) – The type of mask to return. Must be either ‘illum’ for the illumination flat mask or ‘pixel’ for the pixel flat mask.- Returns:
The selected mask. If neither the illumination flat or pixel flat mask exist, the returned array is fully unmasked (all values are False).
- Return type:
- get_procflat(frametype='pixel')[source]
Get the processed flat data.
- Parameters:
frametype (
str
, optional) – The type of flat to return. Must be either ‘illum’ for the illumination flat or ‘pixel’ for the pixel flat.- Returns:
The selected flat. Can be None if the flat has not been instantiated/processed.
- Return type:
- get_spat_bsplines(frametype='illum', finecorr=False)[source]
Grab a list of bspline fits
- Parameters:
- Returns:
The selected list of spatial bsplines. Can be None if the requested data (or the fall-back) do not exist.
- Return type:
- is_synced(slits)[source]
Confirm the slits in WaveTilts are aligned to that in SlitTraceSet
Barfs if not
- Parameters:
slits (
SlitTraceSet
)
- property shape
Shape of the image arrays.
- show(frametype='all', slits=None, wcs_match=True, chk_version=True)[source]
Simple wrapper to
show_flats()
.- Parameters:
frametype (str, optional) – String used to select the flats to be displayed. The frame type should be ‘illum’ to show the illumflat version, ‘pixel’ to show the pixelflat version, or ‘all’ to show both.
slits (
SlitTraceSet
) – Definition of the slit edgeswcs_match (
bool
, optional) – (Attempt to) Match the WCS coordinates of the output images in the ginga viewer.chk_version (
bool
, optional) – When reading in existing files written by PypeIt, perform strict version checking to ensure a valid file. If False, the code will try to keep going, but this may lead to faults and quiet failures. User beware!
- version = '1.1.2'
Provides the string representation of the class version.
This is currently put to minimal use so far, but will used for I/O verification in the future.
Each derived class should provide a version to guard against data model changes during development.
- class pypeit.flatfield.SlitlessFlat(fitstbl, slitless_rows, spectrograph, par, qa_path=None)[source]
Bases:
object
Class to generate a slitless pixel flat-field calibration image.
- Parameters:
fitstbl (
PypeItMetaData
) – The class holding the metadata for all the frames.slitless_rows (numpy.ndarray) – Boolean array selecting the rows in the fitstbl that correspond to the slitless frames.
spectrograph (
Spectrograph
) – The spectrograph object.par (
CalibrationsPar
) – Parameter set defining optional parameters of PypeIt’s algorithms for Calibrationsqa_path (Path) – Path for the QA diagnostics.
- make_slitless_pixflat(msbias=None, msdark=None, calib_dir=None, write_qa=False, show=False)[source]
Generate and save to disc a slitless pixel flat-field calibration images. The pixel flat file will have one extension per detector, even in the case of a mosaic. Contrary to the regular calibration flow, the slitless pixel flat is created for all detectors of the current spectrograph at once, and not only the one for the current detector. Since the slitless pixel flat images are saved to disc, this approach helps with the I/O This is a method is used in ~pypeit.calibrations.get_flats().
Note: par[‘flatfield’][‘pixelflat_file’] is updated in this method.
- Parameters:
msbias (
BiasImage
, optional) – Bias image for bias subtraction; passed tobuildimage_fromlist()
msdark (
DarkImage
, optional) – Dark-current image; passed tobuildimage_fromlist()
calib_dir (Path) – Path for the processed calibration files.
write_qa (
bool
, optional) – Write QA plots to disk?show (
bool
, optional) – Show the diagnostic plots?
- Returns:
The name of the slitless pixel flat file that was generated.
- Return type:
- pypeit.flatfield.detector_structure_qa(det_resp, det_resp_model, outfile=None, title='Detector Structure Correction')[source]
Plot the QA for the fine correction fits to the spatial illumination profile
- Parameters:
det_resp (numpy.ndarray) – Image data showing the detector structure, generated with extract_structure
det_resp_model (numpy.ndarray) – Image containing the structure correction model
outfile (str, optional) – Output file name
title (str, optional) – A title to be printed on the QA
- pypeit.flatfield.illum_profile_spectral(rawimg, waveimg, slits, slit_illum_ref_idx=0, smooth_npix=None, polydeg=None, model=None, gpmask=None, skymask=None, trim=3, flexure=None, maxiter=5, debug=False)[source]
Determine the relative spectral illumination of all slits. Currently only used for image slicer IFUs.
- Parameters:
rawimg (numpy.ndarray) – Image data that will be used to estimate the spectral relative sensitivity
waveimg (numpy.ndarray) – Wavelength image
slits (
SlitTraceSet
) – Information stored about the slitsslit_illum_ref_idx (int) – Index of slit that is used as the reference.
smooth_npix (int, optional) – smoothing used for determining smoothly varying relative weights by sn_weights
polydeg (int, optional) – Degree of polynomial to be used for determining relative spectral sensitivity. If None, coadd.smooth_weights will be used, with the smoothing length set to smooth_npix.
model (numpy.ndarray, None) – A model of the rawimg data. If None, rawimg will be used.
gpmask (numpy.ndarray, None) – Good pixel mask
skymask (numpy.ndarray, None) – Sky mask
trim (int) – Number of pixels to trim from the edges of the slit when deriving the spectral illumination
flexure (float, None) – Spatial flexure
maxiter (
int
) – Maximum number of iterations to performdebug (
bool
) – Show the results of the relative spectral illumination correction
- Returns:
scale_model – An image containing the appropriate scaling
- Return type:
- pypeit.flatfield.load_pixflat(pixel_flat_file, spectrograph, det, flatimages, calib_dir=None, chk_version=False)[source]
Load a pixel flat from a file and add it to the flatimages object. The pixel flat file has one detector per extension, even in the case of a mosaic. Therefore, if this is a mosaic reduction, this script will construct a pixel flat mosaic. The Edges file needs to exist in the Calibration Folder, since the mosaic parameters are pulled from it. This is used in ~pypeit.calibrations.get_flats().
- Parameters:
pixel_flat_file (
str
) – Name of the pixel flat file.spectrograph (
Spectrograph
) – The spectrograph object.det (
int
,tuple
) – The single detector or set of detectors in a mosaic to process.flatimages (
FlatImages
) – The flat field images object.calib_dir (
str
, optional) – The path to the calibration directory.chk_version (
bool
, optional) – Check the version of the file.
- Returns:
The flat images object with the pixel flat added.
- Return type:
- pypeit.flatfield.merge(init_cls, merge_cls)[source]
Merge
merge_cls
intoinit_cls
, and return a mergedFlatImages
class.If
init_cls
is None, the returned value ismerge_cls
, and vice versa. If an element exists in both init_cls and merge_cls, the merge_cls value is taken- Parameters:
init_cls (
FlatImages
) – Initial class (the elements of this class will be considered the default). Can be None.merge_cls (
FlatImages
) – The non-zero elements will be merged into init_cls. Can be None.
- Returns:
flats – A new instance of the FlatImages class with merged properties.
- Return type:
- pypeit.flatfield.show_flats(image_list, wcs_match=True, slits=None, waveimg=None)[source]
Show the flat-field images
- Parameters:
image_list (list) – List of tuples containing the image data, image name and the cut levels
wcs_match (bool, optional) – Match the WCS of the images
slits (
pypeit.slittrace.SlitTraceSet
, optional) – Slit traces to be overplotted on the imageswaveimg (numpy.ndarray, optional) – Wavelength image to be overplotted on the images
- pypeit.flatfield.spatillum_finecorr_qa(normed, finecorr, left, right, ypos, cut, outfile=None, title=None, half_slen=50)[source]
Plot the QA for the fine correction fits to the spatial illumination profile
- Parameters:
normed (numpy.ndarray) – Image data with the coarse spatial illumination profile divided out (normalised in the spectral direction)
finecorr (numpy.ndarray) – Image containing the fine correction
left (numpy.ndarray) – Left slit edge
right (numpy.ndarray) – Right slit edge
ypos (numpy.ndarray) – Spectral coordinate (from 0 to 1, where 0=blue wavelength, 1=red wavelength)
cut (tuple) – A 2-tuple, containing the (x, y) coordinates of the pixels on the slit.
outfile (str, optional) – Output file name
title (str, optional) – A title to be printed on the QA
half_slen (int, optional) – The sampling size for the spatial profile. Should be about half the slit length. In this case, each output pixel shown contains about 2 detector pixels.
- pypeit.flatfield.write_pixflat_to_fits(pixflat_norm_list, detname_list, spec_name, outdir, pixelflat_name, to_cache=True)[source]
Write the pixel-to-pixel flat-field images to a FITS file. The FITS file will have an extension for each detector (never a mosaic). The load_pixflat() method read this file and transform it into a mosaic if needed. This image is generally used as a user-provided pixel flat-field image and ingested in the reduction using the pixelflat_file parameter in the PypeIt file.
- Parameters:
pixflat_norm_list (
list
) – List of 2D numpy.ndarray arrays containing the pixel-to-pixel flat-field images.detname_list (
list
) – List of detector names.spec_name (
str
) – Name of the spectrograph.outdir (
pathlib.Path
) – Path to the output directory.pixelflat_name (
str
) – Name of the output file to be written.to_cache (
bool
, optional) – If True, the file will be written to the cache directory pypeit/data/pixflats.