mangadap.drpfits module

Defines a class used to interface with files produced in the 3D phase of the MaNGA Data Reduction Pipeline (DRP).

License:
Copyright (c) 2016, SDSS-IV/MaNGA Pipeline Group
Licensed under BSD 3-clause license - see LICENSE.rst
Source location:
$MANGADAP_DIR/python/mangadap/drpfits.py
Imports and python version compliance:
from __future__ import division
from __future__ import print_function
from __future__ import absolute_import

import sys
if sys.version > '3':
    long = int

import time
import os.path
import numpy
import warnings

from scipy import sparse
from astropy.io import fits
from astropy.wcs import WCS

from .util.parser import arginp_to_list
from .util.covariance import Covariance
from .config.defaults import default_idlutils_dir, default_drp_version
from .config.defaults import default_redux_path, default_drp_directory_path
from .config.defaults import default_cube_pixelscale, default_cube_width_buffer
from .config.defaults import default_cube_recenter, default_regrid_rlim
from .config.defaults import default_regrid_sigma
from .config.defaults import default_manga_fits_root
from .util.bitmask import BitMask

Class usage examples:

This class provides a basic interface with DRP LOG* files. It provides a number of
Revision history:
20 Nov 2014: Original implementation by K. Westfall (KBW)
12 Feb 2014: (KBW) Added DRPFits.directory_path()
20 Feb 2015: (KBW) Add covariance calculation to DRPFits
19 Mar 2015: (KBW) Added redux_path to DRPFits. Re-arranged arguments in drpfits_list(), made drpver optional, and added redux_path
22 May 2015: (KBW) Sphinx documentation. Changed DRPFits.w to DRPFits.wcs.
26 May 2015: (KBW) Added checksum=True when opening the DRP file.
04 Jun 2015: (KBW) Moved parse_drp_file_name to mangadap.util.parser.parse_drp_file_name()
15 Jun 2015: (KBW) Moved functions that return default values (like DRPFits._default_pixelscale()) to mangadap.config.defaults
05 Aug 2015: (KBW) Changed mode testing to be more robust. Added directory_path keyword to drpfits_list(). Changed how directory path is set; previously required drpver and redux_path defaults, even if directory_path was provided directly. May need to add checks to other code to make sure drpver and redux_path are not None when directory_path has been directly defined.
28 Aug 2015: (KBW) Added usage of mangadap.config.defaults.default_manga_fits_root()
15 Feb 2016: (KBW) Added DRPFits.__getitem__() function
17 Feb 2016: (KBW) Converted drpfile class name to DRPFits
23 Mar 2016: (KBW) Added functionality that abstracts the difference between the RSS and CUBE file formats at the user level. CUBE files are now restructured to matched the intended orientation provided by the DRP; i.e., [x,y,lambda]. RSS files are left the same, which is [fiber, lambda]. Documentation. Testing, particularly of x,y order.
10 Nov 2016: (KBW) Included ‘DISP’ in spectral arrays.
30 Nov 2016: (KBW) Include DRPFits.spectral_resolution(), which returns spectral resolution cube or vector independent of MPL
01 Dec 2016: (KBW) Added DRPFits.spectral_resolution_header().
06 Dec 2016: (KBW) Removed wavelength_mask function, now uses mangadap.util.pixelmask.SpectralPixelMask. Moved the main functionality of DRPFits.copy_to_array() and DRPFits.copy_to_masked_array() to mangadap.util.fitsutil.DAPFitsUtil, what’s left are wrapper functions for the more general functions in mangadap.util.fitsutil.DAPFitsUtil.
17 Feb 2017: (KBW) Return nominal inverse variance in DRPFits.regrid_wavelength_plane() if requested.
17 May 2017: (KBW) Include a response function in DRPFits.flux_stats() and DAPFits.mean_sky_coordinates().
21 Aug 2017: (KBW) In spectral resolution function, select pre- vs. post-pixelized Gaussian respresentation.

Todo

  • Calculation in DRPFits._cube_dimensions() will only be correct if the WCS coordinates have no rotation.
  • Further optimize calculation of transfer matrix
  • Make DRP file class flexible to linear or log-linear wavelength sampling? Incorporate into MODE?
  • Reconstructed broad-band images and PSFs are not restructured in the CUBE files! This is why the are transposed in mangadap.drpfits.DRPFits.gri_composite().
  • Image reconstruction has transpose sense wrt DRP output!
  • Add logging
  • Need to be clear about which functions use the RSS spectra to create CUBE related data, like the covariance matrix and instrumental dispersion calculations.
  • Computing the approximate covariance cube is currently not possible with only the CUBE on disk. There’s a logic problem that needs to be fixed.
class mangadap.drpfits.DRPFits(plate, ifudesign, mode, drpver=None, redux_path=None, directory_path=None, read=False, checksum=False)[source]

Bases: object

A general purpose class used to interface with a MaNGA DRP file.

Parameters:
  • plate (int) – Plate number
  • ifudesign (int) – IFU design
  • mode (str) – 3D mode of the DRP file; must be either ‘RSS’ or ‘CUBE’
  • drpver (str) – (Optional) DRP version, which is used to define the default DRP redux path. Default is defined by mangadap.config.defaults.default_drp_version()
  • redux_path (str) – (Optional) The path to the top level directory containing the DRP output files for a given DRP version. Default is defined by mangadap.config.defaults.default_redux_path().
  • directory_path (str) – (Optional) The exact path to the DRP file. Default is defined by mangadap.config.defaults.default_drp_directory_path().
  • read (bool) – (Optional) Read the DRP file upon instantiation of the object.
  • checksum (bool) – (Optional) Check for file corruption.
Raises:

ValueError – Raised if mode is not ‘RSS’ or ‘CUBE’.

plate, ifudesign

Plate and IFU designation

Type:int
mode

3D mode of the DRP file, see above

Type:str
drpver

DRP version, which is used to define the default DRP redux path, see above.

Type:str
redux_path

The path to the top level directory containing the DRP output files for a given DRP version, see above.

Type:str
directory_path

The exact path to the DRP file, see above.

Type:str
pixelscale

Pixel scale used during the CUBE reconstruction.

Type:float
recenter

If False, the coordinates in the XPOS and YPOS extensions of the DRP file are assumed to be centered at 0,0. If True, the XPOS and YPOS can have any center, and the center of the CUBE is set to be approximately the center of the range in XPOS,YPOS.

Type:bool
width_buffer

The number of pixels to add to the width of the cube in addition to the range needed to cover XPOS and YPOS.

Type:float
xs, ys

The starting on-sky coordinate of the reconstructed image defined by the bottom corner of the first pixel, not its center!

Type:float
nx, ny

The size (number of pixel in x and y) of the reconstructed image.

Type:int
regrid_T

Transfer matrix \({\mathbf T}\) such that:

\[{\mathbf T} \times {\mathbf F} = {\mathbf I}\]

where \({\mathbf F}\) is the vector of fluxes in a single wavelength channel for all the fiber measurements in the field of view and \({\mathbf I}\) is the pre-formatted reconstructed image of that wavelength channel. Saved such that repeat calls to create T for a given wavelength channel do not result in repeat calculations.

