pypeit.spectrographs.ldt_deveny module

Module for LDT/DeVeny specific methods.

The DeVeny spectrograph was built at Kitt Peak National Observatory (KPNO) and known as the White Spectrograph. It had a long career at the #1 36-inch and 84-inch telescopes there before being retired; Lowell Observatory acquired the spectrograph from KPNO on indefinite loan in 1998. A new CCD camera was built for it, and the instrument was further modified for installation on the 72-inch Perkins telescope in 2005. Following 8 years of service there, it was removed in 2013 for upgrades for installation on the Lowell Discovery Telescope (LDT) instrument cube. It has been in service since February 2015. The spectrograph was designed for f/7.5 telescope optics, and new re-imaging optics were designed and fabricated to match the spectrograph with LDT’s f/6.1 beam.

class pypeit.spectrographs.ldt_deveny.LDTDeVenySpectrograph[source]

Bases: Spectrograph

Child to handle LDT/DeVeny specific code

camera = 'DeVeny'

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

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

Check for frames of the provided type.

Parameters
  • ftype (str) – Type of frame to check. Must be a valid frame type; see frame-type Definitions.

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

  • exprng (list, optional) – Range in the allowed exposure time for a frame of type ftype. See check_frame_exptime().

Returns

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

Return type

numpy.ndarray

comment = 'LDT DeVeny Optical Spectrograph, 2015 - present'

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

compound_meta(headarr: list, meta_key: str) object[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_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

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_detector_par(det, hdu=None)[source]

Return metadata for the selected detector.

Warning

Many of the necessary detector parameters are read from the file header, meaning the hdu argument is effectively required for LTD/DeVeny. The optional use of hdu is only viable for automatically generated documentation.

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

  • hdu (astropy.io.fits.HDUList, optional) – The open fits file with the raw image of interest.

Returns

Object with the detector metadata.

Return type

DetectorContainer

get_lamps(fitstbl: Table) list[source]

Extract the list of arc lamps used from header

Note

There are some faint Cd and Hg lines in the DV9 spectra that are helpful for nailing down the wavelength calibration for that grating, but these lines are too faint / close to other lines for use with other gratings. This method loads the more detailed lists for DV9, but loads the usual line lists for all other gratings.

Parameters

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

Returns

List of the used arc lamps

Return type

list

get_rawimage(raw_file, det)[source]

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

For LDT/DeVeny, the LOIS control system automatically adjusts the DATASEC and OSCANSEC regions if the CCD is used in a binning other than 1x1. The get_rawimage() method in the base class assumes these sections are fixed and adjusts them based on the binning – an incorrect assumption for this instrument.

This method is a stripped-down version of the base class method and additionally does NOT send the binning to sec2slice().

Parameters
  • raw_file (str) – File to read

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

Returns

  • detector_par (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 in seconds.

  • 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.

header_name = 'Deveny'

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

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.

name = 'ldt_deveny'

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

ndet = 1

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

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

static rotate_trimsections(section_string: str, nspecpix: int)[source]

In order to orient LDT/DeVeny images into the PypeIt-standard configuration, frames are essentially rotated 90º clockwise. As such, \(x' = y\) and \(y' = -x\).

The TRIMSEC / BIASSEC FITS keywords in LDT/DeVeny data specify the proper regions to be trimmed for the data and overscan arrays, respectively, in the native orientation. This method performs the rotation and returns the slices for the Numpy image section required by the PypeIt processing routines.

The LDT/DeVeny FITS header lists the sections as '[SPEC_SEC,SPAT_SEC]'.

Parameters
  • section_string (str) – The FITS keyword string to be parsed / translated

  • nspecpix (int) – The total number of pixels in the spectral direction

Returns

Numpy image section needed by PypeIt

Return type

section (numpy.ndarray)

static scrub_isot_dateobs(dt_str: str)[source]

Scrub the input DATE-OBS for ingestion by AstroPy Time

The main issue this method addresses is that sometimes the LOIS software at LDT has roundoff abnormalities in the time string written to the DATE-OBS header keyword. For example, in one header 2020-01-30T13:17:010.0 was written, where the seconds has 3 digits – presumably the seconds field was constructed with a leading zero because sec was < 10, but when rounded for printing yielded “10.00”, producing a complete seconds field of 010.00.

This abnormality, along with a seconds field equaling 60.00, causes AstroPy’s Time parser to freak out with a ValueError. This method attempts to return the astropy.time.Time object directly, but then scrubs any values that cause a ValueError.

The scrubbing consists of deconstructing the string into its components, then carefully reconstructing it into proper ISO 8601 format. Also, some recursive edge-case catching is done, but at some point you just have to give up and go buy a lottery ticket.

If you have a truly bizarre DATE-OBS string, simply edit that keyword in the FITS header and then re-run PypeIt.

Parameters

dt_str (str) – Input datetime string from the DATE-OBS header keyword

Returns

The AstroPy Time object corresponding to the DATE-OBS input string

Return type

astropy.time.Time

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          LDT                  KECK     str         False    longitude     -111.42251499999998  None     int, float  False    latitude      34.744305000000004   None     int, float  False    elevation     2336.9999999992638   None     int, float  False    fratio        6.18                 None     int, float  False    diameter      4.25                 None     int, float  False    eff_aperture  12.37                None     int, float  False   

Instance of TelescopePar providing telescope-specific metadata.

tweak_standard(wave_in, counts_in, counts_ivar_in, gpm_in, meta_table)[source]

This routine is for performing instrument- and/or disperser-specific tweaks to standard stars so that sensitivity function fits will be well behaved.

These are tweaks needed by LDT/DeVeny for smooth sensfunc sailing.

Parameters
  • wave_in (numpy.ndarray) – Input standard star wavelengths (float, shape = (nspec,))

  • counts_in (numpy.ndarray) – Input standard star counts (float, shape = (nspec,))

  • counts_ivar_in (numpy.ndarray) – Input inverse variance of standard star counts (float, shape = (nspec,))

  • gpm_in (numpy.ndarray) – Input good pixel mask for standard (bool, shape = (nspec,))

  • meta_table (dict) – Table containing meta data that is slupred from the SpecObjs object. See unpack_object() for the contents of this table.

Returns

  • wave_out (numpy.ndarray) – Output standard star wavelengths (float, shape = (nspec,))

  • counts_out (numpy.ndarray) – Output standard star counts (float, shape = (nspec,))

  • counts_ivar_out (numpy.ndarray) – Output inverse variance of standard star counts (float, shape = (nspec,))

  • gpm_out (numpy.ndarray) – Output good pixel mask for standard (bool, shape = (nspec,))

url = 'https://lowell.edu/research/telescopes-and-facilities/ldt/deveny-optical-spectrograph/'

Reference url