mangadap.proc.spatiallybinnedspectra module

A class hierarchy that perform that spatially bins a set of two-dimensional data.

The base class allows for user-defined definitions of binning procedures.

Todo

  • Check binning an RSS file

Revision history

01 Apr 2016: Implementation begun by K. Westfall (KBW)
19 May 2016: (KBW) Include SPECRES and SPECRESD extensions from DRP file in output. Added loggers and quiet keyword arguments to SpatiallyBinnedSpectra, removed verbose
06 Jul 2016: (KBW) Make the application of a reddening correction an input parameter.
25 Aug 2016: (KBW) Fixed error in bin area when calling SpatiallyBinnedSpetra._unbinned_data_table()
01 Dec 2016: (KBW) Include keyword that describes how to handle the spectral resolution.
02 Dec 2016: (KBW) Incorporate mangadap.util.extinction.GalacticExtinction. Revert main file structure to be bin-based instead of cube-based; include convenience functions that construct the cube for each extension as requested.
06 Dec 2016: (KBW) Significantly restructured.
23 Feb 2017: (KBW) Use DAPFitsUtil read and write functions.
21 Aug 2017: (KBW) Use the new PRE-pixelized assessments of the LSF.

License

Copyright © 2019, SDSS-IV/MaNGA Pipeline Group


class mangadap.proc.spatiallybinnedspectra.SpatiallyBinnedSpectra(method_key, drpf, rdxqa, reff=None, method_list=None, dapsrc=None, dapver=None, analysis_path=None, directory_path=None, output_file=None, hardcopy=True, symlink_dir=None, clobber=False, checksum=False, loggers=None, quiet=False)[source]

Bases: object

Class that holds spatially binned spectra.

Parameters:
  • method_key (str) – The keyword that designates which method, provided in method_list, to use for the binning procedure.
  • drpf (mangadap.drpfits.DRPFits) – The DRP datacube with the spectra to bin.
  • rdxqa (mangadap.proc.reductionassessments.ReductionAssessments) – The basic assessments of the DRP data that are needed for the binning procedures.
  • reff (float) – (Optional) The effective radius of the galaxy.
  • method_list (list) – (Optional) List of SpatiallyBinnedSpectraDef objects that define one or more methods to use for the spatial binning. Default is to use the config files in the DAP source directory to construct the available methods using available_spatial_binning_methods().
  • dapsrc (str) – (Optional) Root path to the DAP source directory. If not provided, the default is defined by mangadap.config.defaults.dap_source_dir().
  • dapver (str) – (Optional) The DAP version to use for the analysis, used to override the default defined by mangadap.config.defaults.default_dap_version().
  • analysis_path (str) – (Optional) The top-level path for the DAP output files, used to override the default defined by mangadap.config.defaults.default_analysis_path().
  • directory_path (str) – The exact path to the directory with DAP output that is common to number DAP “methods”. See directory_path.
  • output_file (str) – (Optional) Exact name for the output file. The default is to use mangadap.config.defaults.default_dap_file_name().
  • hardcopy (bool) – (Optional) Flag to write the HDUList attribute to disk. Default is True; if False, the HDUList is only kept in memory and would have to be reconstructed.
  • symlink_dir (str) – (Optional) Create a symbolic link to the created file in the supplied directory. Default is to produce no symbolic link.
  • clobber (bool) – (Optional) Overwrite any existing files. Default is to use any existing file instead of redoing the analysis and overwriting the existing output.
  • checksum (bool) – (Optional) Use the checksum in the fits header to confirm that the data has not been corrupted. The checksum is always written to the fits header when the file is created; this argument does not toggle that functionality.
  • loggers (list) – (Optional) List of logging.Logger objects to log progress; ignored if quiet=True. Logging is done using mangadap.util.log.log_output(). Default is no logging.
  • quiet (bool) – (Optional) Suppress all terminal and logging output. Default is False.
loggers

List of logging.Logger objects to log progress; ignored if quiet=True. Logging is done using mangadap.util.log.log_output().

Type:list
quiet

Suppress all terminal and logging output.

Type:bool

Todo

  • Allow velocity offsets for registration.
_add_method_header(hdr)[source]