Type:scipy.sparse.csr_matrix
regrid_channel

Wavelength channel for which regrid_T has been calculated.

Type:int
regrid_rlim

The limiting radius of the Gaussian interpolation kernel used during image construction in arcseconds.

Type:float
regrid_sigma

The sigma of the Gaussian interpolation kernel used during image construction in arcseconds.

Type:float
sigma_rho

The sigma, \(\sigma_\rho\), of the Gaussian function used to approximate the trend of the correlation coefficient \(\rho\) with pixel separation as stored in cov_rho. That is:

\[\rho_{ij} = \exp\left(\frac{-d^2_{ij}}{2 \sigma^2_\rho}\right)\]

where \(d_{ij}\) is the distance between pixels \(i\) and \(j\).

Type:float
cov_rho

The matrix \({\mathbf R}\) containing the correlation coefficents, \(\rho\), between elements of the covariance matrix, \({\mathbf C}\), as approximated using the parameterization of \(\rho\) with pixel separation. This matrix will be independent of wavelength. In general,

\[\begin{split}{\mathbf C} = \left[ \begin{array}{rrrr} ... & ... & ... & ... \\ ... & \sigma^2_i & \rho_{ij}\sigma_i\sigma_j & ... \\ ... & \rho_{ji}\sigma_j\sigma_i & \sigma^2_j & ... \\ ... & ... & ... & ... \end{array} \right]\end{split}\]

such that

\[\begin{split}{\mathbf R} = \left[ \begin{array}{rrrr} ... & ... & ... & ... \\ ... & 1.0 & \rho_{ij} & ... \\ ... & \rho_{ji} & 1.0 & ... \\ ... & ... & ... & ... \end{array} \right].\end{split}\]
Type:scipy.sparse.csr_matrix
bitmask

Object used to interpret the DRP bit mask values in the MASK extension.

Type:DRPFitsBitMask
hdu

HDUList read from the DRP file

Type:astropy.io.fits.hdu.hdulist.HDUList
ext

List of fits extensions in the file

Type:list
checksum

Flag to check for file corruption when opening the HDU.

Type:bool
wcs

WCS object based on WCS keywords in the header of the FLUX extension.

Type:astropy.wcs.wcs.WCS
shape

Shape of the main data arrays

Type:tuple
spatial_shape

Shape of the spatial axes only. For RSS files, this is a single element with the number of fibers; for CUBE files, this has the x and y dimensions of the data cube. These are transposed w.r.t. the read-in DRP file!

Type:tuple
nspec

Number of spectra in the DRP file; this is just:

self.nspec = numpy.prod(self.spatial_shape)
Type:int
spatial_index

Array with tuples used to select spectra at specific locations within the data array. This is mainly useful in CUBE mode, where this provides the indices in the spatial coordinates. The order is \((x,y)\); i.e. this is different that what you get if you read the DRP CUBE fits file directly using astropy.io.fits. In RSS mode, this is just the index of the spectrum in the 2D array. See: select().

Type:numpy.ndarray
spectral_arrays

List of viable keywords for the data arrays in the DRP file. For CUBE files, these are ‘FLUX’, ‘IVAR’, and ‘MASK’; for RSS files, this also includes ‘DISP’, ‘XPOS’, and ‘YPOS’.

Type:list
dispaxis

Index of the axis with the spectral channels. The internal data structure always has dispaxis as the last axis in the array. So dispaxis is 2 for CUBE files and 1 for RSS files. This means that the internal data array restructures the input fits data for the CUBE files.

Type:int
nwave

The number of wavelength channels; this is just:

self.nwave = self.shape[self.dispaxis]
Type:int
_approximate_covariance_matrix(channel, pixelscale, recenter, width_buffer, rlim, sigma, sigma_rho, csr=False, quiet=False)[source]

Return an approximate calculation of the covariance matrix assuming

\[C_{ij} = \frac{\rho_{ij}}{(V^{-1}_{ii} V^{-1}_{jj})^{1/2}}\]

where \(\rho_{ij}\) is approximated by a Gaussian with a standard deviation defined by the provided sigma_rho. See the description of attributes cov_rho and sigma_rho. For this to work, the variance matrix of the reconstructed image must have been already calculated, meaning that this approach is most appropriately used with the ‘CUBE’ files. For ‘RSS’ files, the covariance matrix is determined in this case by calling covariance_matrix() on the ‘CUBE’ file, which will raise an exception if the parameters defining the dimensions of the reconstructed image and kernel are not the defaults.

In general, this function should not be used with ‘RSS’ files for efficiency sake.

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the transfer matrix.
  • pixelscale (float) – Desired pixel scale in arcsec
  • recenter (bool) – Flag to recenter the coordinate system
  • width_buffer (int) – Number of pixels to use as buffer for the image reconstruction
  • rlim (float) – The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – The sigma of the image reconstruction kernel in arcseconds.
  • sigma_rho (float) – The sigma of the Gaussian function used to approximate the trend of the correlation coefficient with pixel separation.
  • csr (bool) – (Optional) Instead of reaturning a mangadap.util.covariance.Covariance object, return the covariance matrix as a scipy.sparse.csr_matrix object. Primarily used by covariance_cube() for collating the covariance matrix of each wavelength channel before combining them into a single mangadap.util.covariance.Covariance object
  • quiet (bool) – (Optional) Suppress terminal output
Returns:

mangadap.util.covariance.Covariance or scipy.sparse.csr_matrix: The covariance matrix for the designated wavelength channel. The return type depends on csr.

Raises:

ValueError – Raised for ‘CUBE’ files when the input parameters are not the same as the defaults.

_covariance_wavelength(waverange=None, response_func=None, per_pixel=True, flag=None)[source]

Determine the wavelength at which to calculate the covariance matrix.

Parameters:
  • waverange (array-like) – (Optional) Starting and ending wavelength over which to calculate the statistics. Default is to use the full wavelength range.
  • response_func (array-like) – (Optional) A two-column array with the wavelength and transmission of a broad-band response function to use for the calculation.
  • per_pixel (bool) – (Optional) When providing a response function, continue to calculate the statistics per pixel, as opposed to per angstrom. Default is to compute the statistics on a per pixel basis.
  • flag (str or list) – (Optional) (List of) Flag names that are considered when deciding if a pixel should be masked. The names must be a valid bit name as defined by bitmask (see DRPFitsBitMask).
Returns:

The response-weighted center of the wavelength region used to calculate the S/N, which will be where the covariance matrix is calculated.

Return type:

float

_cube_dimensions(pixelscale=None, recenter=None, width_buffer=None, redo=False)[source]

Determine the on-sky dimensions of the reconstructed image for all wavelength channels and save them in xs, ys, nx, and ny.

For CUBE files, these dimensions are drawn directly from the WCS keywords in the header of the FLUX extension of the DRP fits file. In this case, any entered parameters are ignored and the class attributes are set to the default values used by the DRP.

Warning

The calculation for the CUBE files is only valid if the WCS coordinate system has no rotation.

For RSS files, the dimensions are determined using the data in the ‘XPOS’ and ‘YPOS’ extensions and the same algorithm used by the DRP; however, it is possible to provide different parameters that will alter the dimensions.

