pypeit.spectrographs.keck_deimos module

Implements DEIMOS-specific functions, including reading in slitmask design files.

class pypeit.spectrographs.keck_deimos.DEIMOSCameraDistortion[source]

Bases: object

Class to remove or apply DEIMOS camera distortion.

apply_distortion(y)[source]
remove_distortion(x)[source]
class pypeit.spectrographs.keck_deimos.DEIMOSDetectorMap[source]

Bases: DetectorMap

A map of the center coordinates and rotation of each CCD in DEIMOS.

!! PIXEL COORDINATES ARE 1-INDEXED !!

class pypeit.spectrographs.keck_deimos.DEIMOSMosaicLookUp[source]

Bases: object

Provides the geometry required to mosaic Keck DEIMOS data. Similar to GeminiGMOSMosaicLookUp

geometry = {'MSC01': {'blue_det': {'rotation': -0.197, 'shift': (-2.97, -4116.7)}, 'default_shape': (8218, 2064), 'red_det': {'rotation': 0.0, 'shift': (0.0, 0.0)}}, 'MSC02': {'blue_det': {'rotation': 0.11, 'shift': (-1.92, -4114.9)}, 'default_shape': (8218, 2064), 'red_det': {'rotation': 0, 'shift': (0.0, 0.0)}}, 'MSC03': {'blue_det': {'rotation': 0.02, 'shift': (-1.23, -4114.4)}, 'default_shape': (8218, 2064), 'red_det': {'rotation': 0, 'shift': (0.0, 0.0)}}, 'MSC04': {'blue_det': {'rotation': -0.0511, 'shift': (1.05, -4108.7)}, 'default_shape': (8218, 2064), 'red_det': {'rotation': 0.0, 'shift': (0.0, 0.0)}}}
class pypeit.spectrographs.keck_deimos.DEIMOSOpticalModel(grating)[source]

Bases: OpticalModel

Derived class for the DEIMOS optical model.

mask_coo_to_grating_input_vectors(x, y)[source]

Propagate rays from the mask plane to the grating.

Taken from xidl/DEEP2/spec2d/pro/model/pre_grating.pro

Need to override parent class to add tent mirror reflection.

reset_grating(grating)[source]

Reset the grating to the provided input instance.

class pypeit.spectrographs.keck_deimos.KeckDEIMOSSpectrograph[source]

Bases: Spectrograph

Child to handle Keck/DEIMOS specific code

static _grating_orientation(slider, ruling, tilt)[source]

Return the roll, yaw, and tilt of the provided grating.

Numbers are hardwired.

From xidl/DEEP2/spec2d/pro/omodel_params.pro

Parameters:
  • slider (int) – The slider position. Should be 2, 3, or 4. If the slider is 0, the ruling must be 0.

  • ruling (str, int) – The grating ruling number. Should be 600, 831, 900, 1200, or 'other'.

  • tilt (float) – The input grating tilt.

Returns:

The roll, yaw, and tilt of the grating used by the optical model.

Return type:

tuple

property allowed_mosaics

Return the list of allowed detector mosaics.

Returns:

List of tuples, where each tuple provides the 1-indexed detector numbers that can be combined into a mosaic and processed by PypeIt.

Return type:

list

bpm(filename, det, shape=None, msbias=None)[source]

Generate a default bad-pixel mask.

Even though they are both optional, either the precise shape for the image (shape) or an example file that can be read to get the shape (filename using get_image_shape()) must be provided.

Parameters:
  • filename (str or None) – An example file to use to get the image shape.

  • det (int) – 1-indexed detector number to use when getting the image shape from the example file.

  • shape (tuple, optional) – Processed image shape Required if filename is None Ignored if filename is not None

  • msbias (numpy.ndarray, optional) – Processed bias frame used to identify bad pixels

Returns:

An integer array with a masked value set to 1 and an unmasked value set to 0. All values are set to 0.

Return type:

numpy.ndarray

camera = 'DEIMOS'

Name of the spectrograph camera or arm. This is used by specdb, so use that naming convention

