pypeit.core.flexure module

Module for flexure routines

class pypeit.core.flexure.MultiSlitFlexure(s1dfile=None, PYP_SPEC=None, nslits=None, det=None, SN=None, slitid=None, mn_wv=None, fit_slope=None, fit_b=None, fit_los=None, objra=None, objdec=None, maskdef_id=None, rms_arc=None, resid_sky=None, indiv_fit_slope=None, indiv_fit_b=None, indiv_fit_los=None)[source]

Bases: DataContainer

Class to perform multi-detector flexure analysis.

Based on code written by Marla Geha for DEIMOS.

The datamodel attributes are:

Version: 1.1.0

Attribute

Type

Array Type

Description

PYP_SPEC

str

PypeIt spectrograph name

SN

numpy.ndarray

numpy.floating

S/N (ndet, nslits)

det

numpy.ndarray

int, numpy.integer

Integer identifiers for the detector or mosaic (ndet, nslits)

fit_b

numpy.ndarray

numpy.floating

Fitted b value(nslits)

fit_los

numpy.ndarray

numpy.floating

Fitted line width(nslits)

fit_slope

numpy.ndarray

numpy.floating

Fitted slope (nslits)

indiv_fit_b

numpy.ndarray

numpy.floating

Same as above but for b (nslits)

indiv_fit_los

numpy.ndarray

numpy.floating

Same as above but for line width (nslits)

indiv_fit_slope

numpy.ndarray

numpy.floating

Fits to each slit individually (nslits)

is_msc

numpy.ndarray

int, numpy.integer

Flag that the “det” is the mosaic ID (ndet, nslits)

maskdef_id

numpy.ndarray

numpy.integer

Mask ID (nslits)

mn_wv

numpy.ndarray

numpy.floating

Mininum wavelength of the slit [Ang] (nslits)

ndet

int

Number of detectors per spectrum

nslits

int

Number of slits

objdec

numpy.ndarray

numpy.floating

Object DEC (nslits)

objra

numpy.ndarray

numpy.floating

Object RA (nslits)

resid_sky

numpy.ndarray

numpy.floating

Residuals of flexure model on sky lines (nslits)

rms_arc

numpy.ndarray

numpy.floating

RMS of fit (ndet, nslits)

s1dfile

str

spec1d filename

slitid

numpy.ndarray

numpy.floating

Slit ID (nslits)

_bundle()[source]

Override the base class method simply to set the HDU extension name.

datamodel = {'PYP_SPEC': {'descr': 'PypeIt spectrograph name', 'otype': <class 'str'>}, 'SN': {'atype': <class 'numpy.floating'>, 'descr': 'S/N (ndet, nslits)', 'otype': <class 'numpy.ndarray'>}, 'det': {'atype': (<class 'int'>, <class 'numpy.integer'>), 'descr': 'Integer identifiers for the detector or mosaic (ndet, nslits)', 'otype': <class 'numpy.ndarray'>}, 'fit_b': {'atype': <class 'numpy.floating'>, 'descr': 'Fitted b value(nslits)', 'otype': <class 'numpy.ndarray'>}, 'fit_los': {'atype': <class 'numpy.floating'>, 'descr': 'Fitted line width(nslits)', 'otype': <class 'numpy.ndarray'>}, 'fit_slope': {'atype': <class 'numpy.floating'>, 'descr': 'Fitted slope (nslits)', 'otype': <class 'numpy.ndarray'>}, 'indiv_fit_b': {'atype': <class 'numpy.floating'>, 'descr': 'Same as above but for b (nslits)', 'otype': <class 'numpy.ndarray'>}, 'indiv_fit_los': {'atype': <class 'numpy.floating'>, 'descr': 'Same as above but for line width (nslits)', 'otype': <class 'numpy.ndarray'>}, 'indiv_fit_slope': {'atype': <class 'numpy.floating'>, 'descr': 'Fits to each slit individually (nslits)', 'otype': <class 'numpy.ndarray'>}, 'is_msc': {'atype': (<class 'int'>, <class 'numpy.integer'>), 'descr': 'Flag that the "det" is the mosaic ID (ndet, nslits)', 'otype': <class 'numpy.ndarray'>}, 'maskdef_id': {'atype': <class 'numpy.integer'>, 'descr': 'Mask ID (nslits)', 'otype': <class 'numpy.ndarray'>}, 'mn_wv': {'atype': <class 'numpy.floating'>, 'descr': 'Mininum wavelength of the slit [Ang] (nslits)', 'otype': <class 'numpy.ndarray'>}, 'ndet': {'descr': 'Number of detectors per spectrum', 'otype': <class 'int'>}, 'nslits': {'descr': 'Number of slits', 'otype': <class 'int'>}, 'objdec': {'atype': <class 'numpy.floating'>, 'descr': 'Object DEC (nslits)', 'otype': <class 'numpy.ndarray'>}, 'objra': {'atype': <class 'numpy.floating'>, 'descr': 'Object RA (nslits)', 'otype': <class 'numpy.ndarray'>}, 'resid_sky': {'atype': <class 'numpy.floating'>, 'descr': 'Residuals of flexure model on sky lines (nslits)', 'otype': <class 'numpy.ndarray'>}, 'rms_arc': {'atype': <class 'numpy.floating'>, 'descr': 'RMS of fit (ndet, nslits)', 'otype': <class 'numpy.ndarray'>}, 's1dfile': {'descr': 'spec1d filename', 'otype': <class 'str'>}, 'slitid': {'atype': <class 'numpy.floating'>, 'descr': 'Slit ID (nslits)', 'otype': <class 'numpy.ndarray'>}}