See: pixelscale, recenter, width_buffer.

Parameters:
  • pixelscale (float) – (Optional) Desired pixel scale in arcsec
  • recenter (bool) – (Optional) Flag to recenter the coordinate system
  • width_buffer (int) – (Optional) Number of pixels to use as buffer for the image reconstruction
  • redo (bool) – (Optional) Force the recalculation of the cube dimensions if they are already defined.
_cube_dimensions_correct(pixelscale, recenter, width_buffer)[source]

Check that the saved parameters that define the cube dimensions are the same as the desired values.

Parameters:
  • pixelscale (float) – Desired pixel scale in arcsec
  • recenter (bool) – Flag to recenter the coordinate system
  • width_buffer (int) – Number of pixels to use as buffer for the image reconstruction
Returns:

True if the saved and desired values are the same.

Return type:

bool

_cube_dimensions_undefined()[source]

Return True if any of the cube dimensions are None.

_fix_header()[source]

Use of ‘degrees’ in early versions of the DRP did not adhere to the fits standard causing the astropy.wcs.wcs.WCS to fail when initialized; see, e.g., world_mesh(). This changes the units to be ‘deg’ instead.

Note

This function is obsolete as of v1_5_1 of the DRP and is not actively called in this implementation; see source for world_mesh().

_formal_covariance_matrix(channel, pixelscale, recenter, width_buffer, rlim, sigma, csr=False, quiet=False)[source]

Return the formal covariance matrix as defined by

\[{\mathbf C} = {\mathbf T} \times {\mathbf \Sigma} \times {\mathbf T}^{\rm T},\]

where \({\mathbf \Sigma}\) is the covariance matrix for the ‘RSS’ spectra for the specified wavelength channel. For ‘CUBE’ files, the function will attempt to use the ‘RSS’ counterpart of the file to produce the transfer matrix, \({\mathbf T}\). In this case, the input parameters must be the defaults.

Note

The current DRP does not produce spectral covariance matrices for the ‘RSS’ spectra. Here, it is assumed that the spectral covariance matrix is zero everywhere except along the diagonal, which contains the inverse of the values in the ‘IVAR’ extension.

The use of this function should not be used with the ‘CUBE’ files for efficiency sake.

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the transfer matrix.
  • pixelscale (float) – Desired pixel scale in arcsec
  • recenter (bool) – Flag to recenter the coordinate system
  • width_buffer (int) – Number of pixels to use as buffer for the image reconstruction
  • rlim (float) – The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – The sigma of the image reconstruction kernel in arcseconds.
  • csr (bool) – (Optional) Instead of reaturning a mangadap.util.covariance.Covariance object, return the covariance matrix as a scipy.sparse.csr_matrix object. Primarily used by covariance_cube() for collating the covariance matrix of each wavelength channel before combining them into a single mangadap.util.covariance.Covariance object
  • quiet (bool) – (Optional) Suppress terminal output
Returns:

mangadap.util.covariance.Covariance or scipy.sparse.csr_matrix: The covariance matrix for the designated wavelength channel. The return type depends on csr.

Raises:

ValueError – Raised for ‘CUBE’ files when the input parameters are not the same as the defaults.

_generate_spatial_index()[source]

Generate the tuples with the list of original indices in the input DRP file; see spatial_index.

_init_regrid_pars(pixelscale, recenter, width_buffer, rlim, sigma)[source]

Return the regridding parameters. If provided on input, the same value is returned. Otherwise, the returned values are the defaults. See the defaults defined in mangadap.config.defaults.default_cube_pixelscale(), mangadap.config.defaults.default_cube_recenter(), mangadap.config.defaults.default_cube_width_buffer(), mangadap.config.defaults.default_regrid_rlim(), and mangadap.config.defaults.default_regrid_sigma().

Parameters:
  • pixelscale (float) – Desired pixel scale in arcsec
  • recenter (bool) – Flag to recenter the coordinate system
  • width_buffer (int) – Number of pixels to use as buffer for the image reconstruction
  • rlim (float) – The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – The sigma of the image reconstruction kernel in arcseconds.
Returns:

The validated regridding parameters in the same order as listed by the function arguments.

Return type:

float, bool, int, float, float

_interpolated_response_function(response_func)[source]
_regrid_defaults(pixelscale, recenter, width_buffer, rlim, sigma)[source]

Check that the saved regridding parameters are the same as the defaults. See mangadap.config.defaults.default_cube_pixelscale(), mangadap.config.defaults.default_cube_recenter(), mangadap.config.defaults.default_cube_width_buffer(), mangadap.config.defaults.default_regrid_rlim(), and mangadap.config.defaults.default_regrid_sigma().

Parameters:
  • pixelscale (float) – Desired pixel scale in arcsec
  • recenter (bool) – Flag to recenter the coordinate system
  • width_buffer (int) – Number of pixels to use as buffer for the image reconstruction
  • rlim (float) – The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – The sigma of the image reconstruction kernel in arcseconds.
Returns:

True if the saved and default values are the same.

Return type:

bool

_regrid_kernel_correct(pixelscale, rlim, sigma)[source]

Check that the saved parameters used to define the image reconstruction kernel are the same as the desired values.

Parameters:
  • pixelscale (float) – Desired pixel scale in arcsec
  • rlim (float) – The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – The sigma of the image reconstruction kernel in arcseconds.
Returns:

True if the saved and desired values are the same.

Return type:

bool

_regrid_transfer_correct(channel, pixelscale, rlim, sigma)[source]

Check that the saved parameters used to construct the transfer matrix, \({\mathbf T}\), are the same as the desired values.

Parameters:
  • channel (int) – Index of the wavelength channel to reconstruct
  • pixelscale (float) – Desired pixel scale in arcsec
  • rlim (float) – The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – The sigma of the image reconstruction kernel in arcseconds.
Returns:

True if the saved and desired values are the same.

Return type:

bool

_regrid_transfer_undefined()[source]

Return True if regrid_T is None.

_set_spectral_arrays()[source]

Set the list of extensions with spectral data

_set_variance_correlation(sigma_rho, pixelscale=None, recenter=None, width_buffer=None, rlim=None, sigma=None, redo=False)[source]

Produce cov_rho based on the provided sigma_rho.

By default, the details of the cube dimensions should be the same as the DRP used to produce the ‘CUBE’ files from the ‘RSS’ spectra; however, these can be changed.

The resulting cov_rho is independent of wavelength and can be used in combination with the inverse-variance produced by the DRP to yield a wavelength-dependent covariance matrix that is close to the formal calculation.

See: cov_rho, attr:sigma_rho, pixelscale, recenter, width_buffer, rlim, sigma.

Warning

sigma is not actually used by this function.

Todo

Test the relation between sigma, rlim, and sigma_rho. It may be that sigma_rho should/can be determined by sigma.

