prospect.models
Contents
prospect.models#
prospect.models#
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 casetburst
>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 bothzred
andtage
parameters withouttage
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 atzred
.
- 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.