mangadap.dapfits module

Defines a class used to interface with the final maps files produced by the MaNGA Data Analysis Pipeline (DAP).

License:
Copyright (c) 2015, SDSS-IV/MaNGA Pipeline Group
Licensed under BSD 3-clause license - see LICENSE.rst
Class usage examples:
Add some usage comments here!
Revision history:
08 Jun 2016: Original Implementation by K. Westfall (KBW)
Feb 2017: (KBW) Merged construct maps and cube functions here and removed dapmaps.py and dapcube.py.
** 3 May 2017**: (KBW) Allow output maps and cube types to be single-precision (float32) floats, and changed BINID from int (defaults to 64-bit) to int32.
24 Aug 2017: (KBW) Added BADGEOM DAPQUAL bit; signifies invalid photometric geometry parameters used when instantiating the mangadap.par.obsinput.ObsInputPar object.
29 Aug 2017: (KBW) Add S/N metrics to the MAPS primary header
30 Aug 2017: (KBW) Convert some changes from Xihan Ji.
19 Mar 2018: (KBW) Mask spectral indices with incomplete band coverage as unreliable.

Todo

  • Allow DAPFits to read/access both the MAPS and the LOGCUBE files, not just the former.
class mangadap.dapfits.DAPCubeBitMask(dapsrc=None)[source]

Bases: mangadap.util.bitmask.BitMask

Todo

  • Force read IDLUTILS version as opposed to internal one?
class mangadap.dapfits.DAPFits(plate, ifudesign, method, drpver=None, dapver=None, analysis_path=None, directory_path=None, reference_path=None, par_file=None, read=True, checksum=False)[source]

Bases: object

Object used to hold properties of and read data from a DAP-produced file.

Parameters:
plate

Plate number

Type:int
ifudesign

IFU design

Type:int
method

Method type

Type:str
drpver

DRP version

Type:str
dapver

DAP version

Type:str
analysis_path

The path to the top level directory containing the DAP output files for a given DRP and DAP version.

Type:str
directory_path

The exact path to the DAP file.

Type:str
reference_path

The exact path of the DAP reference directory.

Type:str
par_file

SDSS parameter file used to provide input parameters for the DAP.

Type:str
dapqual

Global quality bit

Type:mangadap.dapmaps.DAPQualityBitMask
bitmask

Map mask bits

Type:mangadap.dapmaps.DAPMapsBitMask
hdu

HDUList read from the DRP file

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

List of parameters used by the DAP.

Type:mangadap.par.ObsInputPar
checksum

Boolean to use for checksum in astropy.io.fits.open.

Type:bool
channel_map(ext, channel=None, flag=None)[source]

Return a 2D array with the map for the specified channel.

effective_radius()[source]

Return the effective radius from the parameter file.

file_name()[source]

Return the name of the DAP file

file_path()[source]

Return the full path to the DAP file

guess_cz()[source]

Return the guess redshift (cz) from the parameter file.

guess_inclination(q0=None)[source]

Return a guess inclination based on the isophotal ellipticity from the input parameter file using the equation:

\[\cos^2 i = \frac{ q^2 - q_0^2 }{ 1 - q_0^2 },\]

where \(q\) is the observed oblateness (the semi-minor to semi-major axis ratio, \(q = 1-\epsilon = b/a\)) and \(q_0\) is the intrinsic oblateness.

If \(q < q_0\), the function returns a 90 degree inclination.

Parameters:q0 (float) – (Optional) The intrinsic oblateness of the system. If not provided, the system is assumed to be infinitely thin.
Raises:Exception – Raised if the input intrinsic oblatenss is not less than one and greater than 0.
guess_position_angle(flip=False)[source]

Return the guess position angle from the parameter file.

Parameters:flip (bool) – (Optional) Offset the provided position angle by 180 deg.
list_channels(ext)[source]

Provide a list of the channels in a given extension.

open_hdu(permissions='readonly', checksum=False, quiet=True)[source]

Open the fits file and save it to hdu; if hdu is not None, the function returns without re-reading the data.

Parameters:
  • permissions (string) – (Optional) Open the fits file with these read/write permissions.
  • checksum (bool) – (Optional) Flag for astropy.io.fits.open checksum operation. This overrides the internal checksum attribute for the current operation only.
  • quiet (bool) – (Optional) Suppress terminal output
Raises:

FileNotFoundError – Raised if the DAP file doesn’t exist.