Parameters:
  • sigma_rho (float) – The sigma of the Gaussian function used to approximate the trend of the correlation coefficient with pixel separation.
  • pixelscale (float) – (Optional) Desired pixel scale in arcsec
  • recenter (bool) – (Optional) Flag to recenter the coordinate system
  • width_buffer (int) – (Optional) Number of pixels to use as buffer for the image reconstruction
  • rlim (float) – (Optional) The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – (Optional) The sigma of the image reconstruction kernel in arcseconds.
  • redo (bool) – (Optional) Force the recalculation of the cube dimensions if they are already defined.
_spectral_resolution_extension(ext=None, pre=False)[source]

Determine the spectral resolution channel to use.

_variance_correlation_correct(sigma_rho, pixelscale, rlim, sigma)[source]

Check that the saved parameters used to construct the correlation coefficient matrix, \({\mathbf R}\), are the same as the desired values.

Parameters:
  • sigma_rho (float) – The sigma of the Gaussian function used to approximate the trend of the correlation coefficient with pixel separation.
  • pixelscale (float) – Desired pixel scale in arcsec
  • rlim (float) – The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – The sigma of the image reconstruction kernel in arcseconds.
Returns:

True if the saved and desired values are the same.

Return type:

bool

_variance_correlation_undefined()[source]

Return True if cov_rho is None.

binned_on_sky_area(bin_indx, x=None, y=None)[source]

Compute the on-sky area of a set of binned spectra. For CUBE files, this is just the number of spaxels in the bin times the spaxel area (as given by pixelscale).

For RSS files, this will try to calculate the overlapping area of the fibers using the shapely python package:

  • The fibers “beams” are all renormalized to have an area of pi arcsec^2 by the DRP, so it’s radius 1 arcsec
  • This function will provide the total area, not the integration-weighted effective area.
Parameters:
  • bin_indx (array-like) – A vector with size \(N_{\rm spec}\) the gives which spaxels or fibers were included in each bin. Valid bins have indices of \(\geq 0\).
  • x (array-like) – (Optional) On-sky \(x\) coordinate. Default is to calculate \(x\) and \(y\) using mean_sky_coordinates() with no arguments.
  • y (array-like) – (Optional) On-sky \(y\) coordinate. Default is to calculate \(x\) and \(y\) using mean_sky_coordinates() with no arguments.
Returns:

The on-sky area of each bin.

Return type:

numpy.ndarray

can_compute_covariance
static check_mode(mode)[source]

Check that the mode is valid.

Parameters:
  • mode (str) – Mode value to check. Valid modes are CUBE and
  • RSS.
copy_to_array(ext='FLUX', waverange=None)[source]

Wrapper for mangadap.util.fitsutil.DAPFitsUtil.copy_to_array() specific for DRPFits.

copy_to_masked_array(ext='FLUX', flag=None, waverange=None)[source]

Wrapper for mangadap.util.fitsutil.DAPFitsUtil.copy_to_masked_array() specific for DRPFits.

covariance_cube(channels=None, pixelscale=None, recenter=None, width_buffer=None, rlim=None, sigma=None, sigma_rho=None, csr=False, quiet=False)[source]

Return the covariance matrices for all wavelength channels; see covariance_matrix().

Parameters:
  • pixelscale (float) – (Optional) Desired pixel scale in arcsec
  • recenter (bool) – (Optional) Flag to recenter the coordinate system
  • width_buffer (int) – (Optional) Number of pixels to use as buffer for the image reconstruction
  • rlim (float) – (Optional) The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – (Optional) The sigma of the image reconstruction kernel in arcseconds.
  • sigma_rho (float) – (Optional) The sigma of the Gaussian function used to approximate the trend of the correlation coefficient with pixel separation.
  • csr (bool) – (Optional) Instead of reaturning a mangadap.util.covariance.Covariance object, return a numpy.ndarray of the covariance matrices for each channel, which are scipy.sparse.csr_matrix objects. Primarily used by covariance_cube() for collating the covariance matrix of each wavelength channel before combining them into a single mangadap.util.covariance.Covariance object
  • quiet (bool) – (Optional) Suppress terminal output
Returns:

mangadap.util.covariance.Covariance or numpy.ndarray: The return type depends on csr: if True, the returned object is an ndarray of scipy.sparse.csr_matrix types.

Raises:

Exception – Raised for ‘CUBE’ files when the input parameters are not the same as the defaults.

covariance_matrix(channel, pixelscale=None, recenter=None, width_buffer=None, rlim=None, sigma=None, sigma_rho=None, csr=False, quiet=False)[source]

Return the covariance matrix for the specified wavelength channel.

For a regrided cube image with \(N_x\times N_y\) pixels, the covariance matrix has a size that is \((N_x N_y)\times (N_x N_y)\); however, the majority of these pixels will be zero. Therefore, the covariance matrix is stored as a sparse matrix and interfaced with using the mangadap.util.covariance.Covariance object class.

The value of the covariance matrix at pixel \((i,j)\) is the covariance between pixels \((n_{x,0},n_{y,0})\) and \((n_{x,1},n_{y,1})\) at the specified wavelength channel of the reconstructed CUBE image, where

\[\begin{split}n_{x,0} &= \lfloor i / N_y \rfloor \\ n_{y,0} &= i - n_{x,0} N_y \\ n_{x,1} &= \lfloor j / N_y \rfloor \\ n_{y,1} &= j - n_{x,1} N_y\end{split}\]

and \(\lfloor m\rfloor\) is the “floor” of \(m\). The diagonal of the covariance matrix (\(i=j\)) should directly provide the inverse of the IVAR values provided by the DRP.

Warning

THIS IS IMPORTANT. Because the sense of the pixel coordinates, \((n_x,n_y)\), is flipped with respect to a direct read of the DRP fits file using astropy.io.fits, it is important that you appreciate this ordering in handling the coordinates in the output covariance matrix using the equation above. I.e., if you want to get the covariance between two pixels in the DRP file, make sure you understand which is the \(x\) pixel and which is the \(y\) pixel!

You can compare the CUBE provide inverse variance matrix and the inverse variance matrix provided along the diagonal of a covariance matrix calculated using this function as follows:

from matplotlib import pyplot
from mangadap.drpfits import DRPFits
from astropy.io import fits

drpf = DRPFits(7495, 3701, 'CUBE', read=False)
hdu = fits.open(drpf.file_path())
drpf = DRPFits(7495, 3701, 'RSS', read=True)

channel = 2000
C = drpf.covariance_matrix(channel)

ivar = numpy.diag(C.toarray()).reshape(-1, \
                        numpy.sqrt(C.shape[0])).T.copy()
indx = ivar > 0
ivar[indx] = 1./ivar[indx]
ivar[numpy.invert(indx)] = 0.

fig = pyplot.figure(figsize=(pyplot.figaspect(1)))
ax = fig.add_axes([0.05, 0.3, 0.35, 0.35])
ax.imshow(ivar, origin='lower', interpolation='nearest',
          cmap='inferno', aspect='auto')
ax.set_title('regrid')
ax = fig.add_axes([0.5, 0.3, 0.35, 0.35])
cax = fig.add_axes([0.86, 0.3, 0.01, 0.35])
ax.set_title('DRP')
cs = ax.imshow(hdu['IVAR'].data[channel,:,:], origin='lower',
               interpolation='nearest', cmap='inferno')
pyplot.colorbar(cs, cax=cax)
pyplot.show()

