mangadap.datacube.datacube module

Base class for a datacube


License

Copyright © 2019, SDSS-IV/MaNGA Pipeline Group


class mangadap.datacube.datacube.DataCube(flux, wave=None, ivar=None, mask=None, bitmask=None, sres=None, covar=None, axes=[0, 1, 2], wcs=None, pixelscale=None, log=True, meta=None, prihdr=None, fluxhdr=None, name=None)[source]

Bases: object

Base container class for a rectilinear datacube.

Datacubes have three axes: two spatial coordinates and one spectral. The datacubes are expected to be rectilinear; i.e., each spatial position has the same spectral range and each wavelength channel has the same spatial coordinates.

On input, the ordering of the three dimensions can be arbitrary; however, the axes are re-ordered for the internal attributes such that wavelengths are ordered along the last axis, with the first two axes being the spatial coordinates. Nominally the spatial coordinates are ordered predominantly coincident with right-ascension (E toward smaller pixels numbers in axis 0) and declination (N toward larger pixel numbers in axis 1); however, rotation of the datacube spatial coordinates can be non-zero with respect to the celestial coordinates, as given by the provided world-coordinate system.

The wavelength vector applicable to all spatial positions can either be provided directly or constructed from the provided WCS. Any directly provided wave vector takes precedence over the WCS.

Spatial covariance/correlation can be provided via the covar keyword argument. If provided as a single array, the provided array is used to define the correlation matrix and is assumed to be identical for all wavelength channels. More options are available if the input is provided as a Covariance object. The ordering of the covariance/correlation matrix is expected to provide the correlation between the row-major flattened version of the spatial coordinates. That is, for a datacube with spatial shape \((N_x,N_y)\), the covariance/correlation value at location \(i,j\) is the correlation between pixels at 2D locations \((i_x,i_y)\) and \((j_x,j_y)\), where

\[\begin{split}\begin{array}{rcl} i_x & = & \lfloor i/N_y \rfloor, \\ i_y & = & i - i_x N_y, \\ j_x & = & \lfloor j/N_y \rfloor, {\rm and} \\ j_y & = & j - j_x N_y. \end{array}\end{split}\]

If the values of axes indicates that the provided flux array should have it’s spatial axes transposed (i.e. axes[0] > axes[1]), the provided covariance/correlation data is appropriately restructured to correspond to the wavelength channels in flux (again assuming a flattened row-major memory block); see transpose_raw_shape().

See Creating a new DataCube subclass for instructions on subclassing this object for non-MaNGA datacubes.

Parameters:
  • flux (numpy.ndarray) – Data array with the flux as a function of spatial and spectral position. Shape should be \((N_x, N_y, N_\lambda)\).

  • wave (numpy.ndarray, optional) – The wavelength vector associated with each spatial position. Shape is \((N_\lambda,)\). If None, wcs must be provided such that the wavelength vector can be constructed.

  • ivar (numpy.ndarray, optional) – The inverse variance of the flux array. Shape must be the same as flux. If None, errors are ignored.

  • mask (numpy.ndarray, optional) – Mask array for the flux measurements. Can be a boolean mask (False=unmasked,good; True=masked,bad) or an integer array associated with the provided mangadap.util.bitmask.BitMask object (see bitmask). Shape must be the same as flux.

  • bitmask (mangadap.util.bitmask.BitMask, optional) – Object used to select and toggle masked pixels from mask. If None, any provided mask must be a boolean array or can be successfully converted to one.

  • sres (numpy.ndarray, optional) – The spectral resolution, \(R = \lambda / \Delta\lambda\). Can be a single vector with shape \((N_\lambda,)\), if the spectral resolution is independent of spatial position, or a 3D datacube with a spatially and spectrally dependent resolution with a shape that matches flux. If None, spectral resolution is ignored.

  • covar (array-like, mangadap.util.covariance.Covariance, optional) – The spatial covariance in one or more wavelength channels. If None, no covariance is assumed. See usage and interpretation in class description.

  • axes (list, optional) – The axes with the \(x\), \(y\), and \(\lambda\) in the provided flux array. For example, [2,1,0] means the spectra are organized along the first axis.

  • wcs (astropy.wcs.WCS, optional) – World-coordinate system for the 3D datacube. If None, spatial coordinates are set to be centered with a pixelscale in arcsec (see pixelscale).

  • pixelscale (int, float, optional) – Spatial extent in arcsec of each spaxel in the datacube, assumed to be square. Superseded by wcs, if the latter is provided. If the wcs is not provided and the pixelscale is not provided upon instantiation, it is set to unity.

  • log (bool, optional) – Flag that the datacube spectral pixels are binned logarithmically in wavelength.

  • meta (dict, optional) – A free-form dictionary used to hold metadata relevant to the datacube. Metadata required by analysis modules are indicated where relevant. If None, meta is instantiated as an empty dictionary and populated with the bare minimum needed to execute the DAP. NOTE: Currently this is a required argument because an estimate of the redshift is needed to execute the DAP.

  • prihdr (astropy.io.fits.Header, optional) – Primary header read from datacube fits file. If None, instantiated as an empty astropy.io.fits.Header.

  • fluxhdr (astropy.io.fits.Header, optional) – Header specifically for the flux extension of the datacube fits file. If None, set to be a copy of the primary header.

