mangadap.proc.stellarcontinuummodel module

A class hierarchy that performs the stellar-continuum fitting.

Revision history

14 Apr 2016: Implementation begun by K. Westfall (KBW)
19 Apr 2016: (KBW) First version
19 May 2016: (KBW) Added loggers and quiet keyword arguments to StellarContinuumModel, removed verbose
05 Jul 2016: (KBW) Removed oversample keyword from instantion of mangadap.proc.ppxffit.PPXFFit objects.
08 Nov 2016: (KBW) Moved StellarContinuumModel.reset_continuum_mask_window() from Elric to here. The function allows one to deal with the subtraction of the continuum over the fully viable spectral range, ignoring small spectral regions that were ignored during the stellar continuum fit. Also added wrapper function, StellarContinuumModel.emission_line_continuum_model().
11 Jan 2017: (KBW) Changed StellarContinuumModel.emission_line_continuum_model to StellarContinuumModel.unmasked_continuum_model().
23 Feb 2017: (KBW) Use DAPFitsUtil read and write functions.

License

Copyright © 2019, SDSS-IV/MaNGA Pipeline Group


class mangadap.proc.stellarcontinuummodel.StellarContinuumModel(method_key, binned_spectra, guess_vel, guess_sig=None, method_list=None, dapsrc=None, dapver=None, analysis_path=None, directory_path=None, output_file=None, hardcopy=True, symlink_dir=None, tpl_symlink_dir=None, clobber=False, checksum=False, loggers=None, quiet=False)[source]

Bases: object

Class that holds the stellar-continuum model results.

Parameters:
  • method_key (str) – The keyword that designates which method, provided in method_list, to use for the continuum-fitting procedure.
  • binned_spectra – (mangadap.proc.spatiallybinnedspectra.SpatiallyBinnedSpectra): The binned spectra to fit.
  • guess_vel (float, numpy.ndarray) – A single or spectrum-dependent initial estimate of the redshift in km/s (\(cz\)).
  • guess_sig (float, numpy.ndarray) – (Optional) A single or spectrum-dependent initial estimate of the velocity dispersion in km/s. Default is 100 km/s.
  • method_list (list) – (Optional) List of StellarContinuumModelDef objects that define one or more methods to use for the stellar continuum binning. The default list is provided by the config files in the DAP source directory and compiled into this list using available_stellar_continuum_modeling_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.
  • tpl_symlink_dir (str) – (Optional) Create a symbolic link to the created template library 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 for a prior?
_add_method_header(hdr)[source]

Add fitting method information to the header.

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

_construct_2d_hdu(good_snr, model_flux, model_mask, model_par)[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.

_fill_method_par(dapsrc=None, dapver=None, analysis_path=None, tpl_symlink_dir=None)[source]

Fill in any remaining modeling parameters.

This creates a full array for the guess redshifts and guess velocity dispersions.

It also creates/reads the template library.

_finalize_cube_mask(mask)[source]

Finalize the mask by setting the DIDNOTUSE, FORESTAR, and LOW_SNR masks

Returns:Bitmask array.
Return type:numpy.ndarray
_get_missing_models(debug=False)[source]
_initialize_primary_header(hdr=None)[source]

Initialize the header of hdu.

Returns:Edited header object.
Return type:astropy.io.fits.Header
_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 DAP reduction assessments file. 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().
all_except_emission_flags()[source]
all_spectrum_flags()[source]
construct_3d_hdu()[source]

Reformat the model spectra into a cube matching the shape of the DRP fits file.

construct_models(replacement_templates=None, redshift_only=False, deredshift=False, corrected_dispersion=False)[source]

Reconstruct the best fitting models.

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

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

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

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 default_paths(plate, ifudesign, rdxqa_method, binning_method, method_key, directory_path=None, drpver=None, dapver=None, analysis_path=None, output_file=None)[source]

Set the default directory and file name for the output file.

Parameters:
Returns:

Two strings with the path for the output file and the name of the output file.

Return type:

str

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

Select the method

file_name()[source]

Return the name of the output file.

file_path()[source]

Return the full path to the output file.

fill_to_match(binid, replacement_templates=None, redshift_only=False, deredshift=False, corrected_dispersion=False, missing=None)[source]

Get the stellar-continuum model from this objects that matches the input bin ID matrix. The output is a 2D matrix ordered by the bin ID; any skipped index numbers in the maximum of the union of the unique numbers in the binid and missing input are masked.

All this does is find the stellar continuum from spaxels in this model and match them to the input spaxels. Depending on how differently the data data sets have been binned, this says nothing about whether or not the returned models are a good fit to the data!

replacement_templates only works if the number of templates is the same as those used during the actual fitting process.

use redshift_only if you want to exclude the best-fitting dispersion from the model.

fit(binned_spectra, guess_vel, guess_sig=None, dapsrc=None, dapver=None, analysis_path=None, directory_path=None, output_file=None, hardcopy=True, symlink_dir=None, tpl_symlink_dir=None, clobber=False, loggers=None, quiet=False)[source]

Fit the binned spectra using the specified method.

Parameters:
  • binned_spectra – (mangadap.proc.spatiallybinnedspectra.SpatiallyBinnedSpectra): The binned spectra to fit.
  • guess_vel (float, numpy.ndarray) – A single or spectrum-dependent initial estimate of the redshift in km/s (\(cz\)).
  • guess_sig (float, numpy.ndarray, optional) – A single or spectrum-dependent initial estimate of the velocity dispersion in km/s.
  • 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 paths, used to override the default defined by mangadap.config.defaults.default_dap_version(). This only sets the path names and does not change the version of the code being used.
  • 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, optional) – 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.
  • tpl_symlink_dir (str, optional) – Create a symbolic link to the created template library 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.