Differences between the inverse variance values should be small (\(\approx 10^{-3}\), which is about 1 part in \(10^7\)).

The covariance matrix can be generated in one of two ways:

1. If sigma_rho is None, the returned matrix is the formal covariance matrix defined as

\[{\mathbf C} = {\mathbf T} \times {\mathbf \Sigma} \times {\mathbf T}^{\rm T},\]

where \({\mathbf \Sigma}\) is the covariance matrix for the ‘RSS’ spectra for the specified wavelength channel. For ‘CUBE’ files, the function will attempt to use the ‘RSS’ counterpart of the file to produce the transfer matrix, \({\mathbf T}\). In this case, the input parameters must be the defaults.

Note

The current DRP does not produce spectral covariance matrices for the ‘RSS’ spectra. Here, it is assumed that the spectral covariance matrix is zero everywhere except along the diagonal, which contains the inverse of the values in the ‘IVAR’ extension.

2. If sigma_rho is not None, the returned matrix is an approximation of the covariance matrix determined by sigma_rho. See the description of cov_rho and sigma_rho. For this to work, the variance matrix of the reconstructed image must have been already calculated, meaning that this approach is most appropriately used with the ‘CUBE’ files. For ‘RSS’ files, the covariance matrix is determined in this case by calling covariance_matrix() on the ‘CUBE’ file, which will raise an exception if the parameters defining the dimensions of the reconstructed image and kernel are not the defaults.

In general, use of sigma_rho should not be used with ‘RSS’ files and the use of the formal calculation should not be used with the ‘CUBE’ files for efficiency.

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the transfer matrix.
  • pixelscale (float) – (Optional) Desired pixel scale in arcsec
  • recenter (bool) – (Optional) Flag to recenter the coordinate system
  • width_buffer (int) – (Optional) Number of pixels to use as buffer for the image reconstruction
  • rlim (float) – (Optional) The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – (Optional) The sigma of the image reconstruction kernel in arcseconds.
  • sigma_rho (float) – (Optional) The sigma of the Gaussian function used to approximate the trend of the correlation coefficient with pixel separation.
  • csr (bool) – (Optional) Instead of reaturning a mangadap.util.covariance.Covariance object, return the covariance matrix as a scipy.sparse.csr_matrix object. Primarily used by covariance_cube() for collating the covariance matrix of each wavelength channel before combining them into a single mangadap.util.covariance.Covariance object
  • quiet (bool) – (Optional) Suppress terminal output
Returns:

mangadap.util.covariance.Covariance or scipy.sparse.csr_matrix: The covariance matrix for the designated wavelength channel. The return type depends on csr.

Raises:

ValueError – Raised for ‘CUBE’ files when the input parameters are not the same as the defaults.

Todo

Need to make sure that the correct masks are being used for the RSS files. Should be 3DREJECT but nothing else?

created_today()[source]

Return True if the file was created today based on its time stamp.

Warning

Intended for use at the survey-level for a daily run of the DAP. Likely not to be used though.

static default_paths(plate, ifudesign, mode, drpver=None, redux_path=None, directory_path=None, output_file=None)[source]

Return the primary directory and file name with the DRP fits LOG-binned file.

plate (int): Plate number of the observation. ifudesign (int): IFU design number of the observation. mode (str): 3D mode of the DRP file; must be either ‘RSS’ or

‘CUBE’
drpver (str): (Optional) DRP version. Default set by
mangadap.config.defaults.default_drp_version().
redux_path (str): (Optional) The path to the top level
directory containing the DRP output files for a given DRP version. Default is defined by mangadap.config.defaults.default_redux_path().
directory_path (str): (Optional) The exact path to the
DAP reduction assessments file. Default set by mangadap.config.defaults.default_dap_common_path().
output_file (str): (Optional) The name of the file with
the DRP data. Default set by mangadap.config.defaults.default_manga_fits_root().
Returns:Two strings with the path to and name of the DRP data file.
Return type:str
static do_not_fit_flags()[source]

Return the maskbit names that should not be fit.

static do_not_stack_flags()[source]

Return the maskbit names that should not be stacked.

file_name()[source]

Return the name of the DRP file

file_path()[source]

Return the full path to the DRP file

finding_chart_path()[source]

Return the full path to the PNG finding chart for the targetted object.

Todo

  • Move this to the defaults file?
flux_stats(waverange=None, response_func=None, per_pixel=True, flag=None, covar=False, correlation=False, covar_wave=None)[source]

Compute the mean flux, propagated error in the mean flux, and mean S/N over the specified wavelength range; if the wavelength range is not specified, the quantities are calculated over the full spectral range.

If an RSS file, no covariance is calculated.

If a CUBE file and covar is True, the code will calculate the covariance matrix at the specified wavelength (see covariance_matrix()). If covar_wave is not provided, the covariance is calculated at the center wavelength of the provided wavelength range.

Parameters:
  • waverange (array-like) – (Optional) Starting and ending wavelength over which to calculate the statistics. Default is to use the full wavelength range.
  • response_func (array-like) – (Optional) A two-column array with the wavelength and transmission of a broad-band response function to use for the calculation.
  • per_pixel (bool) – (Optional) When providing a response function, continue to calculate the statistics per pixel, as opposed to per angstrom. Default is to compute the statistics on a per pixel basis.
  • flag (str or list) – (Optional) (List of) Flag names that are considered when deciding if a pixel should be masked. The names must be a valid bit name as defined by bitmask (see DRPFitsBitMask).
  • covar (bool) – (Optional) Flag to calculate covariance matrix.
  • correlation (bool) – (Optional) Flag to convert the covariance matrix to a correlation matrix on ouput.
  • covar_wave (double) – (Optional) Wavelength to use for the covariance calculation.
Returns:

Four objects are returned: the mean flux, the propagated variance in the mean flux, the mean S/N, and the covariance/correlation matrix for a single wavelength channel. If the object is an RSS file or no covariance calculation is requested, the last returned object is None.

Return type:

numpy.ndarray

Raises:

ValueError – Raised of a provided wavelength range object does not have two elements.

classmethod from_file(input_file, plate=0, ifudesign=0, read=False)[source]

Function to read a DRP-like data cube.

Warning

Work in progress. Currently does not work.

info()[source]

Print the HDU info page.

instrumental_dispersion_plane(channel, dispersion_factor=None, pixelscale=None, recenter=None, width_buffer=None, rlim=None, sigma=None, pre=False, quiet=False)[source]

Return the instrumental dispersion for the reconstructed ‘CUBE’ wavelength plane.

For ‘CUBE’ files, the input parameters must be the same as the defaults. If they are, the function must be able to find and access the ‘RSS’ file to construct the instrumental dispersion map because the necessary information is not in the ‘CUBE’ files (before MPL-6).

Todo

  • Need to implement something that will recognize that the ‘DISP’ extension exists in the MPL-6 and later DRP CUBE files.