Raises:
  • ValueError – Raised if the wavelength vector cannot be produced, if the WCS has the wrong dimensionality, or if there are any shape mismatches between the input data arrays.

  • TypeError – Raised if the input metadata are not provided as a dictionary.

original_axes

Original axis order.

Type:

numpy.ndarray

shape

Datacube shape.

Type:

tuple

spatial_shape

Shape of the datacube spatial axes .

Type:

tuple

nwave

Number of wavelength channels.

Type:

int

spatial_index

Array of tuples with the spatial indices of each spectrum in the flattened datacube.

Type:

numpy.ndarray

wcs

Datacube world-coordinate system. Can be None.

Type:

astropy.wcs.WCS

pixelscale

Spatial extent in arcsec of each spaxel in the datacube, assumed to be square. Only used if wcs is None.

Type:

float

bitmask

Object used to select and toggle masked pixels from mask. Can be None.

Type:

mangadap.util.bitmask.BitMask

log

Flag that the datacube spectral pixels are binned logarithmically in wavelength.

Type:

bool

meta

A free-form dictionary used to hold meta data relevant to the datacube. Metadata required by analysis modules are indicated where relevant. If no metadata has been defined, meta is instantiated as an empty dictionary.

Type:

dict

wave

Wavelength vector applicable to all spatial positions.

Type:

numpy.ndarray

flux

Datacube flux array.

Type:

numpy.ndarray

ivar

Datacube flux inverse variance. Can be None.

Type:

numpy.ndarray

mask

Datacube mask. Can be None.

Type:

numpy.ndarray

sres

Datacube spectral resolution. Can be None.

Type:

numpy.ndarray

covar

Datacube spatial covariance. Can be None.

Type:

mangadap.util.covariance.Covariance

prihdr

Primary header for the datacube. If not provided on instantiation, set to an empty astropy.io.fits.Header.

Type:

astropy.io.fits.Header

fluxhdr

Header specifically for the flux array. If not provided on instantiation, set to be a copy of prihdr.

Type:

astropy.io.fits.Header

rss

The source row-stacked spectra used to build the datacube.

Type:

mangadap.spectra.rowstackedspectra.RowStackedSpectra

sigma_rho

The \(\sigma_{\rho}\) of the Gaussian function used to approximate the trend of the correlation coefficient with spaxel separation. Used to construct the approximate correlation matrix (see approximate_correlation_matrix()).

Type:

float

correl_rlim

The limiting radius of the image reconstruction (Gaussian) kernel in arcseconds. Used to construct the approximate correlation matrix (see approximate_correlation_matrix()).

Type:

float

approx_correl

Approximate correlation matrix; see approximate_correlation_matrix().

Type:

Covariance

_get_pixelscale()[source]

Measure the pixel scale using the WCS.