read_par(quiet=True)[source]

Open the parameter file and save it to par; if par is not None, the function returns without re-reading the data.

Parameters:quiet (bool) – Suppress terminal output
thin_disk_polar_coo(xc=None, yc=None, rot=None, pa=None, inc=None, flip=False)[source]

Calculate the disk-plane coordinates for all the binned spectra using mangadap.util.geometry.SemiMajorAxisCoo. See the documentation their for further information regarding the calculations.

If not provided, the default position angle and inclination will be taken from the input parameter file, if available. If the parameter file cannot be read, the default values for the position angle and ellipticity in mangadap.util.geometry.SemiMajorAxisCoo are used.

Parameters:
  • xc,yc (float,float) – (Optional) a reference on-sky position relative to the center of the projected plane (galaxy center).
  • rot (float) – (Optional) Cartesian rotation of the focal-plane relative to the sky-plane (+x toward East; +y toward North).
  • pa (float) – (Optional) On-sky position angle of the major axis of the ellipse traced on the sky by a circle in the projected plane, defined as the angle from North through East.
  • inc (float) – (Optional) Inclination of the plane, angle between the disk normal and the normal of the sky plane (inc=0 is a face-on plane).
  • flip (bool) – (Optional) Offset the provided position angle by 180 deg; e.g., to set the position angle will be defined along the receding major axis.
Returns:

Two numpy arrays with the projected polar

coordinates: \(R, \theta\).

Return type:

numpy.ndarray

unique_indices()[source]

Use the BINID extension to select the spaxels with unique data.

Returns:Returns the list of indices in the flattened array of a map that are unique.
Return type:numpy.ndarray
unique_mask()[source]

Construct a boolean mask for the unique bins.

The returned map is True for spaxels that are part of a unique bin, Fals otherwise.

class mangadap.dapfits.DAPMapsBitMask(dapsrc=None)[source]

Bases: mangadap.util.bitmask.BitMask

Todo

  • Force read IDLUTILS version as opposed to internal one?
class mangadap.dapfits.DAPQualityBitMask(dapsrc=None)[source]

Bases: mangadap.util.bitmask.BitMask

Todo

  • Force read IDLUTILS version as opposed to internal one?
mangadap.dapfits.add_snr_metrics_to_header(hdr, drpf, r_re, dapsrc=None)[source]

For all valid spaxels within 1 Re < R < 1.5 Re calculate the median S/N and combined S/N (including covariance) in the griz bands.

Adds the following keywords (8) to the header:
  • SNR?MED: Median S/N within 1-1.5 Re in each of the griz bands.
  • SNR?RING: The combined S/N of all spaxels within 1-1.5 Re in the griz bands, incuding covariance

The inclusion of covariance is based on the exact calculation of the correlation matrix at the flux-weighted center of each band, and assuming no change in the correlation between the response-weighted correlation within the broad band.

Parameters:
  • hdr (astropy.io.fits.Header) – Header object to edit
  • drpf (mangadap.drpfits.DRPFits) – DRP fits cube
  • r_re (numpy.ndarray) – Flattened array with the semi-major axis radius in units of the effective radius for all spaxels in the datacube. Shape should be (Nx*Ny,), where the shape of the datacube is (Nx,Ny,Nwave).
Returns:

The edited header object.

Return type:

astropy.io.fits.Header

Raises:

FileNotFoundError – Raised if any of the response function files cannot be found.

mangadap.dapfits.combine_binid_extensions(drpf, binned_spectra, stellar_continuum, emission_line_moments, emission_line_model, spectral_indices, dtype=None)[source]

Combine the bin IDs from the different analysis steps into a single multichannel extension.

mangadap.dapfits.confirm_dap_types(drpf, obs, rdxqa, binned_spectra, stellar_continuum, emission_line_moments, emission_line_model, spectral_indices)[source]
class mangadap.dapfits.construct_cube_file(drpf, obs=None, binned_spectra=None, stellar_continuum=None, emission_line_model=None, dapsrc=None, dapver=None, analysis_path=None, directory_path=None, output_file=None, clobber=True, loggers=None, quiet=False, single_precision=False)[source]

Bases: object

Construct a DAP cube file based on the input data.

Set as a class for coherency reasons, but should not be used as an object!

Should force all intermediate objects to be provided.