check_frame_type(ftype, fitstbl, exprng=None)[source]

Check for frames of the provided type.

Parameters:
Returns:

Boolean array with the flags selecting the exposures in fitstbl that are ftype type frames.

Return type:

numpy.ndarray

comment = 'Supported gratings: 600ZD, 830G, 900ZD, 1200B, 1200G; see :doc:`deimos`'

A brief comment or description regarding PypeIt usage with this spectrograph.

compound_meta(headarr, meta_key)[source]

Methods to generate metadata requiring interpretation of the header data, instead of simply reading the value of a header card.

Parameters:
Returns:

Metadata value read from the header(s).

Return type:

object

config_independent_frames()[source]

Define frame types that are independent of the fully defined instrument configuration.

Bias and dark frames are considered independent of a configuration, but the DATE-OBS keyword is used to assign each to the most-relevant configuration frame group. See set_configurations().

Returns:

Dictionary where the keys are the frame types that are configuration independent and the values are the metadata keywords that can be used to assign the frames to a configuration group.

Return type:

dict

config_specific_par(scifile, inp_par=None)[source]

Modify the PypeIt parameters to hard-wired values used for specific instrument configurations.

Parameters:
  • scifile (str) – File to use when determining the configuration and how to adjust the input parameters.

  • inp_par (ParSet, optional) – Parameter set used for the full run of PypeIt. If None, use default_pypeit_par().

Returns:

The PypeIt parameter set adjusted for configuration specific parameter values.

Return type:

ParSet

configuration_keys()[source]

Return the metadata keys that define a unique instrument configuration.

This list is used by PypeItMetaData to identify the unique configurations among the list of frames read for a given reduction.

Returns:

List of keywords of data pulled from file headers and used to constuct the PypeItMetaData object.

Return type:

list

property default_mosaic

Return the default detector mosaic.

For instruments with no allowed detector mosaics, this must be returned as None.

classmethod default_pypeit_par()[source]

Return the default parameters to use for this instrument.

Returns:

Parameters required by all of PypeIt methods.

Return type:

PypeItPar

get_amapbmap(filename)[source]

Select the pre-grating (amap) and post-grating (bmap) maps according to the slider.

Parameters:

filename (str) – The filename to read the slider information from the header.

Returns:

The two attributes amap and bmap, used by the DEIMOS optical model.

Return type:

tuple

get_detector_map()[source]

Return the DEIMOS detector map.

Returns:

The instance describing the detector layout for DEIMOS. The object returned is the same as detector_map. If detector_map is None when the method is called, this method also instantiates it.

Return type:

DEIMOSDetectorMap

get_detector_par(det, hdu=None)[source]

Return metadata for the selected detector.

Parameters:
  • det (int) – 1-indexed detector number.

  • hdu (astropy.io.fits.HDUList, optional) – The open fits file with the raw image of interest. If not provided, frame-dependent parameters are set to a default.

Returns:

Object with the detector metadata.

Return type:

DetectorContainer

get_grating(filename)[source]

Instantiate grating (a ReflectionGrating instance) based on the grating using to collect the data provided by the filename.

Taken from xidl/DEEP2/spec2d/pro/deimos_omodel.pro and xidl/DEEP2/spec2d/pro/deimos_grating.pro

Parameters:

filename (str) – Name of the file with the grating metadata.

Returns:

The grating instance relevant to the data in filename. The returned object is the same as grating.

Return type:

ReflectionGrating

get_lamps(fitstbl)[source]

Extract the list of arc lamps used from header

Parameters:

fitstbl (astropy.table.Table) – The table with the metadata for one or more arc frames.

Returns:

List used arc lamps

Return type:

lamps (list)

get_maskdef_slitedges(ccdnum=None, filename=None, debug=None, trc_path=None, binning=None)[source]

Provides the slit edges positions predicted by the slitmask design using the mask coordinates already converted from mm to pixels by the method mask_to_pixel_coordinates.