For ‘RSS’ files, the transfer matrix, \({\mathbf T}\), is first calculated using regrid_transfer_matrix(). The transfer matrix is used to construct the ‘CUBE’ wavelength plane image, \({\mathbf I}\), by computing \({\mathbf T} \times {\mathbf F} = {\mathbf I}\), where \({\mathbf F}\) is the vector of the fiber fluxes. Under the assumption that the line-spread-function (LSF) is Gaussian, we determine the instrumental dispersion for the data in the wavelength channel of the reconstructed CUBE by calculating the second moment of the weighted sum of Gaussians of the appropriate dispersion. Assuming all the Gaussians have the normalized form:

\[\mathcal{N}(x|\mu=0,\sigma) = \frac{1}{\sqrt{2\pi}\sigma} \exp\left(-\frac{x^2}{2\sigma^2}\right),\]

the combined instrumental dispersion becomes

\[\sigma_{\rm inst,j}^2 = \frac{\sum_i T_{ij} \sigma^2_i}{\sum_i T_{ij}},\]

where \(T_{ij}\) are the elements of the transfer matrix. In terms of matrix multiplications, this can be written as

\[{\mathbf S} = \frac{ {\mathbf T} \times {\mathbf V} }{ {\mathbf T_c} },\]

where \({\mathbf T_c} = {\mathbf T_c} \times {\mathbf 1}\) is the vector with the sum \(\sum_j T_{ij}\), \({\mathbf V}\) is the instrumental variance for all fibers at the designated wavelength plane, and \({\mathbf S}\) is the variance for all the spaxels in the reconstructed wavelength image; the division by \({\mathbf T_c}\) is element-wise.

The returned matrix is the element-wise square-root of \({\mathbf S}\), rearranged into a 2D array of size ny by nx.

Warning

Because the internal data structure is transposed with respect to the DRP file, the result is transposed on output so that it will match the pixel coordinates in the reconstructed image pulled directly from a DRP file. See regrid_wavelength_channel().

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the transfer matrix.
  • dispersion_factor (float) – (Optional) Artificially multiply the dispersion measurements by this factor before calculating the reconstructed dispersion.
  • pixelscale (float) – (Optional) Desired pixel scale in arcsec
  • recenter (bool) – (Optional) Flag to recenter the coordinate system
  • width_buffer (int) – (Optional) Number of pixels to use as buffer for the image reconstruction
  • rlim (float) – (Optional) The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – (Optional) The sigma of the image reconstruction kernel in arcseconds.
  • pre (bool) – (Optional) Read the pre-pixelized version of the spectral resolution, instead of the post-pixelized version. This prepends ‘PRE’ to the extension name.
  • quiet (bool) – (Optional) Suppress terminal output
Returns:

The reconstructed image for the specified wavelength channel.

Return type:

numpy.ndarray

Raises:

ValueError – Raised for ‘CUBE’ files when the input parameters are not the same as the defaults.

mean_sky_coordinates(waverange=None, response_func=None, per_pixel=True, offset=True, flag=None, fluxwgt=False)[source]

Compute the mean sky coordinates for each spectrum.