Add method-specific metadata to the header.

_add_reddening_header(hdr)[source]

Set the relevant reddening information to the header.

Parameters:hdr (astropy.io.fits.Header) – Input header object to edit.
Returns:Edited header object.
Return type:astropy.io.fits.Header
_apply_reddening(flux, ivar, sdev, covar, deredden=True)[source]
_assign_image_arrays()[source]

Set image_arrays, which contains the list of extensions in hdu that are on-sky image data.

_assign_spectral_arrays()[source]

Set spectral_arrays, which contains the list of extensions in hdu that contain spectral data.

_binned_data_table(bin_indx, stack_flux, stack_ivar, per_pixel=True)[source]

Construct the output data table for the binned spectra.

Parameters:
  • bin_indx (numpy.ndarray) – The integer vector with the bin associated with each spectrum in the DRP cube. This is the flattened BINID array.
  • stack_flux (numpy.ndarray) – The stacked spectra, organized as \(N_{\rm spec}\times\N_\lambda}\).
  • stack_flux – The stacked inverse variance, organized as \(N_{\rm spec}\times\N_\lambda}\).
Returns:

The record array that is put in the BINS extension of hdu.

Return type:

numpy.recarray

_check_snr()[source]

Determine which spectra in rdxqa have a S/N greater than the minimum set by method. Only these spectra will be included in the binning.

Returns:Boolean array for the spectra that satisfy the criterion.
Return type:numpy.ndarray
_construct_2d_hdu(bin_indx, good_fgoodpix, good_snr, bin_data, stack_flux, stack_sdev, stack_ivar, stack_npix, stack_mask, stack_sres, stack_covar)[source]

Construct hdu that is held in memory for manipulation of the object. See construct_3d_hdu() if you want to convert the object into a DRP-like datacube.

bin_indx is 2d

stack_sres is expected to be a MaskedArray

_fill_method_par(good_spec)[source]

Finalize the binning parameters, as needed.

For the radial binning, set the ellipticity and position angle. Set scale radius to the effective radius if desired.

For the Voronoi binning, set the signal and noise.

Parameters:good_spec (numpy.ndarray) – List of spectra to include in the binning. See check_fgoodpix() and _check_snr().
_finalize_cube_mask(mask)[source]

Finalize the mask after the 2D mask in self has already been reconstructed into a 3D cube.

This mostly handles the masks for regions outside the IFU field of view.

This propagates the DIDNOTUSE and FORESTAR bits from the DRP file, and sets the LOW_SPECCOV and LOW_SNR bits based on the input boolean arrays.

Returns:Bitmask array.
Return type:numpy.ndarray
static _get_missing_bins(unique_bins)[source]
_initialize_primary_header(hdr=None)[source]

Initialize the header of hdu.

Returns:Edited header object.
Return type:astropy.io.fits.Header
_interpolated_response_function()[source]
_per_bin_dtype()[source]

Construct the record array data type that holds the information for the binned spectra.

Todo

There currently is no mask; however, could add masks for:
  • Insufficient good wavelength channels in the spectrum.
  • Variance in flux in bin is too large.
_set_paths(directory_path, dapver, analysis_path, output_file)[source]

Set the directory_path and output_file. If not provided, the defaults are set using, respectively, mangadap.config.defaults.default_dap_common_path() and mangadap.config.defaults.default_dap_file_name().

Parameters:
  • directory_path (str) – The exact path to the directory with DAP output that is common to number DAP “methods”. See directory_path.
  • dapver (str) – DAP version.
  • analysis_path (str) – The path to the top-level directory containing the DAP output files for a given DRP and DAP version.
  • output_file (str) – The name of the file with the reduction assessments. See compute().
_unbinned_data_table(bin_indx)[source]

Construct the output data table for the unbinned spectra.

Parameters:bin_indx (numpy.ndarray) – The integer vector with the bin associated with each spectrum in the DRP cube. This is the flattened BINID array.
Returns:The record array that is put in the BINS extension of hdu.
Return type:numpy.recarray
above_snr_limit(sn_limit, debug=False)[source]

Flag bins above a provided S/N limit.