If not already instantiated, the slitmask, amap, and bmap attributes are instantiated. If so, a file must be provided.

Parameters:
  • ccdnum (int) – Detector number

  • filename (str, optional) – The filename to use to (re)instantiate the slitmask and grating. Default is None, i.e., to use previously instantiated attributes.

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

Returns:

Three numpy.ndarray and a SlitMask. Two arrays are the predictions of the slit edges from the slitmask design and one contains the indices to order the slits from left to right in the PypeIt orientation

Return type:

tuple

get_mosaic_par(mosaic, hdu=None, msc_order=5)[source]

Return the hard-coded parameters needed to construct detector mosaics from unbinned images.

The parameters expect the images to be trimmed and oriented to follow the PypeIt shape convention of (nspec,nspat). For returned lists, the length of the list is the same as the number of detectors in the mosaic, and they are ordered by the detector number.

Parameters:
  • mosaic (tuple) – Tuple of detector numbers used to construct the mosaic. Must be one among the list of possible mosaics as hard-coded by the allowed_mosaics() function.

  • hdu (astropy.io.fits.HDUList, optional) – The open fits file with the raw image of interest. If not provided, frame-dependent detector parameters are set to a default. BEWARE: If hdu is not provided, the binning is assumed to be 1,1, which will cause faults if applied to binned images!

  • msc_order (int, optional) – Order of the interpolation used to construct the mosaic.

Returns:

Object with the mosaic and detector parameters.

Return type:

Mosaic

get_rawimage(raw_file, det)[source]

Read raw images and generate a few other bits and pieces that are key for image processing.

Data are unpacked from the multi-extension HDU. Function is based on pypeit.spectrographs.keck_lris.read_lris(), which was based on the IDL procedure readmhdufits.pro.

Warning

PypeIt currently cannot reduce images produced by reading the DEIMOS CCDs with the A+B amplifier or those taken in imaging mode. All image handling assumes DEIMOS images have been read with the B or A amplifier in the “Spectral” observing mode. This method will fault if this is not true based on the header keywords MOSMODE and AMPMODE.

Parameters:
  • raw_file (str) – File to read

  • det (int) – 1-indexed detector to read

Returns:

  • detector_par (pypeit.images.detector_container.DetectorContainer) – Detector metadata parameters.

  • raw_img (numpy.ndarray) – Raw image for this detector.

  • hdu (astropy.io.fits.HDUList) – Opened fits file

  • exptime (float) – Exposure time read from the file header

  • rawdatasec_img (numpy.ndarray) – Data (Science) section of the detector as provided by setting the (1-indexed) number of the amplifier used to read each detector pixel. Pixels unassociated with any amplifier are set to 0.

  • oscansec_img (numpy.ndarray) – Overscan section of the detector as provided by setting the (1-indexed) number of the amplifier used to read each detector pixel. Pixels unassociated with any amplifier are set to 0.

get_slitmask(filename: str)[source]

Parse the slitmask data from a DEIMOS file into slitmask, a SlitMask object.

Parameters:

filename (str) – Name of the file to read.

Returns:

The slitmask data read from the file. The returned object is the same as slitmask.

Return type:

SlitMask

get_telescope_offset(file_list)[source]

For a list of frames compute telescope pointing offset w.r.t. the first frame. Note that the object in the slit will appear moving in the opposite direction (=-tel_off)

Parameters:
  • file_list (list) – List of frames (including the path) for which telescope offset is desired.

  • used. (Both raw frames and spec2d files can be) –

Returns:

List of telescope offsets (in arcsec) w.r.t. the first frame

Return type:

numpy.ndarray

header_name = 'DEIMOS'

Name of the spectrograph camera or arm from the Header. Usually the INSTRUME card.

idname(ftype)[source]

Return the idname for the selected frame type for this instrument.

Parameters:

ftype (str) – Frame type, which should be one of the keys in FrameTypeBitMask.

Returns:

The value of idname that should be available in the PypeItMetaData instance that identifies frames of this type.

Return type:

str

init_meta()[source]