The method uses the coordinates of 4 adjacent pixels to compute the pixel area. Assuming the pixels are square and that the pixel scale is constant in square arcsec over the full image, the square-root of the area is the pixel scale.

The method will fault if wcs is None.

Returns:

The estimated pixel scale.

Return type:

float

_get_wavelength_vector()[source]

Use the astropy.wcs.WCS attribute (wcs) to generate the datacube wavelength vector.

wcs cannot be None and the wavelength coordinate system must be defined along its third axis.

Returns:

Vector with wavelengths along the third axis defined by wcs.

Return type:

numpy.ndarray

Raises:

ValueError – Raised if wcs or nwave is not defined.

approximate_correlation_matrix(sigma_rho, rlim, rho_tol=None, redo=False)[source]

Construct a correlation matrix with correlation coefficients that follow a Gaussian in 2D pixel separation.

This method constructs a correlation matrix with correlation coefficients defined as

\[\rho_{ij} = \exp(-D_{ij}^2 / 2\sigma_{\rho}^2)\]

where \(D_{ij}\) is the distance between two spaxels in the spatial dimension of the datacube (in the number spaxels, not arcsec). Any pixels with \(D_{ij} > R_{\rm lim}\) is set to zero, where \(R_{\rm lim}\) is the limiting radius of the kernel used in the datacube rectification; this is provided here as rlim in arcsec, which is converted to spaxels using pixelscale.

We found in Westfall et al. (2019, AJ, 158, 231) that this is a reasonable approximation for the formal covariance that results from Shepard’s rectification method.

There is an unknown relation between the dispersion of the kernel used by Shepard’s method and the value of \(\sigma_{\rho}\). Tests show that it is tantalizingly close to \(\sigma_{\rho} = \sqrt{2}\sigma\), where \(\sigma\) is in pixels instead of arcsec; however, a formal derivation of this hasn’t been done and is complicated by the focal-plane sampling of the row-stacked spectra.

Within the limits of how the focal-plane sampling changes with wavelength, we found the value of \(\sigma_{\rho}\) varies little with wavelength in MaNGA. Here, we assume \(\sigma_{\rho}\) is fully wavelength independent. Therefore, once this method is run once, it doesn’t need to be run again for different wavelength channels. To force the correlation matrix to be recreated, use redo or change the provided \(\sigma_{\rho}\).

Parameters:
  • sigma_rho (float) – The \(\sigma_{\rho}\) of the Gaussian function used to approximate the trend of the correlation coefficient with spaxel separation.

  • rlim (float) – The limiting radius of the image reconstruction (Gaussian) kernel in arcseconds.

  • rho_tol (float, optional) – Any correlation coefficient less than this is assumed to be equivalent to (and set to) 0.

  • redo (bool, optional) – Force the recalculation of the cube dimensions if they are already defined and \(\sigma_{\rho}\) has not changed.

Returns:

Correlation matrix

Return type:

mangadap.util.covariance.Covariance

approximate_covariance_cube(channels=None, sigma_rho=None, rlim=None, rho_tol=None, csr=False, quiet=False)[source]

Return the approximate covariance matrices for many wavelength channels.

This is a simple wrapper for approximate_covariance_matrix() that iteratively builds the covariance matrix for each wavelength channel.

If approx_correl is not yet constructed (see approximate_correlation_matrix()), sigma_rho and rlim must be provided.

Parameters:
  • channels (int, array-like, optional) – Indices of the spectral channels for which to calculate the covariance matrix. If None, the covariance matrix is calculated for all channels.

  • sigma_rho (float, optional) – The \(\sigma_{\rho}\) of the Gaussian function used to approximate the trend of the correlation coefficient with spaxel separation.

  • rlim (float, optional) – The limiting radius of the image reconstruction (Gaussian) kernel in arcseconds.

  • rho_tol (float, optional) – Any correlation coefficient less than this is assumed to be equivalent to (and set to) 0.

  • csr (bool, optional) – Instead of returning a mangadap.util.covariance.Covariance object, return the covariance matrix as a scipy.sparse.csr_matrix object. Primarily used by approximate_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, numpy.ndarray: If csr is True, the returned object is an numpy.ndarray of scipy.sparse.csr_matrix types.