get_template_library(dapsrc=None, dapver=None, analysis_path=None, tpl_symlink_dir=None, velocity_offset=None, match_to_drp_resolution=False, resolution_fwhm=None, hardcopy=False)[source]
info()[source]
matched_kinematics(binid, redshift=None, dispersion=100.0, constant=False, cz=False, corrected=False, min_dispersion=None, nearest=False, missing=None)[source]

Return the guess redshift and velocity dispersion for all the spectra based on the stellar kinematics.

For spectra the were not fit, use the median (unmasked) kinematic measurement if no default value is provided as an argument.

Parameters:
  • binid (numpy.ndarray) – 2D array with the bin ID associate with each spaxel in the datacube. Shape must be the same as spatial_shape.
  • redshift (float) – (Optional) The default redshift to use for spectra without a stellar-continuum fit. Default is to use the median of the unmasked measurements
  • dispersion (float) – (Optional) The default velocity dispersion to use for spectra without a stellar-continuum fit. Default is 100 km/s.
  • constant (bool) – (Optional) Force the function to return a constant redshift and dispersion for each spectrum, regardless of any fitted kinematics.
  • cz (bool) – (Optional) Return the redshift as cz in km/s, as opposed to unitless z.
  • corrected (bool) – (Optional) By default, the returned velocity dispersions are the measured values, which will include any resolution difference between the templates and the galaxy data. Setting corrected to True will return the corrected dispersions.
  • min_dispersion (float) – (Optional) Impose a minimum dispersion.
  • nearest (bool) – (Optional) Instead of the median of the results for the spectra that were not fit, use the value from the nearest bin.
  • missing (list) – (Optional) Any bin ID numbers missing from the input binid image needed for constructing the output matched data.
Returns:

Returns (unmasked!) arrays with a redshift (or cz) and dispersion for each binned spectrum.

Return type:

numpy.ndarray

Raises:

TypeError – Raised if the input redshift or dispersion values are not single numbers.

matched_template_flags(binned_spectra)[source]

TODO: Function out of date!

Return the templates used during the fit to each spectrum, matched to the spectra in the binned_spectra object. For spectra with no models, the flags are all set to true.

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

Read an existing file with an existing set of stellar-continuum models.

static reset_continuum_mask_window(continuum, dispaxis=1, quiet=False)[source]

Reset the mask of the stellar continuum to a continuous window from the minimum to maximum valid wavelength.

Todo

  • Allow continuum to be n-dimensional and use numpy.apply_along_axis.
unmasked_continuum_model(replacement_templates=None, redshift_only=False, deredshift=False, corrected_dispersion=False)[source]

Return the stellar continuum unmasked over the full continuous fitting region. Models are reconstructed based on the model parameters if new template fluxes are provided or if the returned models should not include the LOSVD convolution (redshift_only=True).

write(match_DRP=False, clobber=False)[source]

Write the hdu object to the file.

class mangadap.proc.stellarcontinuummodel.StellarContinuumModelBitMask[source]

Bases: mangadap.util.dapbitmask.DAPBitMask

Derived class that specifies the mask bits for the stellar-continuum modeling. The maskbits defined are:

Key Bit Description
DIDNOTUSE 0 Pixel was ignored because it was flagged in the binning step; see the DRPFits.do_not_stack_flags() and SpatiallyBinnedSpectra.do_not_fit_flags().
FORESTAR 1 Pixel was ignored because it was flagged as FORESTAR by the DRP.
LOW_SNR 2 Pixel was ignored because the S/N estimate of the spectrum was below the set threshold; see header keyword SCMINSN.
ARTIFACT 3 Pixel was ignored during the stellar-continuum fit because it was designated as containing an artifact.
OUTSIDE_RANGE 4 Pixel was ignored during the stellar-continuum fit because it was outside the designated spectral range of the fit; see the header keyword FITWAVE.
EML_REGION 5 Pixel was ignored during the stellar-continuum fit because it contains an emission-line.
TPL_PIXELS 6 These pixels were removed to ensure that the number of template spectral pixels was >= the number of fitted object pixels during the stellar-continuum fit.
TRUNCATED 7 This region was truncated to avoid convolution errors at the edges of the spectral range during the stellar-continuum fit.
PPXF_REJECT 8 Pixel was rejected during the pPXF fit via the clean parameter.
INVALID_ERROR 9 Pixel ignored because the flux error was invalid.
NO_FIT 10 Spectrum not fit such that not parameters are provided.
MAXITER 11 The fit optimizer reached the maximum number of iterations during the fit, which may imply failure to converge.
LARGE_CHI2 12 This pixel has a large chi^2 (definition of large TBD).
LARGE_RESID 13 This pixel has a large residual (definition of large TBD).
INSUFFICIENT_DATA 14 There were insufficient data in the fitting window to fit the line profile(s).
FIT_FAILED 15 Stellar-continuum fit failed according to status returned by scipy.optimize.least_squares.
NEAR_BOUND 16 Fit parameters are within at or near the imposed boundaries in parameter space.
NEGATIVE_WEIGHTS 17 Weights for templates were negative.
BAD_SIGMA 18 Corrected velocity dispersion is below 50 km/s or above 400 km/s
MIN_SIGMA 19 The fitted velocity dispersion is at the minimum allowed by pPXF (1/100th of a pixel)
BAD_SIGMACORR_SRES 20 The instrumental dispersion correction for the velocity dispersion, based on the quadrature difference in the spectral resolution, is invalid.
BAD_SIGMACORR_EMP 21 The instrumental dispersion correction for the velocity dispersion, based on fitting the native resolution template to the resolution matched template, is invalid.
cfg_root = 'stellar_continuum_model_bits'
class mangadap.proc.stellarcontinuummodel.StellarContinuumModelDef(key=None, minimum_snr=None, fit_type=None, fitpar=None, fitclass=None, fitfunc=None)[source]

Bases: mangadap.par.parset.KeywordParSet

A class that holds the parameters necessary to perform the stellar-continuum fitting.

The provided fitfunc must have the form:

model_wave, model_flux, model_mask, model_par                 = fitfunc(self, binned_spectra, par=None)

where binned_spectra has type mangadap.proc.spatiallybinnedspetra.SpatiallyBinnedSpectra, and the returned objects are the wavelength vector, the fitted model flux, a bitmask for the model, and the set of model parameters. For example, see mangadap.proc.ppxffit.PPXFFit.fit_SpatiallyBinnedSpectra().

The defined parameters are:

Key Type Options Default Description
key str Keyword used to distinguish between different spatial binning schemes.
minimum_snr int, float Minimum S/N of spectrum to fit
fit_type str Currently can be anything. In the DAP, this is used to identify the primary purpose of the fit as either producing stellar kinematics, composition measurements, or emission line fits; see mangadap.proc.spectralfitting. The purpose for this type is to isolate the expected format of the binary table data; see, e.g., mangadap.proc.spectralfitting._per_stellar_kinematics_dtype().
fitpar ParSet, dict Any additional parameters, aside from the spectra themselves, required by the fitting function.
fitclass Undefined Instance of class object to use for the model fitting. Needed in case fitfunc is a non-static member function of a class.
fitfunc Undefined The function that models the spectra
mangadap.proc.stellarcontinuummodel.available_stellar_continuum_modeling_methods(dapsrc=None)[source]

Return the list of available stellar-continuum modeling methods.

pPXF methods:

Todo

Fill in

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.stellarcontinuummodel.StellarContinuumModelDef() objects, each defining a stellar-continuum modeling approach.

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.stellarcontinuummodel.validate_stellar_continuum_modeling_method_config(cnfg)[source]

Validate the mangadap.util.parser.DefaultConfig object with the stellar-continuum-modeling method parameters.

Parameters:

cnfg (mangadap.util.parser.DefaultConfig) – Object with the stellar-continuum-modeling method parameters to validate.

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.