prospect.sources#

Classes in the prospect.sources module are used to instantiate sps objects. They are defined by the presence of a get_spectrum() method that takes a wavelength array, a list of filter objects, and a parameter dictionary and return a spectrum, a set of broadband fluxes, and a blob of ancillary information.

Most of these classes are a wrapper on fsps.StellarPopulation objects, and as such have a significant memory footprint. The parameter dictionary can include any fsps parameter, as well as parameters used by these classes to control redshifting, spectral smoothing, wavelength calibration, and other aspects of the model.

class prospect.sources.CSPSpecBasis(zcontinuous=1, reserved_params=['sigma_smooth'], vactoair_flag=False, compute_vega_mags=False, **kwargs)#

A subclass of SSPBasis for combinations of N composite stellar populations (including single-age populations). The number of composite stellar populations is given by the length of the "mass" parameter. Other population properties can also be vectors of the same length as "mass" if they are independent for each component.

update(**params)#

Update the params attribute, making parameters scalar if possible.

update_component(component_index)#

Pass params that correspond to a single component through to the fsps.StellarPopulation object.

Parameters:

component_index – The index of the component for which to pull out individual parameters that are passed to the fsps.StellarPopulation object.

get_galaxy_spectrum(**params)#

Update parameters, then loop over each component getting a spectrum for each and sum with appropriate weights.

Parameters:

params – A parameter dictionary that gets passed to the self.update method and will generally include physical parameters that control the stellar population and output spectrum or SED.

Returns wave:

Wavelength in angstroms.

Returns spectrum:

Spectrum in units of Lsun/Hz/solar masses formed.

Returns mass_fraction:

Fraction of the formed stellar mass that still exists.

class prospect.sources.SSPBasis(zcontinuous=1, reserved_params=['tage', 'sigma_smooth'], interp_type='logarithmic', flux_interp='linear', mint_log=-3, compute_vega_mags=False, **kwargs)#

This is a class that wraps the fsps.StellarPopulation object, which is used for producing SSPs. The fsps.StellarPopulation object is accessed as SSPBasis().ssp.

This class allows for the custom calculation of relative SSP weights (by overriding all_ssp_weights) to produce spectra from arbitrary composite SFHs. Alternatively, the entire get_galaxy_spectrum method can be overridden to produce a galaxy spectrum in some other way, for example taking advantage of weight calculations within FSPS for tabular SFHs or for parameteric SFHs.

The base implementation here produces an SSP interpolated to the age given by tage, with initial mass given by mass. However, this is much slower than letting FSPS calculate the weights, as implemented in FastSSPBasis.

Furthermore, smoothing, redshifting, and filter projections are handled outside of FSPS, allowing for fast and more flexible algorithms.

Parameters:

reserved_params – These are parameters which have names like the FSPS parameters but will not be passed to the StellarPopulation object because we are overriding their functionality using (hopefully more efficient) custom algorithms.

update(**params)#

Update the parameters, passing the unreserved FSPS parameters through to the fsps.StellarPopulation object.

Parameters:

params – A parameter dictionary.

get_galaxy_spectrum(**params)#

Update parameters, then multiply SSP weights by SSP spectra and stellar masses, and sum.

Returns wave:

Wavelength in angstroms.

Returns spectrum:

Spectrum in units of Lsun/Hz/solar masses formed.

Returns mass_fraction:

Fraction of the formed stellar mass that still exists.

get_galaxy_elines()#

Get the wavelengths and specific emission line luminosity of the nebular emission lines predicted by FSPS. These lines are in units of Lsun/solar mass formed. This assumes that get_galaxy_spectrum has already been called.

Returns ewave:

The restframe wavelengths of the emission lines, AA

Returns elum:

Specific luminosities of the nebular emission lines, Lsun/stellar mass formed

get_spectrum(outwave=None, filters=None, peraa=False, **params)#

Get a spectrum and SED for the given params.

Parameters:
  • outwave – (default: None) Desired vacuum wavelengths. Defaults to the values in sps.wavelength.

  • peraa – (default: False) If True, return the spectrum in erg/s/cm^2/AA instead of AB maggies.

  • filters – (default: None) A list of filter objects for which you’d like photometry to be calculated.

  • params – Optional keywords giving parameter values that will be used to generate the predicted spectrum.

Returns spec:

Observed frame spectrum in AB maggies, unless peraa=True in which case the units are erg/s/cm^2/AA.

Returns phot:

Observed frame photometry in AB maggies.

Returns mass_frac:

The ratio of the surviving stellar mass to the total mass formed.

property all_ssp_weights#

Weights for a single age population. This is a slow way to do this!

class prospect.sources.FastStepBasis(zcontinuous=1, reserved_params=['tage', 'sigma_smooth'], interp_type='logarithmic', flux_interp='linear', mint_log=-3, compute_vega_mags=False, **kwargs)#

Subclass of SSPBasis that implements a “nonparameteric” (i.e. binned) SFH. This is accomplished by generating a tabular SFH with the proper form to be passed to FSPS. The key parameters for this SFH are:

  • agebins - array of shape (nbin, 2) giving the younger and older (in lookback time) edges of each bin in log10(years)

  • mass - array of shape (nbin,) giving the total stellar mass (in solar masses) formed in each bin.

get_galaxy_spectrum(**params)#

Construct the tabular SFH and feed it to the ssp.

convert_sfh(agebins, mformed, epsilon=0.0001, maxage=None)#

Given arrays of agebins and formed masses with each bin, calculate a tabular SFH. The resulting time vector has time points either side of each bin edge with a “closeness” defined by a parameter epsilon.

Parameters:
  • agebins – An array of bin edges, log(yrs). This method assumes that the upper edge of one bin is the same as the lower edge of another bin. ndarray of shape (nbin, 2)

  • mformed – The stellar mass formed in each bin. ndarray of shape (nbin,)

  • epsilon – (optional, default 1e-4) A small number used to define the fraction time separation of adjacent points at the bin edges.

  • maxage – (optional, default: None) A maximum age of stars in the population, in yrs. If None then the maximum value of agebins is used. Note that an error will occur if maxage < the maximum age in agebins.

Returns time:

The output time array for use with sfh=3, in Gyr. ndarray of shape (2*N)

Returns sfr:

The output sfr array for use with sfh=3, in M_sun/yr. ndarray of shape (2*N)

Returns maxage:

The maximum valid age in the returned isochrone.

class prospect.sources.BlackBodyDustBasis(**kwargs)#
get_spectrum(outwave=None, filters=None, **params)#

Given a params dictionary, generate spectroscopy, photometry and any extras (e.g. stellar mass).

Parameters:
  • outwave – The output wavelength vector.

  • filters – A list of sedpy filter objects.

  • params – Keywords forming the parameter set.

Returns spec:

The restframe spectrum in units of erg/s/cm^2/AA

Returns phot:

The apparent (redshifted) maggies in each of the filters.

Returns extras:

A list of None type objects, only included for consistency with the SedModel class.

one_sed(icomp=0, wave=None, filters=None, **extras)#

Pull out individual component parameters from the param dictionary and generate spectra for those components

normalization()#

This method computes the normalization (due do distance dimming, unit conversions, etc.) based on the content of the params dictionary.