approximate_covariance_matrix(channel, sigma_rho=None, rlim=None, rho_tol=None, csr=False, quiet=False)[source]

Return an approximate calculation of the covariance matrix assuming

\[C_{ij} = \rho_{ij}(V_{i} V_{j})^{1/2}\]

where \(\rho_{ij}\) is approximated by approximate_correlation_matrix() and \(V_i\equiv C_{ii}\) are the variances provided by the inverse of ivar.

The method first calculates \(\rho_{ij}\) if it hasn’t been yet or the provided sigma_rho and/or rlim values are different to previous calls, which builds approx_correl. If approx_correl is not yet constructed, sigma_rho and rlim must be provided.

The returned covariance matrix is the correlation matrix rescaled by the variance provided by ivar for the specified channel.

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the covariance matrix.

  • sigma_rho (float, optional) – The \(\sigma_{\rho}\) of the Gaussian function used to approximate the trend of the correlation coefficient with spaxel separation.

  • rlim (float, optional) – The limiting radius of the image reconstruction (Gaussian) kernel in arcseconds.

  • rho_tol (float, optional) – Any correlation coefficient less than this is assumed to be equivalent to (and set to) 0.

  • csr (bool, optional) – Instead of returning a mangadap.util.covariance.Covariance object, return the covariance matrix as a scipy.sparse.csr_matrix object. Primarily used by approximate_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, scipy.sparse.csr_matrix: The approximate covariance matrix for the designated wavelength channel. The return type depends on csr.

binned_on_sky_area(bin_indx)[source]

Compute the on-sky area of a set of binned spectra.

For each bin, this is just the number of spaxels in the bin multiplied by the spaxel area.

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\).

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

property can_analyze

Confirm that the DAP can analyze the datacube.

The only requirement is:

  • the \(cz\) velocity ('vel' in meta) must be greater than -500 km/s.

Returns:

Flag whether or not the DAP can analyze the datacube.

Return type:

bool

property can_compute_covariance

Determine if the object can be used to compute the spatial covariance.

If covar is currently not defined, the method tries to load the row-stacked spectra used to build the datacube; see load_rss(). If that is successful or if covar is already defined, the method returns True. If covar is None and the row-stacked spectra cannot be loaded, the method returns False.

Returns:

Flag that the object can be used to calculate the spatial covariance.

Return type:

bool

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

Determine the mean central wavelength for all spectra under various conditions.

The wavelength channel is set to be the center of the bandpass selected by waverange, weighted by the response function and the flux (if requested/provided). The mask is also incorporated in these calculations. By default (i.e., no wavelength limits or weighting) the wavelength channel is just the central wavelength of the full spectral range. (Note that if the spectra are binned logarithmically, this isn’t necessarily the central wavelength channel.)

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 as a weighting function for the calculation.

  • per_pixel (bool, optional) – When providing a response function, base the calculation on per pixel measurements, instead of per angstrom. Set to False for a per-angstrom calculation.

  • flag (str, list, optional) – One or more 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. If bitmask is None, these are ignored.

  • fluxwgt (bool, optional) – Flag to weight by the flux when determining the mean coordinates.

Returns:

The mean central wavelength of all spectra.

Return type:

float

copy_to_array(attr='flux', waverange=None, nbins=None, select_bins=None, missing_bins=None, unique_bins=None)[source]

Return a copy of the selected data array with a flattened spatial axis.

The array size is always \(N_{\rm spec} \times N_{\rm wavelength}\). The spatial positions within the original datacube for each spectrum are given by tuples in spatial_index.

See copy_to_masked_array() for argument descriptions.

Warning

Any masking is ignored in this function call. To incorporate the mask use copy_to_masked_array().

Returns:

A 2D array with a copy of the data from the selected attribute.

Return type:

numpy.ndarray

copy_to_masked_array(attr='flux', use_mask=True, flag=None, waverange=None, nbins=None, select_bins=None, missing_bins=None, unique_bins=None)[source]