Provides the class data model. This is undefined in the abstract class and should be overwritten in the derived classes.

The format of the datamodel needed for each implementation of a DataContainer derived class is as follows.

The datamodel itself is a class attribute (i.e., it is a member of the class, not just of an instance of the class). The datamodel is a dictionary of dictionaries: Each key of the datamodel dictionary provides the name of a given datamodel element, and the associated item (dictionary) for the datamodel element provides the type and description information for that datamodel element. For each datamodel element, the dictionary item must provide:

  • otype: This is the type of the object for this datamodel item. E.g., for a float or a numpy.ndarray, you would set otype=float and otype=np.ndarray, respectively.

  • descr: This provides a text description of the datamodel element. This is used to construct the datamodel tables in the pypeit documentation.

If the object type is a numpy.ndarray, you should also provide the atype keyword that sets the type of the data contained within the array. E.g., for a floating point array containing an image, your datamodel could be simply:

datamodel = {'image' : dict(otype=np.ndarray, atype=float, descr='My image')}

More advanced examples are given in the top-level module documentation.

Currently, datamodel components are restricted to have otype that are tuple, int, float, numpy.integer, numpy.floating, numpy.ndarray, or astropy.table.Table objects. E.g., datamodel values for otype cannot be dict.

fit_mask_surfaces()[source]

Fit 2D model to linear flexure models from each slit as a function of RA, DEC.

init(spectrograph, par)[source]

Initialize this and that about the slits, par, spectrograph e.g. RA, DEC, S/N

Parameters:
internals = ['flex_par', 'spectrograph', 'specobjs', 'sobj_idx', 'sky_table', 'pmodel_m', 'pmodel_b', 'pmodel_l']

A list of strings identifying a set of internal attributes that are not part of the datamodel.

measure_sky_lines()[source]

Main method to analyze the sky lines for all the slits

qa_plots(plot_dir: str, root: str)[source]

Generate QA plots

Parameters:
  • plot_dir (str) – Top-lvel folder for QA QA/ is generated beneath this, as needed

  • root (str) – Root for output files

update_fit()[source]

Update fits for each slit based on 2D model

version = '1.1.0'

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.

pypeit.core.flexure.calculate_image_offset(im_ref, image, nfit=3)[source]

Calculate the x,y offset between two images

Parameters:
  • im_ref (numpy.ndarray) – Reference image

  • image (numpy.ndarray) – Image that we want to measure the shift of (relative to im_ref)

  • nfit (int, optional) – Number of pixels (left and right of the maximum) to include in fitting the peak of the cross correlation.

Returns:

Returns two floats, the x and y offset of the image.
  • ra_diff – Relative shift (in pixels) of image relative to im_ref (x direction). In order to align image with im_ref, ra_diff should be added to the x-coordinates of image

  • dec_diff – Relative shift (in pixels) of image relative to im_ref (y direction). In order to align image with im_ref, dec_diff should be added to the y-coordinates of image

