mangadap.datacube.manga module¶
Implement the derived class for MaNGA datacubes.
Copyright © 2019, SDSS-IV/MaNGA Pipeline Group
- class mangadap.datacube.manga.MaNGADataCube(ifile, z=None, vdisp=None, ell=None, pa=None, reff=None, sres_ext=None, sres_fill=True, covar_ext=None)[source]¶
Bases:
mangadap.util.drpfits.DRPFits
,mangadap.datacube.datacube.DataCube
Container class for a MaNGA datacube.
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.Note that
z
,vdisp
,ell
,pa
, andreff
are saved to themeta
dictionary.- Parameters
ifile (
str
) – The file with the MaNGA datacube. The name is expected to follow the nominal naming convention, which is used to parse the plate, ifudesign, and whether or not the datacube has been binned logarithmically in wavelength.z (
float
, optional) – Estimated bulk redshift. If None, some of the DAP analysis modules will fault.vdisp (
float
, optional) – Estimated velocity dispersion. If None, some of the DAP analysis modules will assume an initial guess of 100 km/s.ell (
float
, optional) – Characteristic isophotal ellipticity (1-b/a). If None, some of the DAP modules will issue a warning and continue by assumingell=0
.pa (
float
, optional) – Characteristic isophotal position angle (through E from N). If None, some of the DAP modules will issue a warning and continue by assumingpa=0
.reff (
float
, optional) – Effective (half-light) radius in arcsec. If None, some of the DAP modules will issue a warning and continue by assumingreff=1
.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.covar_ext (
str
, optional) – Extension to use as the single spatial correlation matrix for all wavelength channels, read from the DRP file. For generating the covariance matrix directly for an arbitrary wavelength channel using the RSS file, seemangadap.datacube.datacube.DataCube.covariance_matrix()
.
- approximate_correlation_matrix(sigma_rho=1.92, rlim=None, rho_tol=None, redo=False)[source]¶
Construct an approximate correlation matrix for the datacube.
This is a simple wrapper for
mangadap.datacube.datacube.DataCube.approximate_correlation_matrix()
that defaults to the expected values for a MaNGA datacube. See that method for a description of the arguments. The defaultsigma_rho
is from Westfall et al. (2019, AJ, 158, 231). The defaultrlim
is provided bymangadap.config.defaults.regrid_rlim()
.
- approximate_covariance_cube(channels=None, sigma_rho=1.92, rlim=None, rho_tol=None, csr=False, quiet=False)[source]¶
Construct approximate covariance matrices for multiple channels.
This is a simple wrapper for
mangadap.datacube.datacube.DataCube.approximate_covariance_cube()
that defaults to the expected values for a MaNGA datacube. See that method for a description of the arguments. The defaultsigma_rho
is from Westfall et al. (2019, AJ, 158, 231). The defaultrlim
is provided bymangadap.config.defaults.regrid_rlim()
.
- approximate_covariance_matrix(channel, sigma_rho=1.92, rlim=None, rho_tol=None, csr=False, quiet=False)[source]¶
Construct an approximate covariance matrix.
This is a simple wrapper for
mangadap.datacube.datacube.DataCube.approximate_covariance_matrix()
that defaults to the expected values for a MaNGA datacube. See that method for a description of the arguments. The defaultsigma_rho
is from Westfall et al. (2019, AJ, 158, 231). The defaultrlim
is provided bymangadap.config.defaults.regrid_rlim()
.
- static build_file_name(plate, ifudesign, log=True)[source]¶
Return the name of the DRP datacube file.
This is a simple wrapper for
mangadap.util.drpfits.DRPFits.build_file_name()
, specific to the datacube.- 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
- 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 datacube.
This is a simple wrapper for
mangadap.util.drpfits.DRPFits.default_paths()
, specific to the datacube 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.MaNGADataCube
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 datacube that is 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.datacube.manga.MaNGADataCube
ormangadap.spectra.manga.MaNGARSS
.
- load_rss(force=False)[source]¶
Try to load the source row-stacked spectra for this datacube.
If
rss
is not None, this method does not attempt to reload it, unlessforce
is True.The method first looks for the relevant file in the same directory with the datacube. If the file is not there, it uses the DRP version in the header of the datacube file and the paths defined by the keyword arguments to construct the nominal path to the RSS file. If the file is also not there, the method faults.
Nothing is returned. If successful, the method initializes the row-stacked spectra object (
mangadap.spectra.manga.MaNGARSS
) torss
.- Parameters
force (
bool
, optional) – Reload the row-stacked spectra ifrss
is not None.
- mean_sky_coordinates(center_coo=None, offset='OBJ')[source]¶
Calculate the sky coordinates for each spaxel.
This is a simple wrapper for
mangadap.datacube.datacube.DataCube.mean_sky_coordinates()
that passes the relevant set of center coordinates; see that function for the default behavior. Any provided coordinates (seecenter_coo
) take precedence.- Parameters
center_coo (
tuple
, optional) – A two-tuple with the coordinates in right-ascension and declination for the coordinate-frame origin. If None, no offset is performed.offset (
str
, optional) – The offset coordinates to use from the primary header. Must be either ‘obj’ or ‘ifu’. If None, no header coordinates are used.
- Returns
Two numpy.ndarray objects with the RA and declination of each pixel in degrees, or its offset from the center in arcseconds. In both cases the shape of the returned arrays matches the spatial shape of the datacube.
- Return type
tuple