bin_spectra(drpf, rdxqa, reff=None, dapver=None, analysis_path=None, directory_path=None, output_file=None, hardcopy=True, symlink_dir=None, clobber=False, loggers=None, quiet=False)[source]

Bin and stack the spectra.

Parameters:
  • drpf (mangadap.drpfits.DRPFits) – The DRP datacube with the spectra to bin.
  • rdxqa (mangadap.proc.reductionassessments.ReductionAssessments) – The basic assessments of the DRP data that are needed for the binning procedures.
  • reff (float) – (Optional) The effective radius of the galaxy.
  • dapver (str) – (Optional) The DAP version to use for the analysis, used to override the default defined by mangadap.config.defaults.default_dap_version().
  • analysis_path (str) – (Optional) The top-level path for the DAP output files, used to override the default defined by mangadap.config.defaults.default_analysis_path().
  • directory_path (str) – The exact path to the directory with DAP output that is common to number DAP “methods”. See directory_path.
  • output_file (str) – (Optional) Exact name for the output file. The default is to use mangadap.config.defaults.default_dap_file_name().
  • hardcopy (bool) – (Optional) Flag to write the HDUList attribute to disk. Default is True; if False, the HDUList is only kept in memory and would have to be reconstructed.
  • symlink_dir (str) – (Optional) Create a symbolic link to the created file in the supplied directory. Default is to produce no symbolic link.
  • clobber (bool) – (Optional) Overwrite any existing files. Default is to use any existing file instead of redoing the analysis and overwriting the existing output.
  • loggers (list) – (Optional) List of logging.Logger objects to log progress; ignored if quiet=True. Logging is done using mangadap.util.log.log_output(). Default is no logging.
  • quiet (bool) – (Optional) Suppress all terminal and logging output. Default is False.
check_fgoodpix(minimum_fraction=0.8)[source]

Determine which spaxels in rdxqa have a fractional spectral coverage of greater than the provided minimum fraction. Only these spectra will be included in the binning.

Parameters:minimum_fraction (float, optional) – The minimum fraction of the spectrum that must be valid for the spectrum to be included in any bin. Default is 0.8.
Returns:Boolean array for the spectra that satisfy the criterion. The shape of the array is \((N_{\rm spaxel},)\).
Return type:numpy.ndarray
construct_3d_hdu()[source]

Reformat the binned spectra into a cube matching the shape of the DRP fits file from which it was derived.

copy_to_array(ext='FLUX', waverange=None, include_missing=False)[source]

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

Return a copy of the selected data array. The array size is always \(N_{\rm models} \times N_{\rm wavelength}\); i.e., the data is always flattened to two dimensions and the unique spectra are selected.

Parameters:
  • ext (str) – (Optional) Name of the extension from which to draw the data. Must be allowed for the current mode; see data_arrays. Default is 'FLUX'.
  • 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.
  • include_missing (bool) – (Optional) Create an array with a size that accommodates the missing models.
Returns:

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

Return type:

numpy.ndarray

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

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

Return a copy of the selected data array as a masked array. 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 flag.

Parameters:
  • ext (str) – (Optional) Name of the extension from which to draw the data. Must be allowed for the current mode; see data_arrays. Default is ‘FLUX’.
  • 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. If not provided, ANY non-zero mask bit is omitted.
  • 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.
  • include_missing (bool) – (Optional) Create an array with a size that accommodates the missing models.
Returns:

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

Return type:

numpy.ndarray

static define_method(method_key, method_list=None, dapsrc=None)[source]

Select the method

static do_not_fit_flags()[source]
file_name()[source]

Return the name of the output file.

file_path()[source]

Return the full path to the output file.

find_nearest_bin(input_bins, weighted=False, indices=False)[source]

Use a KDTree to find the bins nearest to and excluding a list of input bins.

Parameters:
  • input_bins (array-like) – One or more bin ID numbers to use to locate the bin nearest to it based on its on-sky coordinates. The list must be a unique set.
  • weighted (bool) – (Optional) Use the weighted coordinates (LW_SKY_COO) instead of the unweighted coordinates (SKY_COO) when finding the nearest bin.
  • indices (bool) – (Optional) Return the indices of the nearest bins instead of their ID number (default).