Return type:

tuple

pypeit.core.flexure.calculate_image_phase(imref, imshift, gpm_ref=None, gpm_shift=None, maskval=None)[source]

Perform a masked cross-correlation and optical flow calculation to robustly estimate the subpixel shifts of two images.

If gpm_ref, gpm_shift, and maskval are all None, no pixels will be masked

This routine (optionally) requires skimage to calculate the image phase. If skimage is not installed, a standard (unmasked) cross-correlation is used.

Parameters:
  • im_ref (numpy.ndarray) – Reference image

  • imshift (numpy.ndarray) – Image that we want to measure the shift of (relative to im_ref)

  • gpm_ref (numpy.ndarray) – Mask of good pixels (True = good) in the reference image

  • gpm_shift (numpy.ndarray) – Mask of good pixels (True = good) in the shifted image

  • maskval (float, optional) – If gpm_ref and gpm_shift are both None, a single value can be specified and this value will be masked in both images.

Returns:

  • ra_diff (float) – Relative shift (in pixels) of image relative to im_ref (x direction). In order to align image with im_ref, ra_diff should be added to the x-coordinates of image

  • dec_diff (float) – Relative shift (in pixels) of image relative to im_ref (y direction). In order to align image with im_ref, dec_diff should be added to the y-coordinates of image

pypeit.core.flexure.flexure_interp(shift, wave)[source]

Perform interpolation on wave given a shift in pixels

Parameters:
  • shift (float) – Shift in pixels

  • wave (numpy.ndarray) – extracted wave of size nspec

  • wavein (numpy.ndarray, optional) – Apply the shift to this array of wavelengths

Returns:

Wavelength scale corrected for spectral flexure

Return type:

numpy.ndarray

pypeit.core.flexure.get_archive_spectrum(sky_file, obj_skyspec=None, spec_fwhm_pix=None)[source]

Load an archival sky spectrum

Parameters:
  • sky_file (str) – Name of the archival sky file. If equal to ‘model’, instead, a model sky spectrum will be generated using nearIR_modelsky() and the spectral resolution of obj_skyspec. If obj_skyspec is None, then sky_file cannot be ‘model’.

  • obj_skyspec (linetools.spectra.xspectrum1d.XSpectrum1d, optional) – Sky spectrum associated with the science target. This must be provided if sky_file is ‘model’.

  • spec_fwhm_pix (float, optional) – Spectral FWHM (in pixels) of the sky spectrum related to our object.

Returns:

The sky spectrum (linetools.spectra.xspectrum1d.XSpectrum1D) and the FWHM (float) of the sky lines in pixels.

Return type:

tuple

pypeit.core.flexure.get_fwhm_gauss_smooth(arx_skyspec, obj_skyspec, arx_fwhm_pix, spec_fwhm_pix=None)[source]
Parameters:
Returns:

FWHM of the smoothing Gaussian in pixels.

Return type:

float

pypeit.core.flexure.get_percentile_clipping(arr, percent=90.0)[source]

Get the values for clipping based on a percentile

Parameters:
  • arr (numpy.ndarray) – Array to clip.

  • percent (float) – Percentile to clip at. Default is 90.0

Returns:

Lower value for clipping float: Upper value for clipping

Return type:

float

pypeit.core.flexure.get_sky_spectrum(sciimg, ivar, waveimg, thismask, global_sky, box_radius, slits, trace_spat, pypeline, det)[source]

Obtain a boxcar extraction of the sky spectrum

Parameters:
  • sciimg (numpy.ndarray) – Science image - shape (nspec, nspat)

  • ivar (numpy.ndarray) – Inverse variance of the science image - shape (nspec, nspat)

  • waveimg (numpy.ndarray) – Wavelength image - shape (nspec, nspat)

  • thismask (numpy.ndarray) – Good pixel mask (True=good) that indicates the pixels that should be included in the boxcar extraction

  • global_sky (numpy.ndarray) – 2D array of the global_sky fit - shape (nspec, nspat)

  • box_radius (float) – Radius of the boxcar extraction (in pixels)

  • slits (SlitTraceSet) – Slit trace set

  • trace_spat (numpy.ndarray) – Spatial pixel values (usually the center of each slit) where the sky spectrum will be extracted. The shape of this array should be (nspec, nslits)

  • pypeline (str) – Name of the PypeIt pipeline method. Allowed options are MultiSlit, Echelle, or IFU.

  • det (str) – The name of the detector or mosaic from which the spectrum will be extracted. For example, DET01.