For CUBE files, this just returns the spaxel coordinates in arcseconds relative to the object center, where the object center is \((\alpha_0,\delta_0) = ({\rm OBJRA},{\rm OBJDEC)\) as provided in the primary header, and

\[\begin{split}x &= (\alpha - \alpha_0) \cos \delta_0 \\ y &= (\delta - \delta_0)\end{split}\]

The coordinate grid, \((\alpha, \delta)\) is based on the WCS coordinates in the header as returned by world_mesh().

For RSS files, this returns either the unweighted or flux-weighted XPOS and YPOS values. For observations where the pointing center is different from the object center, the returned coordinates are relative to the object center if offset=True (default) and relative to the pointing center if offset=False.

Warning

Flux-weighting the coordinates can produce spurious results in low-flux regimes.

Parameters:
  • waverange (array-like) – (Optional) Two-element array with the first and last wavelength to include in the computation. Default is to use the full wavelength range.
  • response_func (array-like) – (Optional) A two-column array with the wavelength and transmission of a broad-band response function to use for the calculation.
  • per_pixel (bool) – (Optional) When providing a response function, continue to calculate the statistics per pixel, as opposed to per angstrom. Default is to
  • offset (bool) – Offset the coordinates to the object coordinates.
  • flag (str or list) – (Optional) (List of) Flag names that are considered when deciding if a pixel should be masked. The names must be a valid bit name as defined by bitmask (see DRPFitsBitMask).
  • fluxwgt (bool) – (Optional) Flag to weight by the flux when determining the mean coordinates for the RSS spectra.
Returns:

A 1D array with the coordinates of each spectrum. The index of the vector matches the index in spatial_index in case you want to recreate the a map for the CUBE files.

Return type:

numpy.ndarray

static mode_options()[source]

Return the allowed modes.

Returns:List of the allowed DRP fits file modes.
Return type:list
object_coo()[source]

Return the object right ascension (‘OBJRA’) and declination (‘OBJDEC’) taken from the DRP file primary header. The HDUList is opened if hasn’t been already.

Returns:Object RA and DEC in degrees.
Return type:float
object_data()[source]

Return the MaNGA ID (‘MANGAID’), object right ascension (‘OBJRA’), and object declination (‘OBJDEC’) taken from the DRP file primary header. The HDUList is opened if hasn’t been already.

Returns:MaNGA ID and object RA and DEC.
Return type:str, float, float
open_hdu(checksum=False)[source]

Read the fits file.

For security purposes, the fits file is always opened as read only. If hdu is not None, this function assumes the data has already been read and returns to the calling process.

Parameters:checksum (bool) – (Optional) Check for file corruption.
Raises:FileNotFoundError – Raised if the DRP file does not exist.
pix_mesh(pixelscale=None, recenter=None, width_buffer=None, extent=False)[source]

Uses numpy.meshgrid to create and return the I and J pixel coordinates for an nx*ny mesh.

For CUBE files, these dimensions are drawn directly from the WCS keywords in the header of the FLUX extension of the DRP fits file. In this case, any entered parameters are ignored and the class attributes are set to the default values used by the DRP.

Warning

The calculation for the CUBE files is only valid if the WCS coordinate system has no rotation.

For RSS files, the dimensions are determined using the data in the ‘XPOS’ and ‘YPOS’ extensions and the same algorithm used by the DRP; however, it is possible to provide different parameters that will alter the dimensions.

If extent is True, the returned mesh is of the pixel edges. That is, for a grid of size \(N\), there are \(N\) pixel centers that run from 1 to \(N\); however, there are \(N+1\) pixel edges that run from from \(0.5\) to \(N+0.5\).

See: pixelscale, recenter, width_buffer.

Parameters:
  • pixelscale (float) – (Optional) Desired pixel scale in arcsec
  • recenter (bool) – (Optional) Flag to recenter the coordinate system
  • width_buffer (int) – (Optional) Number of pixels to use as buffer for the image reconstruction
  • extent (bool) – (Optional) Return a grid of the pixel edges instead of the coordinates of the pixel centers.
Returns:

Two arrays with the pixel indices; numpy.meshgrid is run using:

indexing='ij'

Return type:

numpy.ndarray

pix_mesh_range(pixelscale=None, recenter=None, width_buffer=None)[source]

Return the range in x and y of the reconstructed image pixels, including the size of the pixel. Coordinate \((1,1)\) is the center of the first pixel, so its bottom corner is at \((0.5,0.5)\).

For CUBE files, these dimensions are drawn directly from the WCS keywords in the header of the FLUX extension of the DRP fits file. In this case, any entered parameters are ignored and the class attributes are set to the default values used by the DRP.

Warning

The calculation for the CUBE files is only valid if the WCS coordinate system has no rotation.

For RSS files, the dimensions are determined using the data in the ‘XPOS’ and ‘YPOS’ extensions and the same algorithm used by the DRP; however, it is possible to provide different parameters that will alter the dimensions.

See: pixelscale, recenter, width_buffer.

Parameters:
  • pixelscale (float) – (Optional) Desired pixel scale in arcsec
  • recenter (bool) – (Optional) Flag to recenter the coordinate system
  • width_buffer (int) – (Optional) Number of pixels to use as buffer for the image reconstruction
Returns:

Two arrays with, respectively, the lower and upper x range and the lower and upper y range.

Return type:

numpy.ndarray

pointing_offset()[source]

Return the offsets in RA and DEC between the pointing coordinates (IFURA, IFUDEC) and the designated object center coordinates (OBJRA, OBJDEC), drawn from the primary header of the DRP fits file.

Returns:Sky-right arcsecond offsets in RA and DEC.
Return type:float
regrid_transfer_matrix(channel, pixelscale=None, recenter=None, width_buffer=None, rlim=None, sigma=None, quiet=False, rej_flag='3DREJECT')[source]

Calculate the transfer matrix used to produce a reconstructed image of the fiber data at the specified wavelength channel. See: regrid_T.

For ‘RSS’ files, this is done directly using the available on-sky x and y coordinates of each fiber as a function of wavelength taken from the ‘XPOS’ and ‘YPOS’ extensions.

For ‘CUBE’ files, the function will attempt to use the ‘RSS’ counterpart of the file to produce transfer matrix. In this case, the input parameters must be the defaults.

The default behavior is to use the same parameters as used by the DRP. For the defaults, see mangadap.config.defaults.default_cube_pixelscale(), mangadap.config.defaults.default_cube_recenter(), mangadap.config.defaults.default_cube_width_buffer(), mangadap.config.defaults.default_regrid_rlim(), and mangadap.config.defaults.default_regrid_sigma(). However, the code allows the parameters to be freely chosen by the user.

See: pixelscale, recenter, width_buffer, regrid_rlim, regrid_sigma, regrid_T.

Todo

  • Give more detail on the pixels at which the radius is calculated.

  • The details of this calculation need to be checked against what is done by the DRP:

    • Are the positions in XPOS sky-right in the sense that positive XPOS is toward positive RA or not?
    • Does the DRP rebin the data using a kernel that defines the coordinate of a pixel at the pixel center or at the pixel edge?

Warning

The internal data structure of the RSS data is not different from simply reading the fits file using astropy.io.fits. This means you can use this function and apply the result to a separate read of the RSS file (not that you would ever want to do that). But after doing that, the output will be transposed wrt the CUBE file spatial orientation! See regrid_wavelength_plane().

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the transfer matrix.
  • pixelscale (float) – (Optional) Desired pixel scale in arcsec
  • recenter (bool) – (Optional) Flag to recenter the coordinate system
  • width_buffer (int) – (Optional) Number of pixels to use as buffer for the image reconstruction
  • rlim (float) – (Optional) The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – (Optional) The sigma of the image reconstruction kernel in arcseconds.
  • quiet (bool) – (Optional) Suppress terminal output
Returns:

Transfer matrix \({\mathbf T}\)

Return type:

scipy.sparse.csr_matrix

Raises:

ValueError – Raised for ‘CUBE’ files if the input parameters are not the defaults.

regrid_wavelength_plane(channel, pixelscale=None, recenter=None, width_buffer=None, rlim=None, sigma=None, quiet=False, return_ivar=False, return_covar=False)[source]

Return the reconstructed image for the specified wavelength channel.

For ‘CUBE’ files, the input parameters must be the same as the defaults. If they are, the function simply returns the selected 2D array from the ‘FLUX’ extension.

For ‘RSS’ files, the transfer matrix is first calculated using regrid_transfer_matrix() and then used to calculate:

\[{\mathbf T} \times {\mathbf F} = {\mathbf I}\]

where \({\mathbf F}\) is the vector of fluxes in the selected wavelength channel for all the fiber measurements in the field of view and \({\mathbf I}\) is the pre-formatted (flattened) reconstructed image of that wavelength channel.

On output, :math:`{mathbf I}` is rearranged into a 2D array of size :attr:`nx` by :attr:`ny`.

Warning

Because the internal data structure is transposed with respect to the DRP file, the result is transposed on output so that it will match the similar operation done directly on a DRP file.

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the transfer matrix.
  • pixelscale (float) – (Optional) Desired pixel scale in arcsec
  • recenter (bool) – (Optional) Flag to recenter the coordinate system
  • width_buffer (int) – (Optional) Number of pixels to use as buffer for the image reconstruction
  • rlim (float) – (Optional) The limiting radius of the image reconstruction kernel in arcseconds.
  • sigma (float) – (Optional) The sigma of the image reconstruction kernel in arcseconds.
  • quiet (bool) – (Optional) Suppress terminal output
  • return_ivar (bool) – (Optional) Return the nominal inverse variance image (i.e. the diagonal of the covariance matrix).
Returns:

The reconstructed image for the specified wavelength channel.

Return type:

numpy.ndarray

Raises:

Exception – Raised for ‘CUBE’ files when the input parameters are not the same as the defaults.

static sampling_options()[source]

Return the allowed wavelength sampling modes.

Returns:List of the allowed DRP fits wavelength sampling modes.
Return type:list
select(t, ext='FLUX', order='xy')[source]

Select a specific vector from the fits data arrays. The spatial position within the original file accessed by index i is given by spatial_index. For CUBE files this is a tuple with the relevant pixel indices from the original DRP file. This function is setup such that, e.g.:

from mangadap.drpfits import DRPFits
from astropy.io import fits
import numpy
drpf = DRPFits(7495, 1901, 'CUBE', read=True)
hdu = fits.open(drpf.file_path())

# should be no difference
flux_direct = hdu['FLUX'].data[:,15,16]
flux_class = drpf.select( (15,16), order='yx' )
assert not numpy.any( numpy.absolute(flux_direct - flux_class) > 0 ), \
    'Selection error!'

# Or to use more natural x y order
flux_class = drpf.select( (16,15) )
assert not numpy.any( numpy.absolute(flux_direct - flux_class) > 0 ), \
    'Selection error!'

# Or select the spectrum based on its index
for i,t in enumerate(drpf.spatial_index):
    if t == (15,16):
        break
flux_class = drpf.select(i)
assert not numpy.any( numpy.absolute(flux_direct - flux_class) > 0 ), \
    'Selection error!'

# Or use index and the provided spatial_index tuple
i = 550
t = drpf.spatial_index[i]
flux_direct = hdu['FLUX'].data[:,t[0],t[1]]
flux_class = drpf.select(i)
assert not numpy.any( numpy.absolute(flux_direct - flux_class) > 0 ), \
    'Selection error!'

flux_class = drpf.select(t)
assert not numpy.any( numpy.absolute(flux_direct - flux_class) > 0 ), \
    'Selection error!'
Parameters:
  • t (tuple or int) – If an integer, this is the flattened index of the spectrum to return. If a tuple, this is the position within the original DRP file for the non-spectral axes. See the examples above.
  • ext (str) – (Optional) Name of the extension from which to draw the data. Must be allowed for the current mode; see spectral_arrays. Default is 'FLUX'.
Returns:

Vector with the data values as a function of wavelength.

Return type:

numpy.ndarray

Raises:
  • KeyError – Raised if the selected extension cannot be used.
  • ValueError – Raised if the input selection object is not an int or a tuple.

Todo

Add select_near function that uses on-sky coordinates.

spectral_resolution(ext=None, toarray=False, fill=False, pre=False, median=False)[source]

Return the spectral resolution at each spatial and spectral position.

Select the extension ‘DISP’ or ‘SPECRES’. To get the pre-pixelized versions, set pre=True. If you set the extension to ‘PREDISP’ and pre=True, it will try to find the extension ‘PREPREDISP’ and fault.

Parameters:
  • ext (str, optional) – Specify the extension with the spectral estimate to use. Should be in [ None, ‘DISP’, ‘SPECRES’]. The default is None, which means it will return, in order of precedence, the data in ‘DISP’, ‘SPECRES’, or a None value if neither are present.
  • toarray (bool, optiional) – Return the spectral resolution as a 2D array: Nspec x Nwave, even if the DRP file is a CUBE object, similar to DRPFits.copy_to_array(). Default is to return an object with the same shape as the flux array.
  • fill (bool, optional) – Fill masked values by interpolation. Default is to leave masked pixels in returned array.
  • pre (bool, optional) – Read the pre-pixelized version of the spectral resolution, instead of the post-pixelized version. This prepends ‘PRE’ to the extension name.
  • median (bool, optional) – Return a single vector with the median spectral resolution instead of a per spectrum array. When using the SPECRES extension, this just returns the vector provided by the DRP file; when using the DISP extension, this performs a masked median across the array and then interpolates any wavelengths that were masked in all vectors.
Returns:

Even if interpolated such that there should be not masked values, the function returns a masked array. Array contains the spectral resolution (\(R = \lambda/\Delta\lambda\)) pulled from the DRP file.

Return type:

numpy.ma.MaskedArray

spectral_resolution_header(ext=None, pre=False)[source]

Return a fits header for the spectral resolution array. Copies the basic header from the relevant extension in the DRP file.

Parameters:
  • ext (str) – (Optional) Specify the extension with the spectral estimate to use. Should be in [ None, ‘DISP’, ‘SPECRES’]. The default is None, which means it will return, in order of precedence, the header for ‘DISP’, ‘SPECRES’, or an empty header if neither are present.
  • pre (bool) – (Optional) Read the pre-pixelized version of the spectral resolution, instead of the post-pixelized version. This prepends ‘PRE’ to the extension name.
world_mesh(skyright=True)[source]

Return the world X and Y coordinate for the nx by ny mesh of the reconstructed image. The pixel coordinates are first determined using pix_mesh() and then converted to the world coordinate system.

For ‘CUBE’ files, this is done using the WCS information in the header of the ‘FLUX’ extension. See wcs.

Note

Prior to v1_5_1 of the DRP, this function required a correction to the DRP header because of its use of ‘degrees’, which does not adhere to the fits standard causing the astropy.wcs.wcs.WCS to fail when initialized. This call is not longer made.

Todo

  • Make check of DRP version ‘VERSDRP’ to determine if header fix is necessary.

For an ‘RSS’ file, the returned mesh is based on the offset in arcseconds as determined by the dimensions of the reconstructed image. The skyright option will force the mesh to have a decreasing x coordinate as a function of increasing index. I.e.; the front edge of the first pixel is set to:

x0 = self.xs+self.nx*self.pixelscale if skyright else self.xs

and the offset per pixel is set to:

dx = -self.pixelscale if skyright else self.pixelscale

Recall that xs and ys are defined at the bottom-left edge of the first image pixel.

Parameters:skyright (bool) – If True, return the mesh with x coordinates decreasing as a function of increase index.
Returns:Two arrays with the world x and y coordinates.
Return type:numpy.ndarray
world_mesh_range(skyright=True)[source]

Return the range in the world X and Y coordinates including the size of the pixel.

For ‘RSS’ files, converts the result of pix_mesh_range() to world coordinates using the same calculation as in world_mesh().

For ‘CUBE’ files, converts the result of pix_mesh() with extent=True to world coordinates using the same calculation as in world_mesh(). The returned arrays give the range in X and Y respectively.

Warning

May not be accurate if the reconstructed image has an on-sky rotation.

Parameters:skyright (bool) – Return the range in X such that the order is [max(X), min(X)] if True; otherwise, values are returned as [min(X), max(X)].
Returns:Two arrays with, respectively, the range in X and Y.
Return type:numpy.ndarray

Note

Prior to v1_5_1 of the DRP, this function required a correction to the DRP header because of its use of ‘degrees’, which does not adhere to the fits standard causing the astropy.wcs.wcs.WCS to fail when initialized. This call is not longer made.

Todo

  • Make check of DRP version ‘VERSDRP’ to determine if header fix is necessary.
class mangadap.drpfits.DRPFitsBitMask(sdss_maskbits=None, mode='CUBE')[source]

Bases: mangadap.util.bitmask.BitMask

Structure with the DRP mask bits.

class mangadap.drpfits.DRPQuality3DBitMask(sdss_maskbits=None)[source]

Bases: mangadap.util.bitmask.BitMask

Structure with the definition of the DRP3QUAL mask bits.

mangadap.drpfits.drpfits_list(platelist, ifudesignlist, modelist, combinatorics=False, drpver=None, redux_path=None, directory_path=None)[source]

Provided a list of plates, ifudesigns, and modes, return a list of DRP files to be analyzed by the DAP.

If the number of elements in each list is the same, the matching is assumed to be finished unless combinatorics is True. If the number of elements is not the same, or cominatorics is True, the matched list is all the combinations of the provided elements.

Parameters:
  • platelist (str or list) – List of plates to use.
  • ifudesignlist (str or list) – List of IFU designs to use.
  • modelist (str or list) – List of DRP output modes (‘CUBE’ or ‘RSS’)
  • combinatorics (bool) – (Optional) Based on the input platelist and ifudesignlist, create a list with all possible combinations.
  • drpver (str) – (Optional) The DRP version, which must be a single string value used for all DRP files.
  • redux_path (str) – (Optional) The path to the top level directory containing the DRP output files; this is the same as the redux_path in the DRPFits class.
  • directory_path (str) – (Optional) The exact path to the DRP file. Default is defined by mangadap.config.defaults.default_drp_directory_path().
Returns:

A list of DRP file objects

Return type:

list

Raises:
  • Exception – Raised if drpver or redux_path are not strings.
  • ValueError – Raised if the platelist, ifudesignlist, or modelist are None.