Returns:

The bin IDs, one per input bin, of the bin closest to each input bin. Any bins that are in missing_bins have a return value of -1; there are no coordinates for these bins.

Return type:

numpy.ndarray

Raises:

ValueError – Raised if the set of input bins is not unique or contains bin IDs that are not present.

get_bin_indices(bins)[source]

Return the indices of the bins in the BIN table.

Parameters:bins (array-like) – The bin ID numbers to find
Returns:Integer array with the index of each bin ID in the BINID columns of the BINS extension.
Return type:numpy.ndarray
get_mask_array_from_drpfits(drpf, select=None)[source]

Convert the DRP mask for individual spaxels to the binned spectra mask.

Any pixel with bit flags in the list returned by mangadap.drpfits.DRPFits.do_not_stack_flags() are consolidated into the DIDNOTUSE flag. Any FORESTAR flags are propagated.

info()[source]
is_unbinned

Determine if the spectra are unbinned.

read(ifile=None, strict=True, checksum=False)[source]

Read an existing file with a previously binned set of spectra.

Parameters:
  • ifile (str) – (Optional) Name of the file with the data. Default is to use the name provided by file_path().
  • strict (bool) – (Optional) Force a strict reading of the file to make sure that it adheres to the expected format. At the moment, this only checks to make sure the method keyword printed in the header matches the expected value in method.
  • checksum (bool) – (Optional) Use the checksum in the header to check for corruption of the data. Default is False.
Raises:
  • FileNotFoundError – Raised if the file does not exist.
  • ValueError – Raised if strict is true and the header keyword ‘BINKEY’ does not match the method keyword.
replace_with_data_from_nearest_bin(data, bad_bins)[source]

Replace data in the list of provided bad bins with the data from the nearest good bin.

Parameters:
  • data (array-like) – Data for each bin. The length must be the same as nbins.
  • bad_bins (array-like) – The list of indices (must not be a boolean array) with bad values to be replaced.
Returns:

A new array with the bad data filled with the data from the nearest bin.

Return type:

numpy.ndarray

Raises:

ValueError – Raised if the input array doesn’t have the correct shape or if the list of bad bins has numbers outside the viable range (0,self.nbins-1).

static spectral_resolution_options()[source]

Return the allowed options for treating the spectral resolution.

Options are:

‘spaxel’: If available, use the spectral resolution determined for each spaxel. This is pulled from the ‘DISP’ extension in the DRP file; an exception will be raised if this extension does not exist!

‘cube’: Only consider the median spectral resolution determine for the entire datacube. This is pulled from the ‘SPECRES’ extension in the DRP file; an exception will be raised if this extension does not exist!

Returns:List of the available method keywords.
Return type:list
write(match_DRP=False, clobber=False)[source]

Write the hdu object to the file.

class mangadap.proc.spatiallybinnedspectra.SpatiallyBinnedSpectraBitMask[source]

Bases: mangadap.util.dapbitmask.DAPBitMask

Derived class that specifies the mask bits for the spatially binned spectra. The maskbits defined are:

Key Bit Description
DIDNOTUSE 0 Pixel was ignored because it was flagged as either DONOTUSE or FORESTAR by the data reduction pipeline (DRP).
FORESTAR 1 Pixel was ignored because it was flagged as FORESTAR by the data reduction pipeline (DRP).
LOW_SPECCOV 2 Pixel was ignored because the fraction of valid spectral channels limited the spectral coverage below the set threshold; see header keyword FSPECCOV.
LOW_SNR 3 Pixel was ignored because the S/N estimate of the spectrum was below the set threshold; see header keyword BINMINSN.
NONE_IN_STACK 4 No valid pixels were available or valid for the stacked spectrum.
NO_STDDEV 5 Insufficient pixels to calculate the standard deviation in the stacked spectrum.
IVARINVALID 6 Pixel ignored because inverse variance invalid.
cfg_root = 'spatially_binned_spectra_bits'
class mangadap.proc.spatiallybinnedspectra.SpatiallyBinnedSpectraDef(key=None, galactic_reddening=None, galactic_rv=None, minimum_snr=None, binpar=None, binclass=None, binfunc=None, stackpar=None, stackclass=None, stackfunc=None, spec_res=None, prepixel_sres=None)[source]

