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 vectortheta
.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. -
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.
-
configure
(reset=False, **kwargs)¶ Use the
config_dict
to generate atheta_index
mapping, and propogate the initial parameters into theparams
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
- kwargs – Keyword parameters can be used to override or add to the initial
parameter values specified in
-
fixed_params
¶ A list of the names fixed model parameters that are specified in the
config_dict
.
-
free_params
¶ A list of the names of the free model parameters.
-
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.
-
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
- theta – The parameter vector for which you want to calculate the
prior. ndarray of shape
-
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.
-
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,)
-
theta
¶ The current value of the theta vector, pulled from the
params
state dictionary.
-
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_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,)
-
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,)
- initial_disp – (default: 0.1)
The default dispersion to use in case the
-
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.
-
class
prospect.models.
SpecModel
(configuration, verbose=True, param_order=None, **kwargs)¶ A subclass of
ProspectorParams
that passes the models through to ansps
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.
-
cache_eline_parameters
(obs, nsigma=5)¶ 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_elines_to_fit
- If fitting and marginalizing over emission lines, this stores indices of the lines to actually fit, as a boolean array. Only lines that are withinnsigma
of an observed wavelength points are included._eline_wavelength_mask
- A mask of the _outwave vector that indicates which pixels to use in the emission line fitting. Only pixels withinnsigma
of an emission line are used.
Can be subclassed to add more sophistication redshift - first looks for
eline_delta_zred
, and defaults tozred
sigma - first looks foreline_sigma
, defaults to 100 km/sParameters: 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.
-
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.
-
get_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 observer-frame spectrum in the same units as the
observed spectrum, ndarray of shape
(n_wave,)
- 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)
wheren_wave_neb
is the number of wavelength elements withinnsigma
of a line, andn_fitted_lines
is the number of lines that fall withinnsigma
of a wavelength pixel. Units are same ascalibrated_spec
- obs – A dictionary containing the
-
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)
-
get_eline_spec
(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: 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
-
mean_model
(theta, obs, sps=None, sigma=None, **extras)¶ Legacy wrapper around predict()
-
nebline_photometry
(filters)¶ Compute the emission line contribution to photometry. This requires several cached attributes:
_ewave_obs
_eline_lum
Parameters: filters – List of sedpy.observate.Filter
objectsReturns nebflux: The flux of the emission line through the filters, in units of maggies. ndarray of shape (len(filters),)
-
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
-
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 maggiesReturns 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.
- theta – ndarray of parameter values, of shape
-
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._eline_wave
and_eline_lum
- emission line parameters from the SPS model
Parameters: filters – List of sedpy.observate.Filter
objects. If there is no photometry,None
should be suppliedReturns phot: Observed frame photometry of the model SED through the given filters. ndarray of shape (len(filters),)
, in units of maggies. Iffilters
is None, this returns 0.0
-
predict_spec
(obs, sigma_spec, **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_elinespec
- emission line spectrum_sed
- Intrinsic spectrum (before cilbration vector applied)
And if emission line marginalization is being performed, numerous quantities related to the emission lines are also cached (see
get_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.
-
smoothspec
(wave, spec)¶ Smooth the spectrum. See
prospect.utils.smoothing.smoothspec()
for details.
-
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
-
-
class
prospect.models.
SedModel
(configuration, verbose=True, param_order=None, **kwargs)¶ A subclass of
ProspectorParams
that passes the models through to ansps
object and returns spectra and photometry, including optional spectroscopic calibration and sky emission.-
mean_model
(theta, obs, sps=None, sigma_spec=None, **extras)¶ Legacy wrapper around predict()
-
predict
(theta, obs=None, sps=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_spectrum()
method defined. - sigma_spec – (optional, unused) 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 maggiesReturns 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.
- theta – ndarray of parameter values, of shape
-
sed
(theta, obs=None, sps=None, **kwargs)¶ Given a vector of parameters
theta
, generate a spectrum, photometry, and any extras (e.g. surviving mass fraction), *not including any instrument calibration effects. The intrinsic spectrum thus produced is cached in _spec attributeParameters: - theta – ndarray of parameter values.
- 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_spectrum()
method defined.
Returns spec: The model spectrum for these parameters, at the wavelengths specified by
obs['wavelength']
. Default units are maggies, and the calibration vector is not applied.Returns phot: The model photometry for these parameters, for the filters specified in
obs['filters']
. Units are 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 steallr mass formed.
-
sky
(obs)¶ Model for the additive sky emission/absorption
-
spec_calibration
(theta=None, obs=None, **kwargs)¶ Implements an overall scaling of the spectrum, given by the parameter
'spec_norm'
Returns cal: (float) A scalar multiplicative factor that gives the ratio between the true spectrum and the observed spectrum
-
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
Parameters: - wavelength – The input wavelengths. ndarray of shape
(nwave,)
- mask – optional
The mask. slice or boolean array with
True
for unmasked elements. The interval (-1, 1) will be defined only by unmasked wavelength points
Returns x: The wavelength vector, remapped to the interval (-1, 1). ndarray of same shape as
wavelength
- wavelength – The input wavelengths. ndarray of shape
-
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()(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 – 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…. -
inverse_unit_transform
(x, **kwargs)¶ Go from the parameter value to the unit coordinate using the cdf.
-
loc
¶ This should be overridden.
-
sample
(nsample=None, **kwargs)¶ Draw a sample from the prior distribution.
Parameters: nsample – (optional) Unused
-
scale
¶ This should be overridden.
-
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.
-
update
(**kwargs)¶ Update params values using alias.
-
-
class
prospect.models.priors.
Uniform
(parnames=[], name='', **kwargs)¶ A simple uniform prior, described by two parameters
Parameters: - mini – Minimum of the distribution
- maxi – Maximum of the distribution
-
distribution
= <scipy.stats._continuous_distns.uniform_gen object>¶
-
loc
¶ This should be overridden.
-
scale
¶ This should be overridden.
-
class
prospect.models.priors.
TopHat
(parnames=[], name='', **kwargs)¶ Uniform distribution between two bounds, renamed for backwards compatibility :param mini:
Minimum of the distributionParameters: maxi – Maximum of the distribution
-
class
prospect.models.priors.
Normal
(parnames=[], name='', **kwargs)¶ A simple gaussian prior.
Parameters: - mean – Mean of the distribution
- sigma – Standard deviation of the distribution
-
distribution
= <scipy.stats._continuous_distns.norm_gen object>¶
-
loc
¶ This should be overridden.
-
scale
¶ This should be overridden.
-
class
prospect.models.priors.
ClippedNormal
(parnames=[], name='', **kwargs)¶ A Gaussian prior clipped to some range.
Parameters: - mean – Mean of the normal distribution
- sigma – Standard deviation of the normal distribution
- mini – Minimum of the distribution
- maxi – Maximum of the distribution
-
distribution
= <scipy.stats._continuous_distns.truncnorm_gen object>¶
-
loc
¶ This should be overridden.
-
scale
¶ This should be overridden.
-
class
prospect.models.priors.
LogNormal
(parnames=[], name='', **kwargs)¶ A log-normal prior, where the natural log of the variable is distributed normally. Useful for parameters that cannot be less than zero.
Note that
LogNormal(np.exp(mode) / f) == LogNormal(np.exp(mode) * f)
andf = np.exp(sigma)
corresponds to “one sigma” from the peak.Parameters: - mode – Natural log of the variable value at which the probability density is highest.
- sigma – Standard deviation of the distribution of the natural log of the variable.
-
distribution
= <scipy.stats._continuous_distns.lognorm_gen object>¶
-
loc
¶ This should be overridden.
-
scale
¶ This should be overridden.
-
class
prospect.models.priors.
LogUniform
(parnames=[], name='', **kwargs)¶ Like log-normal, but the distribution of natural log of the variable is distributed uniformly instead of normally.
Parameters: - mini – Minimum of the distribution
- maxi – Maximum of the distribution
-
distribution
= <scipy.stats._continuous_distns.reciprocal_gen object>¶
-
class
prospect.models.priors.
Beta
(parnames=[], name='', **kwargs)¶ A Beta distribution.
Parameters: - mini – Minimum of the distribution
- maxi – Maximum of the distribution
- alpha –
- beta –
-
distribution
= <scipy.stats._continuous_distns.beta_gen object>¶
-
loc
¶ This should be overridden.
-
scale
¶ This should be overridden.
-
class
prospect.models.priors.
StudentT
(parnames=[], name='', **kwargs)¶ A Student’s T distribution
Parameters: - mean – Mean of the distribution
- scale – Size of the distribution, analogous to the standard deviation
- df – Number of degrees of freedom
-
distribution
= <scipy.stats._continuous_distns.t_gen object>¶
-
loc
¶ This should be overridden.
-
scale
¶ This should be overridden.
-
class
prospect.models.priors.
SkewNormal
(parnames=[], name='', **kwargs)¶ A normal distribution including a skew parameter
Parameters: - location – Center (not mean, mode, or median) of the distribution. The center will approach the mean as skew approaches zero.
- sigma – Standard deviation of the distribution
- skew – Skewness of the distribution
-
distribution
= <scipy.stats._continuous_distns.skew_norm_gen object>¶
-
loc
¶ This should be overridden.
-
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 – FSPS stellar metaliicity parameter. Returns logzsol: The same.
-
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 – The log10(mass) Returns mass: The mass in linear units
-
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 casetburst
>tage
.Parameters: - tage – The age of the host galaxy (Gyr)
- fage_burst – 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)
-
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 bothzred
andtage
parameters withouttage
exceeding the age of the universe.Parameters: - zred – Cosmological redshift.
- tage_tuniv – The ratio of
tage
to the age of the universe atzred
.
Returns tage: The stellar population age, in Gyr
-
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 – Cosmological redshift. This sets the age of the universe.
- agebins – The SFH bin edges in log10(years). ndarray of shape
(nbin, 2)
.
Returns agebins: The new SFH bin edges.
-
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 – The diffuse dust V-band optical depth (the FSPS
dust2
parameter.) - dust_ratio – 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.)- dust2 – The diffuse dust V-band optical depth (the FSPS
-
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
(logsfr_ratios=None, agebins=None, **extras)¶ This transforms from SFR ratios to agebins by assuming a constant amount of mass forms in each bin agebins = np.array([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 – The total mass formed over all bins in the SFH.
- z_fraction – latent variables drawn form a specific set of Beta distributions. (see Betancourt 2010)
Returns masses: The stellar mass formed in each age bin.
-
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 – latent variables drawn form a specific set of Beta distributions. (see Betancourt 2010) Returns sfrac: The star formation fractions (See Leja et al. 2017 for definition).
-
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: The total mass Returns zfrac: The dimensionless z variables used for sfr fraction parameterization.