prospect.models

Contents

prospect.models#

prospect.models#

This module includes objects that store parameter specfications and efficiently convert between parameter dictionaries and parameter vectors necessary for fitting algorithms. There are submodules for parameter priors, common parameter transformations, and pre-defined sets of parameter specifications.

class prospect.models.ProspectorParams(configuration, verbose=True, param_order=None, **kwargs)#

This is the base model class that holds model parameters and information about them (e.g. priors, bounds, transforms, free vs fixed state). In addition to the documented methods, it contains several important attributes:

  • params: model parameter state dictionary.

  • theta_index: A dictionary that maps parameter names to indices (or rather slices) of the parameter vector theta.

  • config_dict: Information about each parameter as a dictionary keyed by parameter name for easy access.

  • config_list: Information about each parameter stored as a list.

Intitialization is via, e.g.,

model_dict = {"mass": {"N": 1, "isfree": False, "init": 1e10}}
model = ProspectorParams(model_dict, param_order=None)
Parameters:

configuration – A list or dictionary of model parameters specifications.

configure(reset=False, **kwargs)#

Use the config_dict to generate a theta_index mapping, and propogate the initial parameters into the params state dictionary, and store the intital theta vector thus implied.

Parameters:
  • kwargs – Keyword parameters can be used to override or add to the initial parameter values specified in config_dict

  • reset – (default: False) If true, empty the params dictionary before re-reading the config_dict

map_theta()#

Construct the mapping from parameter name to the index in the theta vector corresponding to the first element of that parameter. Called during configuration.

set_parameters(theta)#

Propagate theta into the model parameters params dictionary.

Parameters:

theta – A theta parameter vector containing the desired parameters. ndarray of shape (ndim,)

prior_product(theta, nested=False, **extras)#

Public version of _prior_product to be overridden by subclasses.

Parameters:
  • theta – The parameter vector for which you want to calculate the prior. ndarray of shape (..., ndim)

  • nested – If using nested sampling, this will only return 0 (or -inf). This behavior can be overridden if you want to include complicated priors that are not included in the unit prior cube based proposals (e.g. something that is difficult to transform from the unit cube.)

Returns lnp_prior:

The natural log of the prior probability at theta

prior_transform(unit_coords)#

Go from unit cube to parameter space, for nested sampling.

Parameters:

unit_coords – Coordinates in the unit hyper-cube. ndarray of shape (ndim,).

Returns theta:

The parameter vector corresponding to the location in prior CDF corresponding to unit_coords. ndarray of shape (ndim,)

propagate_parameter_dependencies()#

Propogate any parameter dependecies. That is, for parameters whose value depends on another parameter, calculate those values and store them in the self.params dictionary.

rectify_theta(theta, epsilon=1e-10)#

Replace zeros in a given theta vector with a small number epsilon.

property theta#

The current value of the theta vector, pulled from the params state dictionary.

property free_params#

A list of the names of the free model parameters.

property fixed_params#

A list of the names fixed model parameters that are specified in the config_dict.

theta_labels(name_map={})#

Using the theta_index parameter map, return a list of the model parameter names that has the same order as the sampling chain array.

Parameters:

name_map – A dictionary mapping model parameter names to output label names.

Returns labels:

A list of labels of the same length and order as the theta vector.

theta_bounds()#

Get the bounds on each parameter from the prior.

Returns bounds:

A list of length ndim of tuples (lo, hi) giving the parameter bounds.

theta_disps(default_disp=0.1, fractional_disp=False)#

Get a vector of absolute dispersions for each parameter to use in generating sampler balls for emcee’s Ensemble sampler. This can be overridden by subclasses if fractional dispersions are desired.

Parameters:
  • initial_disp – (default: 0.1) The default dispersion to use in case the "init_disp" key is not provided in the parameter configuration.

  • fractional_disp – (default: False) Treat the dispersion values as fractional dispersions.

Returns disp:

The dispersion in the parameters to use for generating clouds of walkers (or minimizers.) ndarray of shape (ndim,)

theta_disp_floor(thetas=None)#

Get a vector of dispersions for each parameter to use as a floor for the emcee walker-calculated dispersions. This can be overridden by subclasses.

Returns disp_floor:

The minimum dispersion in the parameters to use for generating clouds of walkers (or minimizers.) ndarray of shape (ndim,)

clip_to_bounds(thetas)#

Clip a set of parameters theta to within the priors.

Parameters:

thetas – The parameter vector, ndarray of shape (ndim,).

Returns thetas:

The input vector, clipped to the bounds of the priors.

class prospect.models.SpecModel(*args, **kwargs)#