Bases: mangadap.par.parset.KeywordParSet

A class that holds the two parameter sets and the key designator for the binning scheme.

The provided binfunc method must have a call with the following form:

id = bin(x, y, par=par)

where id is the index of the bin to place each spectrum and x and y are the on-sky Cartesian coordinates; e.g.:

x = ReductionAssessment.hdu['SPECTRUM'].data['SKY_COO'][:,0]
y = ReductionAssessment.hdu['SPECTRUM'].data['SKY_COO'][:,1]

The provided stackfunc method must have a call with the following form:

stack_wave, stack_flux, stack_sdev, stack_npix, stack_ivar,                 stack_sres, stack_covar = stack(drpf, id, par=par)

where drpf is a mangadap.drpfits.DRPFits object. Note that the wavelengths are not returned because the input and output wavelength ranges are expected to be the same!

As long as they are mutable, the values in par can change, meaning that some products of the bin algorithm can be passed to the stack algorithm. For example, if you want to weight the inclusion of the spectrum in the bin, you’d have to provide both the binning and stacking routines. Actually, if that’s the case, you’re better off providing a class object that will both bin and stack the spectra!

The defined parameters are:

Key Type Options Default Description
key str Keyword used to distinguish between different spatial binning schemes.
galactic_reddening str The string identifier for the Galactic extinction curve to use. See mangadap.util.extinction.GalacticExtinction.valid_forms() for the available curves. Default is ODonnell.
galactic_rv int, float Ratio of V-band extinction to the B-V reddening. Default is 3.1.
minimum_snr int, float Minimum S/N of spectra to include in any bin.
binpar ParSet, dict The parameter set defining how to place each spectrum in a bin.
binclass Undefined Instance of class object to use for the binning. Needed in case binfunc is a non-static member function of the class.
binfunc Undefined The function that determines which spectra go into each bin.
stackpar ParSet, dict The parameter set defining how to stack the spectra in each bin.
stackclass Undefined Instance of class object to used to stack the spectra. Needed in case stackfunc is a non-static member function of the class.
stackfunc Undefined The function that stacks the spectra in a given bin.
spec_res str spaxel, cube Keyword defining the treatment of the spectral resolution. See SpatiallyBinnedSpectra.spectral_resolution_options() for a list of the options.
prepixel_sres bool Use the prepixelized version of the LSF measurements.
mangadap.proc.spatiallybinnedspectra.available_spatial_binning_methods(dapsrc=None)[source]

Return the list of available binning schemes.

The list of available binning schemes is defined by the set of configuration files at:

config_path = os.path.join(dapsrc, 'python', 'mangadap', 'config', 'spatial_binning')
ini_files = glob.glob(os.path.join(config_path, '/*.ini'))
Parameters:

dapsrc (str) – (Optional) Root path to the DAP source directory. If not provided, the default is defined by mangadap.config.defaults.dap_source_dir().

Returns:

A list of mangadap.proc.spatiallybinnedspectra.SpatiallyBinnedSpectraDef() objects, each defining a separate binning method.

Return type:

list

Raises:
  • NotADirectoryError – Raised if the provided or default dapsrc is not a directory.
  • OSError/IOError – Raised if no binning scheme configuration files could be found.
  • KeyError – Raised if the binning method keywords are not all unique.

Todo

  • Somehow add a python call that reads the databases and constructs the table for presentation in sphinx so that the text above doesn’t have to be edited with changes in the available databases.
mangadap.proc.spatiallybinnedspectra.validate_spatial_binning_scheme_config(cnfg)[source]

Validate the mangadap.util.parser.DefaultConfig object with spatial-binning scheme parameters.

Parameters:

cnfg (mangadap.util.parser.DefaultConfig) – Object with the spatial-binning method parameters as needed by mangadap.proc.spatiallybinnedspectra.SpatiallyBinnedSpectraDef.

Raises:
  • KeyError – Raised if any required keywords do not exist.
  • ValueError – Raised if keys have unacceptable values.
  • FileNotFoundError – Raised if a file is specified but could not be found.