pypeit.spectrographs.mmt_binospec module
Module for MMT/BINOSPEC specific methods.
- class pypeit.spectrographs.mmt_binospec.MMTBINOSPECSpectrograph[source]
Bases:
SpectrographChild to handle MMT/BINOSPEC specific code
- bino_get_slit_region(filename, det=None, Nx=4096, Ny=4112, pady=0)[source]
Compute the pixel-space rectangular regions for each slit in a Binospec mask.
This function reads the slitmask design from a FITS file (or an already-loaded SlitMask object), converts slit and object positions from mask coordinates to pixel coordinates, and determines the x/y pixel boundaries for each slit on the detector. It returns these boundaries along with the updated slitmask object.
- Parameters:
filename (
str) – Path to the slitmask FITS file. Must be provided unless the slitmask is already loaded via self.get_slitmask.det (
int, optional) – Detector number (1 or 2). Must be specified.Nx (
int, optional) – Detector size in the x-direction (default: 4096 pixels).Ny (
int, optional) – Detector size in the y-direction (default: 4112 pixels).pady (
float, optional) – Additional padding (in pixels) applied to the slit boundaries (default: 0).
- Returns:
region (
list) – A list containing: - slit_x_range : array of x-boundaries for each slit [Nslits, 2] - slit_y_range : array of y-boundaries for each slit [Nslits, 2] - x_slitobj_pix : array of x pixel positions for slit objects - y_slitobj_pix : array of y pixel positions for slit objectsslitmask (
SlitMask) – The updated SlitMask object containing slit geometry and metadata.
Notes
Converts mask coordinates to pixel coordinates using the appropriate scale factor.
Handles detector 2 by reversing slit order and applying a vertical flip.
Slit boundaries are clipped to remain within detector dimensions.
- bpm(filename, det, shape=None, msbias=None)[source]
Generate a default bad-pixel mask.
Loads a pre-built static BPM from the IDL pipeline calibration data (
badpix_binospec.fits+ hard-coded bad columns and detector trap regions frombino_mosaic.pro).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 (filenameusingget_image_shape()) must be provided.- Parameters:
filename (
stror 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:
- camera = 'BINOSPEC'
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:
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 typeftype. Seepypeit.core.framematch.check_frame_exptime().
- Returns:
Boolean array with the flags selecting the exposures in
fitstblthat areftypetype frames.- Return type:
- 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:
headarr (
list) – List of astropy.io.fits.Header objects.meta_key (
str) – Metadata keyword to construct.
- Returns:
Metadata value read from the header(s).
- Return type:
- config_specific_par(inp, inp_par=None)[source]
Modify the PypeIt parameters to hard-wired values used for specific instrument configurations.
- Parameters:
inp (
str,list, Path, astropy.io.fits.Header, astropy.table.Table) – Input filename, an astropy.io.fits.Header object, or a list of astropy.io.fits.Header objects. Or a row from the metadata table.inp_par (
ParSet, optional) – Parameter set used for the full run of PypeIt. If None, usedefault_pypeit_par().
- Returns:
The PypeIt parameter set adjusted for configuration specific parameter values.
- Return type:
- configuration_keys()[source]
Return the metadata keys that define a unique instrument configuration.
This list is used by
PypeItMetaDatato 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
PypeItMetaDataobject.- Return type:
- classmethod default_pypeit_par()[source]
Return the default parameters to use for this instrument.
- Returns:
Parameters required by all of PypeIt methods.
- Return type:
- 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:
- get_maskdef_slitedges(filename=None, det=1, debug=None, binning=None, trc_path=None)[source]
Provides the slit edges positions predicted by the slitmask design.
This method is not defined for all spectrographs. This base-class method raises an exception. This may be because
use_maskdesignhas been set to True for a spectrograph that does not support it.- Parameters:
filename (
str,list, optional:) – Name of the file holding the mask design info or the maskfile and wcs_file in that orderdet (
int, optional) – Detector numberdebug (
bool, optional) – Flag to run in debugging modetrc_path (str, optional) – Path to the first trace file used to generate the trace flat
binning (str, optional) – String with the comma-separated number of pixels binned in each dimension of the flat-field image. Order must be spectral then spatial.
- Returns:
top_edges (
numpy.ndarray) – Predicted locations of the top edges of the slits in spatial pixel coordinates.bot_edges (
numpy.ndarray) – Predicted locations of the bottom edges of the slits in spatial pixel coordinates.sortindx (
numpy.ndarray) – Indices of the slits in the providedslitmaskobject that orders the slits from left to right, in the PypeIt orientation.slitmask (
SlitMask) – Slit mask metadata read from the provided input file(s).
Notes
Edges are sorted by bottom edge y-coordinate to order slits spatially.
- get_rawimage(raw_file, det)[source]
Read raw images and generate a few other bits and pieces that are key for image processing.
- Parameters:
- 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 headerrawdatasec_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, det=1)[source]
Parse the slitmask data from a raw file into
slitmask, aSlitMaskobject.- Parameters:
- Returns:
The slitmask data read from the file. The returned object is the same as
slitmask.- Return type:
Notes
Target-slit alignment is characterized via distances from slit edges.
Slit corners and on-sky positions are stored for each target.
- header_name = 'Binospec'
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 = 'mmt_binospec'
The name of the spectrograph. See Spectrographs for the currently supported spectrographs.
- ndet = 2
Number of detectors for this instrument.
- nonlinearity_coeffs = array([[ 0.00000000e+00, 1.00400089e+00, -1.39235362e-06, 8.31711824e-12, -1.20653479e-17], [ 0.00000000e+00, 1.00361458e+00, -1.29223833e-06, 6.93723177e-12, -9.67406255e-18], [ 0.00000000e+00, 1.00269542e+00, -9.29361806e-07, 5.97902827e-12, -2.30257302e-17], [ 0.00000000e+00, 1.00339616e+00, -8.47134521e-07, 7.92441693e-12, -4.46542834e-17], [ 0.00000000e+00, 1.00727205e+00, -1.69093388e-06, 2.07225055e-11, -1.62655178e-16], [ 0.00000000e+00, 1.00858745e+00, -2.35668901e-06, 2.40641019e-11, -1.50286358e-16], [ 0.00000000e+00, 1.00728526e+00, -1.80779473e-06, 1.73427719e-11, -1.01685780e-16], [ 0.00000000e+00, 1.00845168e+00, -2.02050567e-06, 2.97587091e-11, -2.65508521e-16]])
- plot_mask(filename, det=None, save_dir=None)[source]
Plot the slit mask layout and target positions for one or both detectors.
This function retrieves slit region data for a given Binospec mask and plots the rectangular slit outlines and target positions for detector 1, detector 2, or both. It is useful for visually validating mask design and target alignment.
- Parameters:
- Returns:
- 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 inconfig_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:
- 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 MMT KECK str False longitude -110.87750000000003 None int, float False latitude 31.68094444444444 None int, float False elevation 2319.9999999995903 None int, float False fratio None None int, float False diameter 6.5 None int, float False eff_aperture None None int, float False
Instance of
TelescopeParproviding 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.EdgeTraceParThe modified parameters used to guide slit tracing.
- url = 'https://lweb.cfa.harvard.edu/mmti/binospec.html'
Reference url
- pypeit.spectrographs.mmt_binospec.binospec_read_amp(inp, ext)[source]
Read one amplifier of an MMT BINOSPEC multi-extension FITS image
- Parameters:
inp (str,
astropy.io.fits.HDUList) – The input FITS file name or already opened HDU list.ext (
int) – FITS extension to read
- Returns:
data (
numpy.ndarray) – Array with data from the data section of the image.overscan (
numpy.ndarray) – Array with the overscan section of the image.datasec (
str) – String with the data section in IRAF format, e.g. ‘[x1:x2,y1:y2]’.biassec (
str) – String with the bias section in IRAF format, e.g. ‘[x1:x2,y1:y2]’.