Return a copy of the selected data array as a masked array with a flattened spatial axis.

The array size is always \(N_{\rm spec} \times N_{\rm wavelength}\). The spatial positions within the original datacube for each spectrum are given by tuples in spatial_index.

This is functionally identical to copy_to_array(), except the output format is a numpy.ma.MaskedArray. The pixels that are considered to be masked can be specified using the flag option.

Parameters:
  • attr (str, optional) – The attribute for the returned array. Can be ‘flux’, ‘ivar’, ‘sres’, or ‘mask’. For ‘mask’, you’re likely better off using copy_to_array(). Strings are always set to be lower-case, so capitalization shouldn’t matter.

  • use_mask (bool, optional) – Use the internal mask to mask the data. This is largely here to allow for copy_to_array() to wrap this function while not applying the internal mask.

  • 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.

  • flag (str, list, optional) – One or more bitmask flags used to select bad pixels. The names must be a valid bit name as defined by bitmask (see mangadap.util.bitmask.BitMask). If not provided, any non-zero mask bit is omitted.

  • nbins (int, optional) – The total number of defined bins. Default is to find the maximum number in the unique_bins list.

  • select_bins (array-like, optional) – A boolean array selecting spectra in the flattened cube (with shape \(N_{\rm spec} \times N_\lambda\)) to return.

  • missing_bins (list, optional) – A list of bin numbers that are missing from the output selection and that will be replaced with masked spectra. If specified, must also provide unique_bins.

  • unique_bins (array-like, optional) – The indices of the bins that have valid spectra. The length of this array should match the number of selected spectra from select_bins.

Returns:

A 2D array with a copy of the data from the selected extension, masked where mask is either True or with the selected bitmask flags.

Return type:

numpy.ma.MaskedArray

Raises:
  • AttributeError – Raised if attr is not a valid attribute.

  • ValueError – Raised if the unique bins are not specified, but the missing bins are, or if the number of unique bins does not match the number of spectra in the datacube.

covariance_cube(**kwargs)[source]

Construct the formal covariance matrix for one or more channels.

This is a simple wrapper for a call to mangadap.spectra.rowstackedspectra.RowStackedSpectra.covariance_cube() executed using rss, which cannot be None. See that method for the argument description.

Warning

The provided rectifications parameters should be the same as used to construct the datacube. If not, the calculated covariance matrix will not be correct.

covariance_matrix(channel, **kwargs)[source]

Construct the formal covariance matrix for the provided channel.

This is a simple wrapper for a call to mangadap.spectra.rowstackedspectra.RowStackedSpectra.covariance_matrix() executed using rss, which cannot be None. See that method for the argument description.

Warning

The provided rectifications parameters should be the same as used to construct the datacube. If not, the calculated covariance matrix will not be correct.

static do_not_fit_flags()[source]

Return the maskbit names that should not be fit.

The base class returns None.

static do_not_stack_flags()[source]

Return the maskbit names that should not be stacked.

The base class returns None.

static do_not_use_flags()[source]

Return the maskbit names that should not be used.

The base class returns None.

property file_path

Return the full path to the input file.

Warning

The required attributes for this method are not defined by the base class. If the attributes directory_path or file_name do not exist or if either of them are None, this returns None.