A subclass of ProspectorParams that passes the models through to an sps object and returns spectra and photometry, including optional spectroscopic calibration, and sky emission.

This class performs most of the conversion from intrinsic model spectrum to observed quantities, and additionally can compute MAP emission line values and penalties for marginalization over emission line amplitudes.

predict(theta, obs=None, sps=None, sigma_spec=None, **extras)#

Given a theta vector, generate a spectrum, photometry, and any extras (e.g. stellar mass), including any calibration effects.

Parameters:
  • theta – ndarray of parameter values, of shape (ndim,)

  • obs – An observation dictionary, containing the output wavelength array, the photometric filter lists, and the observed fluxes and uncertainties thereon. Assumed to be the result of utils.obsutils.rectify_obs()

  • sps – An sps object to be used in the model generation. It must have the get_galaxy_spectrum() method defined.

  • sigma_spec – (optional) The covariance matrix for the spectral noise. It is only used for emission line marginalization.

Returns spec:

The model spectrum for these parameters, at the wavelengths specified by obs['wavelength'], including multiplication by the calibration vector. Units of maggies

Returns phot:

The model photometry for these parameters, for the filters specified in obs['filters']. Units of maggies.

Returns extras:

Any extra aspects of the model that are returned. Typically this will be mfrac the ratio of the surviving stellar mass to the stellar mass formed.

predict_spec(obs, sigma_spec=None, **extras)#

Generate a prediction for the observed spectrum. This method assumes that the parameters have been set and that the following attributes are present and correct

  • _wave - The SPS restframe wavelength array

  • _zred - Redshift

  • _norm_spec - Observed frame spectral fluxes, in units of maggies

  • _eline_wave and _eline_lum - emission line parameters from the SPS model

It generates the following attributes

  • _outwave - Wavelength grid (observed frame)

  • _speccal - Calibration vector

  • _sed - Intrinsic spectrum (before cilbration vector applied but including emission lines)

And the following attributes are generated if nebular lines are added

  • _fix_eline_spec - emission line spectrum for fixed lines, intrinsic units

  • _fix_eline_spec - emission line spectrum for fitted lines, with spectroscopic calibration factor included.

Numerous quantities related to the emission lines are also cached (see cache_eline_parameters() and fit_el() for details.)

Parameters:
  • obs – An observation dictionary, containing the output wavelength array, the photometric filter lists, and the observed fluxes and uncertainties thereon. Assumed to be the result of utils.obsutils.rectify_obs()

  • sigma_spec – (optional) The covariance matrix for the spectral noise. It is only used for emission line marginalization.

Returns spec:

The prediction for the observed frame spectral flux these parameters, at the wavelengths specified by obs['wavelength'], including multiplication by the calibration vector. ndarray of shape (nwave,) in units of maggies.

predict_phot(filters)#

Generate a prediction for the observed photometry. This method assumes that the parameters have been set and that the following attributes are present and correct:

  • _wave - The SPS restframe wavelength array

  • _zred - Redshift

  • _norm_spec - Observed frame spectral fluxes, in units of maggies.

  • _ewave_obs and _eline_lum - emission line parameters from the SPS model

Parameters:

filters – Instance of sedpy.observate.FilterSet or list of sedpy.observate.Filter objects. If there is no photometry, None should be supplied.

Returns phot:

Observed frame photometry of the model SED through the given filters. ndarray of shape (len(filters),), in units of maggies. If filters is None, this returns 0.0

flux_norm()#

Compute the scaling required to go from Lsun/Hz/Msun to maggies. Note this includes the (1+z) factor required for flux densities.

Returns norm:

(float) The normalization factor, scalar float.

nebline_photometry(filters, elams=None, elums=None)#

Compute the emission line contribution to photometry. This requires several cached attributes:

  • _ewave_obs

  • _eline_lum

Parameters:
  • filters – Instance of sedpy.observate.FilterSet or list of sedpy.observate.Filter objects

  • elams – (optional) The emission line wavelength in angstroms. If not supplied uses the cached _ewave_obs attribute.

  • elums – (optional) The emission line flux in erg/s/cm^2. If not supplied uses the cached _eline_lum attribute and applies appropriate distance dimming and unit conversion.

Returns nebflux:

The flux of the emission line through the filters, in units of maggies. ndarray of shape (len(filters),)

cache_eline_parameters(obs, nsigma=5, forcelines=False)#