Define how metadata are derived from the spectrograph files.

That is, this associates the PypeIt-specific metadata keywords with the instrument-specific header cards using meta.

list_detectors(mosaic=False)[source]

List the names of the detectors in this spectrograph.

This is primarily used average_maskdef_offset() to measure the mean offset between the measured and expected slit locations.

Detectors separated along the dispersion direction should be ordered along the first axis of the returned array. For example, Keck/DEIMOS returns:

dets = np.array([['DET01', 'DET02', 'DET03', 'DET04'],
                 ['DET05', 'DET06', 'DET07', 'DET08']])

such that all the bluest detectors are in dets[0], and the slits found in detectors 1 and 5 are just from the blue and red counterparts of the same slit.

Parameters:

mosaic (bool, optional) – Is this a mosaic reduction? It is used to determine how to list the detector, i.e., ‘DET’ or ‘MSC’.

Returns:

The list of detectors in a numpy.ndarray. If the array is 2D, there are detectors separated along the dispersion axis.

Return type:

numpy.ndarray

mask_to_pixel_coordinates(x=None, y=None, wave=None, order=1, filename=None, corners=False)[source]

Convert the mask coordinates in mm to pixel coordinates on the DEIMOS detector.

If not already instantiated, the slitmask, grating, optical_model, and detector_map attributes are instantiated. If these are not instantiated, a file must be provided. If no arguments are provided, the function expects these attributes to be set and will output the pixel coordinates for the centers of the slits in the slitmask at the central wavelength of the grating.

Method generally expected to be executed in one of two modes:
  • Use the filename to read the slit mask and determine the detector positions at the central wavelength.

  • Specifically map the provided x, y, and wave values to the detector.

If arrays are provided for both x, y, and wave, the returned objects have the shape \(N_\lambda\times S_x\), where \(S_x\) is the shape of the x and y arrays.

Parameters:
  • x (array-like, optional) – The x coordinates in the slit mask in mm. Default is to use the center of the slits in the slitmask.

  • y (array-like, optional) – The y coordinates in the slit mask in mm. Default is to use the center of the slits in the slitmask.

  • wave (array-like, optional) – The wavelengths in angstroms for the propagated coordinates. If not provided, an array of wavelength covering the full DEIMOS wavelength range will be used.

  • order (int, optional) – The grating order. Default is 1.

  • filename (str, optional) – The filename to use to (re)instantiate the slitmask and grating. Default is to use previously instantiated attributes.

  • corners (bool, optional) – Instead of using the centers of the slits in the slitmask, return the detector pixel coordinates for the corners of all slits.

Returns:

Returns 5 arrays: (1-2) the x and y coordinates in the image plane in mm, (3) the detector (1-indexed) where the slit should land at the provided wavelength(s), and (4-5) the pixel coordinates (1-indexed) in the relevant detector.

Return type:

numpy.ndarray

Raises:

ValueError – Raised if the user provides one but not both of the x and y coordinates, if no coordinates are provided or available within the slitmask, or if the grating, amap or bmap haven’t been defined and not file is provided.

name = 'keck_deimos'

The name of the spectrograph. See Spectrographs for the currently supported spectrographs.

ndet = 8

Number of detectors for this instrument.

pypeit_file_keys()[source]

Define the list of keys to be output into a standard PypeIt file.

Returns:

The list of keywords in the relevant PypeItMetaData instance to print to the PypeIt Reduction File.

Return type:

list

ql_supported = True

Flag that PypeIt code base has been sufficiently tested with data from this spectrograph in quicklook mode that it is officially supported by the development team.

raw_header_cards()[source]

Return additional raw header cards to be propagated in downstream output files for configuration identification.

The list of raw data FITS keywords should be those used to populate the configuration_keys() or are used in config_specific_par() for a particular spectrograph, if different from the name of the PypeIt metadata keyword.

This list is used by subheader_for_spec() to include additional FITS keywords in downstream output files.

Returns:

List of keywords from the raw data files that should be propagated in output files.