flux_stats(waverange=None, response_func=None, per_pixel=True, flag=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.

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 as a weighting function for the calculation.

  • per_pixel (bool, optional) – When providing a response function, continue to calculate the statistics per pixel. Set to False for a per-angstrom calculation.

  • flag (str, list, optional) – One or more 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.

Returns:

Three objects are returned: the mean flux, the propagated variance in the mean flux, and the mean S/N. The shape of each is the same as the shape of a single wavelength channel in the datacube (i.e., spatial_shape).

Return type:

numpy.ndarray

Raises:

ValueError – Raised of a provided wavelength range object does not have two elements or if a response function is provided and has an incorrect shape.

classmethod from_config(cfgfile, **kwargs)[source]

Construct a datacube object using a configuration file.

The method is undefined for the base class. If called, an exception is raised. Derived classes must override this method to allow for a configuration file to be used to instantiate the relevant DataCube subclass.

Parameters:
  • cfgfile (str) – Configuration file. See configparser.ConfigParser.

  • **kwargs – Any other keyword arguments that invoke optional instantiation methods. Note that these arguments will never be used in a command-line level execution of the DAP. They should only be available for custom scripts.

instrument = None

The name of the instrument used to collect the data. This is only used in construction of the root name for output files. See output_root().

interpolate_to_match(func, fill_value=0.0)[source]

Interpolate a function to match the datacube wavelength sampling.

Parameters:
  • func (numpy.ndarray) – Function to linear interpolate. Shape is \((N_{\rm wave}, 2)\), where the first column has the wavelength and the second has the function to interpolate. If None, simply returns a unity vector of the correct length.

  • fill_value (float, optional) – The value to use for spectral regions not sampled by the provided function.

Returns:

The function interpolated (not resampled) to match wave. Any spectral regions not within the wavelength range of the provided function is set to 0.

Return type:

numpy.ndarray

load_rss(**kwargs)[source]

Try to load the source row-stacked spectra for this datacube.

This method is undefined in the base class, and simply passes.

Derived classes should override this method if it’s possible to load the row-stacked spectra. If the load is successful, rss should no longer be None.

mean_sky_coordinates(center_coo=None)[source]

Compute the mean sky coordinates for each spectrum.

If the WCS is available, the coordinates can be returned in RA and declination (degrees). If center coordinates are provided (see center_coo), however, the coordinates are offset set as follows:

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

where \((\alpha_0, \delta_0)\) are the provided coordinates.

If the WCS is not available, the returned coordinates are in arcsec from the center of the image (regardless of the value of center_coo) determined using the pixelscale. At least in this case, the coordinates are assumed to be relative to the pixel center (not, e.g., its edge).

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.

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

metakeys()[source]

Get a list of the keys in the datacube metadata.

property nspec

Number of spectra in the datacube.

property output_root
populate_metadata()[source]

Populate and validate the DataCube metadata (in meta) to ensure the it can be analyzed by the DAP.

The only required metadata keyword is z, which sets the initial guess for the bulk redshift of the galaxy. If this key is not available or its value doesn’t meet the criterion below, this function will raise an exception, meaning the DAP will fault before it starts processing the DataCube.

The metadata provided must meet the following critical criteria or the method will fault:

  • Velocity (\(cz\)) must be greater than -500.

For the remainder of the metadata, if the keyword does not exist, if the value is None, or if the value is outside the accepted range, a default is chosen. The metadata keywords, acceptable ranges, and defaults are provided below.

Keyword

Range

Default

vdisp

\(\sigma > 0\)

100

ell

\(0 \leq \varepsilon < 1\)

None

pa

\(0 \leq \phi_0 < 360\)

None

reff

\(R_{\rm eff} > 0\)

None

ebv

\(E(B-V) \geq 0\)

None

drpcrit

True or False

False

All other metadata in meta are ignored.

If they do not already exist in meta, this method adds the following keywords:

  • z: The bulk redshift of the galaxy, used to calculate \(cz\).

  • vel: The initial guess velocity (\(cz\)) in km/s.

  • vdisp: The initial guess velocity dispersion in km/s.

  • ell: The isophotal ellipticity (\(1-b/a\)) to use when

    calculating semi-major axis coordinates.

  • pa: The isophotal position angle in deg from N through E, used when calculating semi-major axis coordinates.

  • reff: The effective radius in arcsec (DataCube WCS coordinates are expected to be in deg), used as a normalization of the semi-major axis radius in various output data.

  • ebv: The E(B-V) Galactic reddening along the line-of-sight to the galaxy.

  • drpcrit: A flag that indicates critical failures in the data reduction (True implies poor data quality). In MaNGA, the DRP will flag datacubes as having CRITICAL failures, but the survey approach was to run the DAP on these datacubes anyway. This flag is then used to propagate the caution in using the data analysis products to the user.

static propagate_flags()[source]

Flags that should be propagated from the observed data to the analyzed data.

The base class returns None.