mangadap.spectra.manga module¶
Implement the derived class for MaNGA row-stacked spectra.
Copyright © 2019, SDSS-IV/MaNGA Pipeline Group
- class mangadap.spectra.manga.MaNGARSS(ifile, sres_ext=None, sres_fill=True)[source]¶
Bases:
mangadap.util.drpfits.DRPFits
,mangadap.spectra.rowstackedspectra.RowStackedSpectra
Container class for MaNGA row-stacked spectra.
For additional description and attributes, see the two base classes.
See
from_plateifu()
to instantiate using the plate and ifu numbers. Seefrom_config()
to instantiate from a configuration file.- Parameters
ifile (
str
) – The file with the MaNGA row-stacked spectra. The name is expected to follow the nominal naming convention, which is used to parse the plate, ifudesign, and whether or not the spectra are binned logarithmically in wavelength.sres_ext (
str
, optional) – The extension to use when constructing the spectral resolution vectors. Seespectral_resolution()
.sres_fill (
bool
, optional) – Fill masked values by interpolation. Default is to leave masked pixels in returned array.
- static _parse_rectification_parameters(pixelscale, rlim, sigma, recenter, width_buffer)[source]¶
Parse the datacube rectification parameters.
Values that are None are set to their default methods in
mangadap.config.defaults
.- Returns
The default or input values of arguments, in the same order.
- Return type
tuple
- binned_on_sky_area(bin_indx, offset=True, **kwargs)[source]¶
Compute the on-sky area of a set of binned spectra.
This is a simple wrapper for
mangadap.spectra.rowstackedspectra.RowStackedSpectra.binned_on_sky_area()
andmean_sky_coordinates()
.The arguments are passed directly to
mean_sky_coordinates()
, and the results are then passed to the base class function. The area of the fibers “beams” are all renormalized to \(\pi {\rm arcsec}^2\) by the DRP.- Parameters
bin_indx (array-like) – An array with size \(N_{\rm spec}\) that gives which spaxels were included in each bin. Valid bins have indices of \(\geq 0\).
offset (
bool
, optional) – Offset the coordinates inxpos
andypos
by the difference between the IFU and OBJ coordinates in the header (seepointing_offset()
). The coordinates are nominally defined with respect to the IFU pointing center, such that this offsets the coordinates to be centered on the galaxy.**kwargs – Passed directly to
mean_sky_coordinates()
.
- Returns
Two numpy.ndarray objects are returned. The first has the unique (non-negative) bin indices, and the second provides the on-sky area of that bin.
- Return type
tuple
- static build_file_name(plate, ifudesign, log=True)[source]¶
Return the name of the DRP file with the row-stacked spectra.
This is a simple wrapper for
mangadap.util.drpfits.DRPFits.build_file_name()
, specific to the RSS.- Parameters
plate (
int
) – Plate numberifudesign (
int
) – IFU designlog (
bool
, optional) – Use the spectra that are logarithmically sampled in wavelength. If False, sampling is linear in wavelength.
- Returns
The relevant file name.
- Return type
str
- covariance_cube(channels=None, pixelscale=None, rlim=None, sigma=None, recenter=None, width_buffer=None, csr=False, quiet=False, rej_flag='3DREJECT')[source]¶
Construct the covariance matrices for many wavelength channels.
A simple wrapper for
mangadap.spectra.rowstackedspectra.RowStackedSpectra.covariance_cube()
that defaults to the nominal rectification parameters for MaNGA datacubes. See this base class function for the argument descriptions.
- covariance_matrix(channel, pixelscale=None, rlim=None, sigma=None, recenter=None, width_buffer=None, csr=False, quiet=False, rej_flag='3DREJECT')[source]¶
Construct the covariance matrix for a given wavelength channel.
A simple wrapper for
mangadap.spectra.rowstackedspectra.RowStackedSpectra.covariance_matrix()
that defaults to the nominal rectification parameters for MaNGA datacubes. See this base class function for the argument descriptions.
- static default_paths(plate, ifudesign, log=True, drpver=None, redux_path=None, directory_path=None)[source]¶
Construct the default path and file name with the MaNGA row-stacked spectra.
This is a simple wrapper for
mangadap.util.drpfits.DRPFits.default_paths()
, specific to the RSS files.- Parameters
plate (
int
) – Plate numberifudesign (
int
) – IFU designlog (
bool
, optional) – Use the spectra that are logarithmically sampled in wavelength. If False, sampling is linear in wavelength.drpver (
str
, optional) – DRP version, which is used to define the default DRP redux path. Default is defined bymangadap.config.defaults.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 bymangadap.config.defaults.drp_redux_path()
.directory_path (
str
, optional) – The exact path to the DRP file. Default is defined bymangadap.config.defaults.drp_directory_path()
. Providing this ignores anything provided fordrpver
orredux_path
.
- Returns
Two strings with the default path to and name of the DRP data file.
- Return type
tuple
- classmethod from_plateifu(plate, ifudesign, log=True, drpver=None, redux_path=None, directory_path=None, **kwargs)[source]¶
Construct a
mangadap.datacube.manga.MaNGARSS
object based on its plate-ifu designation.This is a simple wrapper function that calls
default_paths()
to construct the file names, and then calls the class instantiation method.- Parameters
plate (
int
) – Plate numberifudesign (
int
) – IFU designlog (
bool
, optional) – Use the spectra that are logarithmically binned in wavelength.drpver (
str
, optional) – DRP version, which is used to define the default DRP redux path. Default is defined bymangadap.config.defaults.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 bymangadap.config.defaults.drp_redux_path()
.directory_path (
str
, optional) – The exact path to the DRP file. Default is defined bymangadap.config.defaults.drp_directory_path()
. Providing this ignores anything provided fordrpver
orredux_path
.**kwargs – Keyword arguments passed directly to the primary instantiation method; see
mangadap.spectra.manga.MaNGARSS
.
- instrumental_dispersion_plane(channel, dispersion_factor=None, pixelscale=None, rlim=None, sigma=None, recenter=None, width_buffer=None, quiet=False, rej_flag='3DREJECT')[source]¶
Construct the instrumental dispersion for a given wavelength channel.
A simple wrapper for
mangadap.spectra.rowstackedspectra.RowStackedSpectra.instrumental_dispersion_plane()
that defaults to the nominal rectification parameters for MaNGA datacubes. See this base class function for the argument descriptions.
- mean_sky_coordinates(offset=True, **kwargs)[source]¶
Compute the weighted or unweighted mean sky coordinates for each spectrum.
This is a simple wrapper for
mangadap.spectra.rowstackedspectra.RowStackedSpectra.mean_sky_coordinates()
; see that method for additional keyword arguments.Warning
Flux-weighting the coordinates can produce spurious results in low-flux regimes.
- Parameters
offset (
bool
, optional) – Offset the coordinates inxpos
andypos
by the difference between the IFU and OBJ coordinates in the header (seepointing_offset()
). The coordinates are nominally defined with respect to the IFU pointing center, such that this offsets the coordinates to be centered on the galaxy.**kwargs – Passed directly to
mangadap.spectra.rowstackedspectra.RowStackedSpectra.mean_sky_coordinates()
;
- Returns
Two numpy.ndarray objects with the fiducial x and y on-sky positions of each spectrum in arcseconds relative to a given center. The shape of the two arrays is \((N_{\rm spec},)\).
- Return type
tuple
- 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
- rectification_transfer_matrix(channel, pixelscale=None, rlim=None, sigma=None, recenter=None, width_buffer=None, quiet=False, rej_flag='3DREJECT')[source]¶
Construct the rectification transfer matrix.
A simple wrapper for
mangadap.spectra.rowstackedspectra.RowStackedSpectra.rectification_transfer_matrix()
that defaults to the nominal rectification parameters for MaNGA datacubes. See this base class function for the argument descriptions.
- rectify_wavelength_plane(channel, pixelscale=None, rlim=None, sigma=None, recenter=None, width_buffer=None, quiet=False, rej_flag='3DREJECT', return_ivar=False, return_covar=False)[source]¶
Rectify a wavelength channel into a monochromatic flux map.
A simple wrapper for
mangadap.spectra.rowstackedspectra.RowStackedSpectra.rectify_wavelength_plane()
that defaults to the nominal rectification parameters for MaNGA datacubes. See this base class function for the argument descriptions.