_get_data_mask(binned_spectra, binned_spectra_3d_hdu)[source]
For the binned spectra:
  • consolidate ANY flag except NO_STDDEV from binned_spectra into IGNORED; done use SpatiallyBinnedSpectra.do_not_fit_flags
  • propagate NONE_IN_STACK from binned_spectra into FLUXINVALID
  • propagate IVARINVALID and FORESTAR
_get_emission_line_model_mask(emission_line_model, emission_line_model_3d_hdu)[source]
For the emission-line models:
  • propagate FORESTAR
  • propagate ARTIFACT
  • consolidate DIDNOTUSE, FORESTAR, LOW_SNR, ARTIFACT, OUTSIDE_RANGE, TPL_PIXELS, TRUNCATED, PPXF_REJECT, INVALID_ERROR into FITIGNORED
_get_stellar_continuum_mask(stellar_continuum, stellar_continuum_3d_hdu)[source]
For the stellar continuum models:
  • propagate FORESTAR
  • propagate INVALID_ERROR into IVARINVALID
  • propagate ARTIFACTs
  • consolidate DIDNOTUSE, FORESTAR, LOW_SNR, ARTIFACT, OUTSIDE_RANGE, EML_REGION, TPL_PIXELS, TRUNCATED, PPXF_REJECT, INVALID_ERROR, NO_FIT into FITIGNORED
  • consolidate FIT_FAILED into FITFAILED
_include_no_model_mask(mask)[source]

Assign the NOMODEL bit to spectral regions outside of the fitted spectral range.

Modeled off of StellarContiuumModel.reset_continuum_mask_window.

_set_no_model(mask)[source]
_set_paths(directory_path, dapver, analysis_path, output_file, binned_spectra, stellar_continuum, emission_line_model)[source]

Set the paths relevant to the map file.

binned_data_cube(prihdr, binned_spectra, binned_spectra_3d_hdu)[source]

Constructs the ‘FLUX’, ‘IVAR’, ‘MASK’, ‘WAVE’, and ‘REDCORR’ model cube extensions, and begins the construction of the ‘MASK’.

Returns the primary header, a list of ImageHDUs, and the mask array.

model_cubes(prihdr, binned_spectra, stellar_continuum, stellar_continuum_3d_hdu, emission_line_model, emission_line_model_3d_hdu)[source]

Constructs the ‘MODEL’, ‘MODEL_MASK’, ‘EMLINE’, ‘STELLAR’, and ‘STELLAR_MASK’ model cube extensions, adds the model information to the header, and appends data to the mask.

class mangadap.dapfits.construct_maps_file(drpf, obs=None, rdxqa=None, binned_spectra=None, stellar_continuum=None, emission_line_moments=None, emission_line_model=None, spectral_indices=None, nsa_redshift=None, dapsrc=None, dapver=None, analysis_path=None, directory_path=None, output_file=None, clobber=True, loggers=None, quiet=False, single_precision=False)[source]

Bases: object

Construct a DAP MAPS file.

Set as a class for coherency reasons, but should not be used as an object!

Should force all intermediate objects to be provided.

_build_binning_mask(binned_spectra)[source]

Consolidate the map mask values into NOVALUE.

_consolidate_donotuse(mask)[source]
_emission_line_model_mask_to_map_mask(emission_line_model, elf_mask, for_dispersion=False)[source]

Propagate and consolidate the emission-line moment masks to the map pixel mask.

INSUFFICIENT_DATA propagated to NOVALUE

FIT_FAILED, UNDEFINED_COVAR propagated to FITFAILED

NEAR_BOUND copied to NEARBOUND

UNDEFINED_SIGMA propagated to UNRELIABLE

_emission_line_moment_mask_to_map_mask(emission_line_moments, elm_mask)[source]

Propagate and consolidate the emission-line moment masks to the map pixel mask.

DIDNOTUSE, FORESTAR propagated from already existing mask (self.bin_mask)

MAIN_EMPTY, BLUE_EMPTY, RED_EMPTY, UNDEFINED_BANDS propagated to NOVALUE

MAIN_JUMP, BLUE_JUMP, RED_JUMP, JUMP_BTWN_SIDEBANDS propagated to FITFAILED

MAIN_INCOMP, BLUE_INCOMP, RED_INCOMP propagated to UNRELIABLE

NO_ABSORPTION_CORRECTION propagated to NOCORRECTION

DIVBYZERO propagated to MATHERROR

Second moments not provided so no need to include UNDEFINED_MOM2