Return type:

list

spec1d_match_spectra(sobjs)[source]

Match up slits in a SpecObjs file based on coords. Specific to DEIMOS.

Parameters:

sobjs (pypeit.specobjs.SpecObjs) – Spec1D objects

Returns:

array of indices for the blue detector, array of indices for the red (matched to the blue).

Return type:

tuple

subheader_for_spec(row_fitstbl, raw_header, extra_header_cards=None, allow_missing=False)[source]

Generate a dict that will be added to the Header of spectra files generated by PypeIt (e.g. SpecObjs). This version overrides the parent version to include KOA specific header cards.

Parameters:
  • row_fitstbl (dict-like) – Typically an astropy.table.Row or astropy.io.fits.Header with keys defined by define_core_meta().

  • raw_header (astropy.io.fits.Header) – Header that defines the instrument and detector, meaning that the header must contain the INSTRUME and DETECTOR header cards. If provided, this must also contain the header cards provided by extra_header_cards.

  • extra_header_cards (list, optional) – Additional header cards from raw_header to include in the output dictionary. Can be an empty list or None.

  • allow_missing (bool, optional) – Ignore any keywords returned by define_core_meta() are not present in row_fitstbl. Otherwise, raise PypeItError.

Returns:

Dictionary with data to include an output fits header file or table downstream.

Return type:

dict

supported = True

Flag that PypeIt code base has been sufficiently tested with data from this spectrograph that it is officially supported by the development team.

telescope = Parameter     Value                Default  Type        Callable ---------------------------------------------------------------- name          KECK                 KECK     str         False    longitude     -155.47833333333335  None     int, float  False    latitude      19.828333333333333   None     int, float  False    elevation     4160.000000000756    None     int, float  False    fratio        15                   None     int, float  False    diameter      10                   None     int, float  False    eff_aperture  72.3674              None     int, float  False   

Instance of TelescopePar providing telescope-specific metadata.

update_edgetracepar(par)[source]

This method is used in pypeit.edgetrace.EdgeTraceSet.maskdesign_matching() to update EdgeTraceSet parameters when the slitmask design matching is not feasible because too few slits are present in the detector.

Parameters:

par (pypeit.par.pypeitpar.EdgeTracePar) – The parameters used to guide slit tracing.

Returns:

pypeit.par.pypeitpar.EdgeTracePar The modified parameters used to guide slit tracing.

url = 'https://www2.keck.hawaii.edu/inst/deimos/'

Reference url

valid_configuration_values()[source]

Return a fixed set of valid values for any/all of the configuration keys.

Returns:

A dictionary with any/all of the configuration keys and their associated discrete set of valid values. If there are no restrictions on configuration values, None is returned.

Return type:

dict

pypeit.spectrographs.keck_deimos.deimos_read_1chip(hdu, chipno)[source]

Read one of the DEIMOS detectors

Parameters:
  • hdu (astropy.io.fits.HDUList) –

  • chipno (int) –

Returns:

data, oscan

Return type:

np.ndarray, np.ndarray

pypeit.spectrographs.keck_deimos.indexing(itt, postpix, det=None)[source]

Some annoying book-keeping for instrument placement.

Parameters:
  • itt (int) –

  • postpix (int) –

  • det (int, optional) –

pypeit.spectrographs.keck_deimos.load_wmko_std_spectrum(fits_file: str, outfile=None, pad=False, split=True)[source]

Load up a Standard spectrum generated by WMKO IDL scripts of the great Greg Wirth

The SpecObjs generated is checked that it is ready for fluxing

Parameters:
  • fits_file (str) – filename

  • outfile ([type], optional) – Write the SpecObjs object to a FITS file. Defaults to None.

  • pad (bool,optional) – True if the resulting SpecObjs should be padded to an even length. Defaults to False

  • split (bool,optional) – True if the resulting SpecObjs should be split into two detectors, False if it should be treated as a mosaic. Defaults to True.

Returns:

object holding the spectra

Return type:

specobjs.SpecObjs