This computes and caches a number of quantities that are relevant for predicting the emission lines, and computing the MAP values thereof, including

  • _ewave_obs - Observed frame wavelengths (AA) of all emission lines.

  • _eline_sigma_kms - Dispersion (in km/s) of all the emission lines

  • _fit_eline - If fitting and marginalizing over emission lines, this stores a boolean mask of the lines to actually fit. Only lines that are within nsigma of an observed wavelength points are included.

  • _fix_eline - this stores a boolean mask of the lines that are to be added with the cloudy amplitudes Only lines that are within nsigma of an observed wavelength point are included.

  • _fit_eline_pixelmask - A mask of the _outwave vector that indicates which pixels to use in the emission line fitting. Only pixels within nsigma of an emission line are used.

  • _fix_eline_pixelmask - A mask of the _outwave vector that indicates which pixels to use in the fixed emission line prediction.

Can be subclassed to add more sophistication redshift - first looks for eline_delta_zred, and defaults to zred sigma - first looks for eline_sigma, defaults to 100 km/s

Parameters:

nsigma – (float, optional, default: 5.) Number of sigma from a line center to use for defining which lines to fit and useful spectral elements for the fitting. float.

parse_elines()#

Create mask arrays to identify the lines that should be fit and the lines that should be fixed based on the content of params[“elines_to_fix”] and params[“elines_to_fit”]

This can probably be cached just once unless you want to change between separate fits.

fit_el(obs, calibrated_spec, sigma_spec=None)#

Compute the maximum likelihood and, optionally, MAP emission line amplitudes for lines that fall within the observed spectral range. Also compute and cache the analytic penalty to log-likelihood from marginalizing over the emission line amplitudes. This is cached as _ln_eline_penalty. The emission line amplitudes (in maggies) at _eline_lums are updated to the ML values for the fitted lines.

Parameters:
  • obs – A dictionary containing the 'spectrum' and 'unc' keys that are observed fluxes and uncertainties, both ndarrays of shape (n_wave,)

  • calibrated_spec – The predicted (so far) observer-frame spectrum in the same units as the observed spectrum, ndarray of shape (n_wave,) Should include pixed lines but not lines to be fit

  • sigma_spec – Spectral covariance matrix, if using a non-trivial noise model.

Returns el:

The maximum likelihood emission line flux densities. ndarray of shape (n_wave_neb, n_fitted_lines) where n_wave_neb is the number of wavelength elements within nsigma of a line, and n_fitted_lines is the number of lines that fall within nsigma of a wavelength pixel. Units are same as calibrated_spec

predict_eline_spec(line_indices=slice(None, None, None), wave=None)#

Compute a complete model emission line spectrum. This should only be run after calling predict(), as it accesses cached information. Relatively slow, useful for display purposes

Parameters:
  • line_indices – optional If given, this should give the indices of the lines to predict.

  • wave – (optional, default: None) The wavelength ndarray on which to compute the emission line spectrum. If not supplied, the _outwave vector is used.

Returns eline_spec:

An (n_line, n_wave) ndarray, units of Lsun/Hz intrinsic (no calibration vector applied)

get_eline_gaussians(lineidx=slice(None, None, None), wave=None)#

Generate a set of unit normals with centers and widths given by the previously cached emission line observed-frame wavelengths and emission line widths.

Parameters:
  • lineidx – (optional) A boolean array or integer array used to subscript the cached lines. Gaussian vectors will only be constructed for the lines thus subscripted.

  • wave – (optional) The wavelength array (in Angstroms) used to construct the gaussian vectors. If not given, the cached _outwave array will be used.

Returns gaussians:

The unit gaussians for each line, in units Lsun/Hz. ndarray of shape (n_wave, n_line)

smoothspec(wave, spec)#

Smooth the spectrum. See prospect.utils.smoothing.smoothspec() for details.

observed_wave(wave, do_wavecal=False)#

Convert the restframe wavelngth grid to the observed frame wavelength grid, optionally including wavelength calibration adjustments. Requires that the _zred attribute is already set.

Parameters:

wave – The wavelength array

wave_to_x(wavelength=None, mask=slice(None, None, None), **extras)#

Map unmasked wavelengths to the interval -1, 1 masked wavelengths may have x>1, x<-1

absolute_rest_maggies(filters)#

Return absolute rest-frame maggies (=10**(-0.4*M)) of the last computed spectrum.

Parameters:

filters (list of sedpy.observate.Filter() instances) – The filters through which you wish to compute the absolute mags

Returns:

maggies – The absolute restframe maggies of the model through the supplied filters, including emission lines. Convert to absolute rest-frame magnitudes as M = -2.5 * log10(maggies)

Return type:

ndarray of shape (nbands,)