Returns:

Sky spectrum

Return type:

(linetools.spectra.xspectrum1d.XSpectrum1D)

pypeit.core.flexure.sky_em_residuals(wave: ndarray, flux: ndarray, ivar: ndarray, sky_waves: ndarray, plot=False, noff=5.0, nfit_min=20)[source]

Calculate residuals and other metrics for a set of input sky emission lines

Parameters:
  • wave (numpy.ndarray) – Wavelengths (in air!)

  • flux (numpy.ndarray) – Fluxes

  • ivar (numpy.ndarray) – Inverse variance

  • sky_waves (numpy.ndarray) – Skyline wavelengths (in air!)

  • plot (bool, optional) – If true, plot the residuals

  • noff (int, optional) – Range in Ang to analyze labout emission line. Defaults to 5.

  • nfit_min (int, optional) – Minimum number of pixels required to do a fit. Defaults to 20.

Returns:

tuple of numpy.ndarray – sky line wavelength of good lines, wavelength offset,

error in wavelength offset, sky line width, error in sky line width

pypeit.core.flexure.spat_flexure_shift(sciimg, slits, debug=False, maxlag=20)[source]

Calculate a rigid flexure shift in the spatial dimension between the slitmask and the science image.

It is important to use original=True when defining the slitmask as everything should be relative to the initial slits

Otherwise, the WaveTilts could get out of sync with science images

Parameters:
Returns:

The spatial flexure shift relative to the initial slits

Return type:

float

pypeit.core.flexure.spec_flex_shift(obj_skyspec, sky_file=None, arx_skyspec=None, arx_fwhm_pix=None, spec_fwhm_pix=None, mxshft=20, excess_shft='crash', method='boxcar', minwave=None, maxwave=None)[source]

Calculate shift between object sky spectrum and archive sky spectrum

Parameters:
  • obj_skyspec (linetools.spectra.xspectrum1d.XSpectrum1d) – Spectrum of the sky related to our object

  • sky_file (str, optional) – Name of the archival sky file. If equal to ‘model’, instead, a model sky spectrum will be generated using nearIR_modelsky() and the spectral resolution of obj_skyspec. If None, arx_skyspec and arx_fwhm_pix must be provided.

  • arx_skyspec (linetools.spectra.xspectrum1d.XSpectrum1d, optional) – Archived sky spectrum. If None, it will be loaded from the sky_file (sky_file must be provided).

  • arx_fwhm_pix (float, optional) – Spectral FWHM (in pixels) of the archived sky spectrum. If None, it will be calculated using sky_file (sky_file must be provided).

  • spec_fwhm_pix (float, optional) – Spectral FWHM (in pixels) of the sky spectrum related to our object/slit.

  • mxshft (int, optional) – Maximum allowed shift from flexure; note there are cases that have been known to exceed even 30 pixels.

  • excess_shft (str, optional) – Behavior of the code when a measured flexure exceeds mxshft. Options are “crash”, “set_to_zero”, and “continue”, where “set_to_zero” sets the shift to zero and moves on, and “continue” simply uses the large flexure shift value.

  • method (str, optional) – Which method is used for the spectral flexure correction. Two methods are available: ‘boxcar’ and ‘slitcen’ (see spec_flexure_slit()). In this routine, ‘method’ is only passed to final dict.

  • minwave (float, optional) – Minimum wavelength to use for the correlation. If None or less than the minumum wavelength of either obj_skyspec or arx_skyspec, this has no effect. Default is None.

  • maxwave (float, optional) – Maximum wavelength to use for the correlation. If None or greater than the maximum wavelength of either obj_skyspec or arx_skyspec, this has no effect. Default is None.

Returns:

Contains flexure info. Keys are:

  • polyfit= fit to the cross-correlation

  • shift= best shift in pixels

  • subpix= subpixelation of input spectrum

  • corr= correlation function

  • sky_spec= object sky spectrum used (rebinned, etc.)

  • arx_spec= archived sky spectrum used

  • corr_cen= center of the correlation function

  • smooth= Degree of smoothing of input spectrum to match archive

Return type:

dict

pypeit.core.flexure.spec_flex_shift_global(slit_specs, islit, sky_file, empty_flex_dict, return_later_slits, flex_list, keys_to_update, spec_fwhm_pix=None, mxshft=20, excess_shft='crash', method='slitcen', minwave=None, maxwave=None)[source]

Calculate flexure shifts using the sky spectrum extracted at the center of the slit

Parameters:
  • slit_specs (list) – A list of linetools.xspectrum1d, one for each slit. The spectra stored in this list are sky spectra, extracted from the center of each slit.

  • islit (int) – Index of the slit where the sky spectrum related to our object is.

  • sky_file (str) – Name of the archival sky file. If equal to ‘model’, instead, a model sky spectrum will be generated using nearIR_modelsky() and the spectral resolution of each spectrum from slit_specs.

  • empty_flex_dict (dict) – Empty dictionary to be filled with flexure results.

  • return_later_slits (list) – List of slit indexes that failed the shift calcultion and we want to come back to to assign a value from a different slit.

  • flex_list (list) – A list of dict objects containing flexure results of each slit.

  • keys_to_update (list) – List of flexure dictionary keys that we need to update.

  • spec_fwhm_pix (float, optional) – Spectral FWHM (in pixels) of the sky spectrum related to our object.

  • mxshft (int, optional) – Maximum allowed shift from flexure. Passed to spec_flex_shift().

  • excess_shft (str, optional) – Behavior of the code when a measured flexure exceeds mxshft. Passed to spec_flex_shift()

  • method (str, optional) – Which method is used for the spectral flexure correction. Two methods are available: ‘boxcar’ and ‘slitcen’ (see spec_flexure_slit()). Passed to spec_flex_shift().

  • minwave (float, optional) – Minimum wavelength to use for the correlation. If None or less than the minumum wavelength of either this sky or sky_spectrum, this has no effect. Default is None.

  • maxwave (float, optional) – Maximum wavelength to use for the correlation. If None or greater than the maximum wavelength of either this sky or sky_spectrum, this has no effect. Default is None.

Returns:

A list of dict objects containing flexure results of each slit. This is filled with a basically empty dict if the shift calculation failed for the relevant slit.

Return type:

list

pypeit.core.flexure.spec_flex_shift_local(slits, slitord, specobjs, islit, sky_file, empty_flex_dict, return_later_slits, flex_list, keys_to_update, spec_fwhm_pix=None, mxshft=20, excess_shft='crash', method='boxcar', minwave=None, maxwave=None)[source]

Calculate flexure shifts using the sky spectrum boxcar-extracted at the location of the detected objects

Parameters:
  • slits (SlitTraceSet) – Slit trace set.

  • slitord (numpy.ndarray) – Array of slit/order numbers.

  • specobjs (SpecObjs, optional) – Spectral extractions.

  • islit (int) – Index of the slit where the sky spectrum related to our object is.

  • sky_file (str) – Name of the archival sky file. If equal to ‘model’, instead, a model sky spectrum will be generated using nearIR_modelsky() and the spectral resolution of each spectrum in specobjs.

  • empty_flex_dict (dict) – Empty dictionary to be filled with flexure results.

  • return_later_slits (list) – List of slit indexes that failed the shift calcultion and we want to come back to to assign a value from a different slit.

  • flex_list (list) – A list of dict objects containing flexure results of each slit.

  • keys_to_update (list) – List of flexure dictionary keys that we need to update.

  • spec_fwhm_pix (float, optional) – Spectral FWHM (in pixels) of the sky spectrum related to our object.

  • mxshft (int, optional) – Maximum allowed shift from flexure. Passed to spec_flex_shift().

  • excess_shft (str, optional) – Behavior of the code when a measured flexure exceeds mxshft. Passed to spec_flex_shift()

  • method (str, optional) – Which method is used for the spectral flexure correction. Two methods are available: ‘boxcar’ and ‘slitcen’ (see spec_flexure_slit()). Passed to spec_flex_shift().

  • minwave (float, optional) – Minimum wavelength to use for the correlation. If None or less than the minumum wavelength of either this sky or sky_spectrum, this has no effect. Default is None.

  • maxwave (float, optional) – Maximum wavelength to use for the correlation. If None or greater than the maximum wavelength of either this sky or sky_spectrum, this has no effect. Default is None.