static _join_maps(a)[source]
_set_paths(directory_path, dapver, analysis_path, output_file, binned_spectra, stellar_continuum, emission_line_model)[source]

Set the paths relevant to the map file.

_spectral_index_mask_to_map_mask(spectral_indices, si_mask)[source]

Propagate and consolidate the spectral index masks to the map pixel mask.

DIDNOTUSE, FORESTAR propagated from already existing mask (self.bin_mask)

MAIN_EMPTY, BLUE_EMPTY, RED_EMPTY, UNDEFINED_BANDS propagated to NOVALUE

MAIN_INCOMP, BLUE_INCOMP, RED_INCOMP propagated to UNRELIABLE

NO_DISPERSION_CORRECTION propagated to NOCORRECTION

DIVBYZERO propagated to MATHERROR

Need to further assess *_INCOMP to see if these lead to bad values (BADVALUE).

_stellar_continuum_mask_to_map_mask(stellar_continuum, sc_mask, for_dispersion=False)[source]

sc_mask constains the reconstructed map of the masks in each stellar continuum fit (using the stellar continuum bitmask).

Propagate and consolidate the stellar-continuum masks to the map pixel mask.

DIDNOTUSE, FORESTAR propagated from already existing mask
(self.bin_mask)

LOW_SNR, NO_FIT, INSUFFICIENT_DATA propagated to NOVALUE

FIT_FAILED, NEAR_BOUND propagated to FITFAILED and NEARBOUND

NEGATIVE_WEIGHTS consolidated into UNRELIABLE; if
dispersion=True, include BAD_SIGMA in this
binned_spectra_maps(prihdr, obs, binned_spectra)[source]

Constructs the ‘BIN_LWSKYCOO’, ‘BIN_LWELLCOO’, ‘BIN_AREA’, ‘BIN_FAREA’, ‘BIN_MFLUX’, ‘BIN_MFLUX_IVAR’, ‘BIN_MFLUX_MASK’, and ‘BIN_SNR’ map extensions.

emission_line_model_maps(prihdr, emission_line_model)[source]

Construct the ‘EMLINE_GFLUX’, ‘EMLINE_GFLUX_IVAR’, ‘EMLINE_GFLUX_MASK’, ‘EMLINE_GEW’, ‘EMLINE_GEW_IVAR’, ‘EMLINE_GEW_MASK’, ‘EMLINE_GVEL’, ‘EMLINE_GVEL_IVAR’, ‘EMLINE_GVEL_MASK’, ‘EMLINE_GSIGMA’, ‘EMLINE_GSIGMA_IVAR’, ‘EMLINE_GSIGMA_MASK’, ‘EMLINE_INSTSIGMA’, ‘EMLINE_TPLSIGMA’, ‘EMLINE_GA’, ‘EMLINE_GANR’, ‘EMLINE_FOM’, ‘EMLINE_LFOM’ map extensions.

emission_line_moment_maps(prihdr, emission_line_moments)[source]

Construct the ‘EMLINE_SFLUX’, ‘EMLINE_SFLUX_IVAR’, ‘EMLINE_SFLUX_MASK’, ‘EMLINE_SEW’, ‘EMLINE_SEW_IVAR’, and ‘EMLINE_SEW_MASK’ maps extensions.

reduction_assessment_maps(prihdr, obs, rdxqa)[source]

Constructs the ‘SPX_SKYCOO’, ‘SPX_ELLCOO’, ‘SPX_MFLUX’, ‘SPX_MFLUX_IVAR’, and ‘SPX_SNR’ map extensions.

Todo

Add the spaxel correlation data.

spectral_index_maps(prihdr, spectral_indices)[source]

Construct the ‘SPECINDEX’, ‘SPECINDEX_IVAR’, ‘SPECINDEX_MASK’, and ‘SPECINDEX_CORR’.

stellar_continuum_maps(prihdr, stellar_continuum)[source]

Construct the ‘STELLAR_VEL’, ‘STELLAR_VEL_IVAR’, ‘STELLAR_VEL_MASK’, ‘STELLAR_SIGMA’, ‘STELLAR_SIGMA_IVAR’, ‘STELLAR_SIGMA_MASK’, ‘STELLAR_SIGMACORR’, ‘STELLAR_FOM’ maps extensions.

mangadap.dapfits.finalize_dap_primary_header(prihdr, drpf, obs, binned_spectra, stellar_continuum, dapsrc=None, loggers=None, quiet=False)[source]