mangadap.proc.emissionlinetemplates module

Class that constructs a set of emission-line templates, primarily for use with Sasuke.

Class usage examples

To construct emission-line templates, you need a wavelength vector, the instrumental resolution at which to construct the templates, and an emission line database (see mangadap.par.emissionlinedb.EmissionLineDB). A simple construction would be:

# Imports
import numpy
from mangadap.par.emissionlinedb import EmissionLineDB
from mangadap.proc.emissionlinetemplates import EmissionLineTemplates

wave = numpy.logspace(*(numpy.log10([3600,10300]), 4563)
sigma_inst = 30                     # Instrumental resolution in km/s

tpl = EmissionLineTemplates(wave, sigma_inst, emldb=EmissionLineDB.from_key('ELPMILES'))

License

Copyright © 2019, SDSS-IV/MaNGA Pipeline Group


class mangadap.proc.emissionlinetemplates.EmissionLineTemplates(wave, sigma_inst, log=True, base=10, emldb=None, flux_density=True, loggers=None, quiet=False)[source]

Bases: object

Construct a set of emission-line templates based on an emission-line database, primarily for use in the Sasuke emission-line fitter.

The templates are constructed based on the constraints provided by the emission-line database. See EmissionLinePar for the structure of each row in the database and an explanation for each of its columns.

Only lines with action=f are included in any template. The array tpli provides the index of the template that contains each line in the emission-line database. Lines that are not assigned to any template — either because they do not have action=f or their center lies outside the wavelength range in wave — are given an index of -1.

Only lines with mode=a (i.e., flux, velocity, and velocity dispersion are all tied) are included in the same template.

Lines with tied velocities are assigned the same velocity component (vgrp) and lines with the tied velocity dispersions are assigned the same sigma component (sgrp).

Warning

The construction of templates for use with Sasuke does not allow one to tie fluxes while leaving the velocities and/or velocity dispersions as independent.

Parameters:
  • wave (array-like) – A single wavelength vector with the wavelengths for the template spectra. The wavelengths are expected to be sample either linearly or geometrically (see log).

  • sigma_inst (float, array-like) – The single value or value as a function of wavelength for the instrumental dispersion (km/s) to use for the template construction.

  • log (bool, optional) – Flag that the wavelengths have been sampled geometrically.

  • base (float, optional) – Base for the geometric sampling.

  • emldb (mangadap.par.emissionlinedb.EmissionLineDB, optional) – Emission-line database that is parsed to setup which lines to include in the same template because they are modeled as having the same velocity, velocity dispersion and flux ratio. If not provided, no templates are constructed in instantiation; to build the templates using an existing instantiation, use build_templates().

  • flux_density (bool, optional) – Return spectrum in units of flux density (flux per angstrom). Default is to return the spectrum in units of flux per pixel.

  • 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. This can be reset in some methods.

  • quiet (bool, optional) – Suppress all terminal and logging output.

wave

Array with the wavelength (angstroms) of each pixel for all the constructed templates. Shape is \((N_{\rm pix},)\).

Type:

numpy.ndarray

sigma_inst

The object used to interpolate the instrumental dispersion (km/s) at the rest wavelength of each spectral line.

Type:

scipy.interpolate.interp1d

log

Flag that the spectrum is sampled geometrically.

Type:

bool

base

Base for the geometric sampling.

Type:

float

emldb

Database with the emission-line parameters. Shape is \((N_{\rm pix},)\).

Type:

mangadap.par.emissionlinedb.EmissionLineDB

ntpl

Total number of templates.

Type:

int

flux

The template spectra. Shape is \((N_{\rm tpl},N_{\rm pix})\).

Type:

numpy.ndarray

tpli

The index of the template containing each emission line. Any emission-lines with tpli==-1 means that the emission line was not included in any template, which should only occur for lines with action==i in the database. Shape is \((N_{\rm eml},)\).

Type:

numpy.ndarray

comp

The kinematic component assigned to each template. Templates with the same component number are forced to have the same velocity and velocity dispersion by pPXF. Shape is \((N_{\rm tpl},)\).

Type:

numpy.ndarray

vgrp

The velocity group assigned to each template. Templates in the same velocity group have their velocity parameters tied in pPXF, but the velocity dispersion parameters are independent. Shape is \((N_{\rm tpl},)\).

Type:

numpy.ndarray

sgrp

The velocity disperison (sigma) group assigned to each template. Templates in the same sigma group have their velocity dispersion parameters tied in pPXF, but the velocity parameters are independent. Shape is \((N_{\rm tpl},)\).

Type:

numpy.ndarray

A_ineq

A matrix used to constrain the kinematic parameters between kinematic components. See the constr_kinem for pPXF (version 7.0 and later). Shape is \((N_{\rm constr}, N_{\rm comp})\); i.e., the number of constraints applied by the number of kinematic components.

Type:

numpy.ndarray

b_ineq

A vector used to constrain the kinematic parameters between kinematic components. See the constr_kinem for pPXF (version 7.0 and later). Shape is \((N_{\rm constr},)\); i.e., the number of constraints applied.

Type:

numpy.ndarray

eml_sigma_inst

The instrumental dispersion (km/s) at the rest wavelength of each emission line. This is mostly used to aid the velocity dispersion corrections determined by Sasuke. Shape is \((N_{\rm eml},)\).

Type:

numpy.ndarray

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

_parse_emission_line_database()[source]

Parse the input emission-line database; see mangadap.par.emissionlinedb.EmissionLinePar.

Only lines with action=f are included in any template. The array tpli provides the index of the template that contains each line in the emission-line database. Lines that are not assigned to any template — either because they do not have action=f or their center lies outside the wavelength range in wave — are given an index of -1.

Only lines with mode=a (i.e., tie flux, velocity, and velocity dispersion) are included in the same template.

Lines with tied velocities are assigned the same velocity component (vgrp) and lines with the tied velocity dispersions are assigned the same sigma component (sgrp).

Warning

The construction of templates for use with Sasuke does not allow one to tie fluxes while leaving the velocities and/or velocity dispersions as independent.

This is an entirely internal procedure, taking no arguments and only assigning results to self.

Raises:

ValueError – Raised if the parsing of the database leaves some templates without an assigned component, velocity group, and/or sigma group. This implies an error in the construction of the emission-line database, not the code itself.

build_templates(emldb, flux_density=True, profile='FFTGaussianLSF', loggers=None, quiet=False)[source]

Build the set of templates for a given emission-line database. The function uses the current values in wave and sigma_inst. Any existing templates from a previous call to build_templates() or from the object instantiation will be overwritten using the provided emission-line database.

See check_database() for the requirements of the provided emission-line database, and see _parse_emission_line_database() for how the database is interpretted when constructing the templates.

The function constructs and returns the following attributes: flux, comp, vgrp, and sgrp

Todo

  • Warn the user if any line is undersampled; i.e., the FWHM of the line is less than 2.1 or \(\sigma < 0.9\).

  • Warn the user if any line grouped in the same template falls outside the spectral range.

Parameters:
  • emldb (mangadap.par.emissionlinedb.EmissionLineDB) – Emission-line database.

  • flux_density (bool, optional) – Return spectrum in units of flux density (flux per angstrom). Default is to return the spectrum in units of flux per pixel.

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

  • quiet (bool, optional) – Suppress all terminal and logging output.

Returns:

Returns 6 numpy.ndarray objects: (1) the set of templates with shape \(N_{\rm tpl}\times N_{\rm wave}\), (2) the kinematic component assignement for each template, (3) the velocity group associated with each template, (4) the sigma group associated with each template, (5) the A matrix used to impose constraints on the kinematics, and (6) the b vector used to impose constraints on the kinematics; for the latter two objects, see the constr_kinem keyword argument for pPXF.

Return type:

tuple

static fill_guess_kinematics(guess_kin, ncomp, lb, ub, A_ineq)[source]

Construct the array with the initial guess kinematics, ensuring that they meet any inequality constraints.

Warning

Always assumes that only 2 moments are fit for each emission-line component.

Parameters:
  • guess_kin (array-like) – Guess kinematics for one or more spectra to be fit. Shape is \((N_{\rm spec}, 2)\) or \((2,)\). The number of parameters will be expanded to the number of components in this template set and adjusted to meet any inequality constraints.

  • ncomp (int) – Number of kinematic components.

  • lb (numpy.ndarray) – Array providing the fractional lower bounds on all the components. Can be None if no lower bounds are set; otherwise, shape must be \((N_{\rm comp}, 2)\).

  • ub (numpy.ndarray) – Array providing the fractional upper bounds on all the components. Can be None if no upper bounds are set; otherwise, shape must be \((N_{\rm comp}, 2)\).

  • A_ineq (numpy.ndarray) – Inequality constraint matrix used by pPXF. Can be None if no constraints are to be applied; otherwise, shape must be \((N_{\rm constr}, 2 N_{\rm comp})\).

Returns:

A 2D or 3D array with the guess parameters. Shape is \((N_{\rm spec}, N_{\rm comp}, 2)\) or \((N_{\rm comp}, 2)\).

Return type:

numpy.ndarray