Returns:

A list of dict objects containing flexure results of each slit. This is filled with a basically empty dict if the shift calculation failed for the relevant slit.

Return type:

list

pypeit.core.flexure.spec_flexure_corrQA(ax, this_flex_dict, cntr, name)[source]
pypeit.core.flexure.spec_flexure_qa(slitords, bpm, basename, flex_list, specobjs=None, out_dir=None)[source]

Generate QA for the spectral flexure calculation

Parameters:
  • slitords (numpy.ndarray) – Array of slit/order numbers

  • bpm (numpy.ndarray) – Boolean mask; True = masked slit

  • basename (str) – Used to generate the output file name

  • flex_list (list) – list of dict objects containing the flexure information

  • specobjs (SpecObjs, optional) – Spectrally extracted objects

  • out_dir (str, optional) – Path to the output directory for the QA plots. If None, the current is used.

pypeit.core.flexure.spec_flexure_slit(slits, slitord, slit_bpm, sky_file, method='boxcar', specobjs=None, slit_specs=None, wv_calib=None, mxshft=None, excess_shft='crash', minwave=None, maxwave=None)[source]

Calculate the spectral flexure for every slit (global) or object (local)

Parameters:
  • slits (SlitTraceSet) – Slit trace set

  • slitord (numpy.ndarray) – Array of slit/order numbers

  • slit_bpm (numpy.ndarray) – True = masked slit

  • sky_file (str) – Name of the archival sky file. If equal to ‘model’, instead, a model sky spectrum will be generated using nearIR_modelsky() and the spectral resolution of each spectrum that we want to correct for flexure.

  • method (str, optional) –

    Two methods are available:
    • ’boxcar’: Recommended for object extractions. This method uses the boxcar extracted sky and wavelength spectra from the input specobjs

    • ’slitcen’: Recommended when no objects are being extracted. This method uses a spectrum (stored in slitspecs) that is extracted from the center of each slit.

  • specobjs (SpecObjs, optional) – Spectral extractions

  • slit_specs (list, optional) – A list of linetools.xspectrum1d, one for each slit. The spectra stored in this list are sky spectra, extracted from the center of each slit.

  • wv_calib (pypeit.wavecalib.WaveCalib) – Wavelength calibration object

  • mxshft (int, optional) – Passed to spec_flex_shift()

  • excess_shft (str, optional) – Passed to spec_flex_shift()

  • minwave (float, optional) – Minimum wavelength to use for the correlation. If None or less than the minumum wavelength of either this sky or sky_spectrum, this has no effect. Default is None.

  • maxwave (float, optional) – Maximum wavelength to use for the correlation. If None or greater than the maximum wavelength of either this sky or sky_spectrum, this has no effect. Default is None.

Returns:

A list of dict objects containing flexure results of each slit. This is filled with a basically empty dict if the slit is skipped.

Return type:

list

pypeit.core.flexure.spec_flexure_slit_global(sciImg, waveimg, global_sky, par, slits, slitmask, trace_spat, gd_slits, wv_calib, pypeline, det)[source]

Calculate the spectral flexure for every slit

Parameters:
  • sciImg (PypeItImage) – Science image.

  • waveimg (numpy.ndarray) – Wavelength image - shape (nspec, nspat)

  • global_sky (numpy.ndarray) – 2D array of the global_sky fit - shape (nspec, nspat)

  • par (PypeItPar) – Parameters of the reduction.

  • slits (SlitTraceSet) – Slit trace set

  • slitmask (numpy.ndarray) – An image with the slit index identified for each pixel (returned from slittrace.slit_img).

  • trace_spat (numpy.ndarray) – Spatial pixel values (usually the center of each slit) where the sky spectrum will be extracted. The shape of this array should be (nspec, nslits)

  • gd_slits (numpy.ndarray) – True = good slit

  • wv_calib (pypeit.wavecalib.WaveCalib) – Wavelength calibration

  • pypeline (str) – Name of the PypeIt pipeline method. Allowed options are MultiSlit, Echelle, or IFU.

  • det (str) – The name of the detector or mosaic from which the spectrum will be extracted. For example, DET01.

Returns:

A list of dict objects containing flexure results of each slit. This is filled with a basically empty dict if the slit is skipped.

Return type:

list