mean_model(theta, obs, sps=None, sigma=None, **extras)#

Legacy wrapper around predict()

class prospect.models.PolySpecModel(*args, **kwargs)#

This is a subclass of SpecModel that generates the multiplicative calibration vector at each model predict call as the maximum likelihood chebyshev polynomial describing the ratio between the observed and the model spectrum.

spec_calibration(theta=None, obs=None, spec=None, **kwargs)#

Implements a Chebyshev polynomial calibration model. This uses least-squares to find the maximum-likelihood Chebyshev polynomial of a certain order describing the ratio of the observed spectrum to the model spectrum, conditional on all other parameters. If emission lines are being marginalized out, they are excluded from the least-squares fit.

Returns cal:

A polynomial given by \(\sum_{m=0}^M a_{m} * T_m(x)\).

class prospect.models.SplineSpecModel(*args, **kwargs)#

This is a subclass of SpecModel that generates the multiplicative calibration vector at each model predict call as the maximum likelihood cubic spline describing the ratio between the observed and the model spectrum.

spec_calibration(theta=None, obs=None, spec=None, **kwargs)#

Implements a spline calibration model. This fits a cubic spline with determined knot locations to the ratio of the observed spectrum to the model spectrum. If emission lines are being marginalized out, they are excluded from the least-squares fit.

The knot locations must be specified as model parameters, either explicitly or via a number of knots or knot spacing (in angstroms)

make_knots(wlo, whi, as_x=True, **params)#

prospect.models.priors#

priors.py – This module contains various objects to be used as priors. When called these return the ln-prior-probability, and they can also be used to construct prior transforms (for nested sampling) and can be sampled from.

class prospect.models.priors.Prior(parnames=[], name='', **kwargs)#

Encapsulate the priors in an object. Each prior should have a distribution name and optional parameters specifying scale and location (e.g. min/max or mean/sigma). These can be aliased at instantiation using the parnames keyword. When called, the argument should be a variable and the object should return the ln-prior-probability of that value.

ln_prior_prob = Prior(param=par)(value)

Should be able to sample from the prior, and to get the gradient of the prior at any variable value. Methods should also be avilable to give a useful plotting range and, if there are bounds, to return them.

Parameters:

parnames (sequence of strings) – A list of names of the parameters, used to alias the intrinsic parameter names. This way different instances of the same Prior can have different parameter names, in case they are being fit for….

params#

The values of the parameters describing the prior distribution.

Type:

dictionary

update(**kwargs)#

Update self.params values using alias.

sample(nsample=None, **kwargs)#

Draw a sample from the prior distribution.

Parameters:

nsample – (optional) Unused

unit_transform(x, **kwargs)#

Go from a value of the CDF (between 0 and 1) to the corresponding parameter value.

Parameters:

x – A scalar or vector of same length as the Prior with values between zero and one corresponding to the value of the CDF.

Returns theta:

The parameter value corresponding to the value of the CDF given by x.

inverse_unit_transform(x, **kwargs)#

Go from the parameter value to the unit coordinate using the cdf.

property loc#

This should be overridden.

property scale#

This should be overridden.

prospect.models.transforms#

transforms.py – This module contains parameter transformations that may be useful to transform from parameters that are easier to _sample_ in to the parameters required for building SED models.

They can be used as "depends_on" entries in parameter specifications.

prospect.models.transforms.stellar_logzsol(logzsol=0.0, **extras)#

Simple function that takes an argument list and returns the value of the logzsol argument (i.e. the stellar metallicity)

Parameters:

logzsol (float) – FSPS stellar metaliicity parameter.

Returns:

logzsol – The same.

Return type:

float

prospect.models.transforms.delogify_mass(logmass=0.0, **extras)#

Simple function that takes an argument list including a logmass parameter and returns the corresponding linear mass.

Parameters:

logmass (float) – The log10(mass)

Returns:

mass – The mass in linear units

Return type:

float

prospect.models.transforms.tburst_from_fage(tage=0.0, fage_burst=0.0, **extras)#

This function transfroms from a fractional age of a burst to an absolute age. With this transformation one can sample in fage_burst without worry about the case tburst > tage.

Parameters:
  • tage (float, Gyr) – The age of the host galaxy.

  • fage_burst (float between 0 and 1) – The fraction of the host age at which the burst occurred.

Returns:

tburst – The age of the host when the burst occurred (i.e. the FSPS tburst parameter)

Return type:

float, Gyr

prospect.models.transforms.tage_from_tuniv(zred=0.0, tage_tuniv=1.0, **extras)#

This function calculates a galaxy age from the age of the universe at zred and the age given as a fraction of the age of the universe. This allows for both zred and tage parameters without tage exceeding the age of the universe.

Parameters:
  • zred (float) – Cosmological redshift.

  • tage_tuniv (float between 0 and 1) – The ratio of tage to the age of the universe at zred.

Returns:

tage – The stellar population age, in Gyr

Return type:

float

prospect.models.transforms.zred_to_agebins(zred=0.0, agebins=[], **extras)#

Set the nonparameteric SFH age bins depending on the age of the universe at zred. The first bin is not altered and the last bin is always 15% of the upper edge of the oldest bin, but the intervening bins are evenly spaced in log(age).

Parameters:
  • zred (float) – Cosmological redshift. This sets the age of the universe.

  • agebins (ndarray of shape (nbin, 2)) – The SFH bin edges in log10(years).

Returns:

agebins – The new SFH bin edges.

Return type:

ndarray of shape (nbin, 2)

prospect.models.transforms.dustratio_to_dust1(dust2=0.0, dust_ratio=0.0, **extras)#

Set the value of dust1 from the value of dust2 and dust_ratio

Parameters:
  • dust2 (float) – The diffuse dust V-band optical depth (the FSPS dust2 parameter.)

  • dust_ratio (float) – The ratio of the extra optical depth towards young stars to the diffuse optical depth affecting all stars.

Returns:

dust1 – The extra optical depth towards young stars (the FSPS dust1 parameter.)

Return type:

float

prospect.models.transforms.logsfr_ratios_to_masses(logmass=None, logsfr_ratios=None, agebins=None, **extras)#

This converts from an array of log_10(SFR_j / SFR_{j+1}) and a value of log10(Sum_i M_i) to values of M_i. j=0 is the most recent bin in lookback time.

prospect.models.transforms.logsfr_ratios_to_sfrs(logmass=None, logsfr_ratios=None, agebins=None, **extras)#

Convenience function

prospect.models.transforms.logsfr_ratios_to_agebins([NBINS, 2])#
use equation:

delta(t1) = tuniv / (1 + SUM(n=1 to n=nbins-1) PROD(j=1 to j=n) Sn) where Sn = SFR(n) / SFR(n+1) and delta(t1) is width of youngest bin

prospect.models.transforms.zfrac_to_masses(total_mass=None, z_fraction=None, agebins=None, **extras)#

This transforms from independent dimensionless z variables to sfr fractions and then to bin mass fractions. The transformation is such that sfr fractions are drawn from a Dirichlet prior. See Betancourt et al. 2010 and Leja et al. 2017

Parameters:
  • total_mass (float) – The total mass formed over all bins in the SFH.

  • z_fraction (ndarray of shape (Nbins-1,)) – latent variables drawn from a specific set of Beta distributions. (see Betancourt 2010)

Returns:

masses – The stellar mass formed in each age bin.

Return type:

ndarray of shape (Nbins,)

prospect.models.transforms.zfrac_to_sfrac(z_fraction=None, **extras)#

This transforms from independent dimensionless z variables to sfr fractions. The transformation is such that sfr fractions are drawn from a Dirichlet prior. See Betancourt et al. 2010 and Leja et al. 2017

Parameters:

z_fraction (ndarray of shape (Nbins-1,)) – latent variables drawn from a specific set of Beta distributions. (see Betancourt 2010)

Returns:

sfrac – The star formation fractions (See Leja et al. 2017 for definition).

Return type:

ndarray of shape (Nbins,)

prospect.models.transforms.zfrac_to_sfr(total_mass=None, z_fraction=None, agebins=None, **extras)#

This transforms from independent dimensionless z variables to SFRs.

Returns sfrs:

The SFR in each age bin (msun/yr).

prospect.models.transforms.masses_to_zfrac(mass=None, agebins=None, **extras)#

The inverse of zfrac_to_masses(), for setting mock parameters based on mock bin masses.

Returns:

  • total_mass (float) – The total mass

  • zfrac (ndarray of shape (Nbins-1,)) – latent variables drawn from a specific set of Beta distributions. (see Betancourt 2010) related to the fraction of mass formed in each bin.

prospect.models.transforms.zred_to_agebins_pbeta(zred=None, agebins=[], **extras)#

New agebin scheme, refined so that none of the bins is overly wide when the universe is young.

Parameters:
  • zred (float) – Cosmological redshift. This sets the age of the universe.

  • agebins (ndarray of shape (nbin, 2)) – The SFH bin edges in log10(years).

Returns:

agebins – The new SFH bin edges.

Return type:

ndarray of shape (nbin, 2)