Modelling Stellar Populations

fsps.StellarPopulation drives FSPS and is generally the only API needed. fsps.find_filter() can also be used interactively to search for a filter (see also Table of FSPS Filters).

Example

Lets start by initializing a simple stellar population with solar metallicity and some dust with a Calzetti et al (2000) extinction curve:

>>> import fsps
>>> sp = fsps.StellarPopulation(compute_vega_mags=False, zcontinuous=1,
                                sfh=0, logzsol=0.0, dust_type=2, dust2=0.2)
>>> sp.libraries
('mist', 'miles', 'DL07')

The last line indicates that we are using the MIST isochrones and MILES spectral library. These can be changed only by reinstalling python-FSPS with appropriate flags.

Let’s get the AB magnitudes in SDSS bands, for a simple stellar population (SSP) that is 13.7 Gyr old:

>>> sdss_bands = fsps.find_filter('sdss')
>>> print(sdss_bands)
['sdss_u', 'sdss_g', 'sdss_i', 'sdss_r', 'sdss_z']
>>> sp.get_mags(tage=13.7, bands=sdss_bands)
array([ 9.85484694,  7.8785663 ,  6.56580511,  6.99828574,  6.16779663])

Now we can change a parameter (say, lower the metallicity) and get a new set of magnitudes:

>>> sp.params['logzsol'] = -1
>>> sp.get_mags(tage=13.7, bands=sdss_bands)
array([ 8.5626572 ,  7.07918435,  6.05304881,  6.38592117,  5.84199   ])

We can also get the spectrum, here in units of \(L_\odot/\mathrm{Hz}\), as well as the total stellar mass formed by 13.7 Gyr and the surviving stellar mass at 13.7 Gyr (both in units of \(M_\odot\)):

>>> wave, spec = sp.get_spectrum(tage=13.7)
>>> sp.formed_mass
1.0
>>> sp.stellar_mass
0.57809244144339011

It is highly recommended that only one instance of fsps.StellarPopulation be used in a given program.

Example using nebular emission

We initialize a simple stellar population and set the flag to include nebular emission:

>>> sp = fsps.StellarPopulation(zcontinuous=1,
                                add_neb_emission=1)

We can change the stellar metallicity, the gas phase metallicity, the gas ionization parameter, and then return the total spectrum at 1 Myr:

>>> sp.params['logzsol'] = -1.0
>>> sp.params['gas_logz'] = -1.0
>>> sp.params['gas_logu'] = -2.5
>>> wave, spec = sp.get_spectrum(tage=0.001)

Note: for the nebular model to be fully self-consistent, the gas phase metallicity and the stellar metallicity should be set to the same value. This effectively adds the emission spectrum to the same stellar spectrum that was used as the ionizing spectrum in Cloudy. If users choose to vary the gas phase metallicity at constant stellar metallicity, expect non-hydrogenic emission lines to be accurate within 1-15%.

Emission line wavelengths and line luminosities can be accessed through the stellar population object:

>>> sp.emline_wavelengths
>>> sp.emline_luminosity

API Reference

class fsps.StellarPopulation(compute_vega_mags=False, vactoair_flag=False, zcontinuous=0, **kwargs)

This is the main interface to use when interacting with FSPS from Python. Most of the Fortran API is exposed through Python hooks with various features added for user friendliness. It is recommended to only instantiate one StellarPopulation object in a given program. When initializing, you can set any of the parameters of the system using keyword arguments. Below, you’ll find a list of the options that you can include (with the comments taken directly from the FSPS docs). Unless otherwise noted, you can change these values later using the params property—which is dict-like. For example:

sp = StellarPopulation(imf_type=2, zcontinuous=1)
sp.params["imf_type"] = 1
sp.params["logzsol"] = -0.3
sp.params["sfh"] = 1
Parameters:
  • compute_vega_mags – (default: False) A switch that sets the zero points of the magnitude system: True uses Vega magnitudes versus AB magnitudes. Can only be changed during initialization.

  • vactoair_flag – (default: False) If True, output wavelengths in air (rather than vac). Can only be changed during initialization.

  • zcontinuous

    (default: 0) Flag specifying how interpolation in metallicity of the simple stellar populations (SSPs) is performed before computing composite stellar population (CSP) models:

    • 0: No interpolation, use the metallicity index specified by zmet.

    • 1: The SSPs are interpolated to the value of logzsol before the spectra and magnitudes are computed, and the value of zmet is ignored.

    • 2: The SSPs are convolved with a metallicity distribution function specified by the logzsol and pmetals parameters. The value of zmet is ignored.

    • 3: Use all available SSP metallicities when computing the composite model, for use exclusively with tabular SFHs where the metallicity evolution as function of age is given (see set_tabular_sfh()). The values of zmet and logzsol are ignored. Furthermore add_neb_emission must be set to False.

    Can only be changed during initialization.

  • add_agb_dust_model – (default: True) Switch to turn on/off the AGB circumstellar dust model presented in Villaume (2014). NB: The AGB dust emission is scaled by the parameter agb_dust.

  • add_dust_emission – (default: True) Switch to turn on/off the Draine & Li 2007 dust emission model.

  • add_igm_absorption – (default: False) Switch to include IGM absorption via Madau (1995). The zred parameter must be non-zero for this switch to have any effect. The optical depth can be scaled using the igm_factor parameter.

  • add_neb_emission – (default: False) Switch to turn on/off a nebular emission model (both continuum and line emission), based on Cloudy models from Nell Byler. Contrary to FSPS, this option is turned off by default.

  • add_neb_continuum – (default: True) Switch to turn on/off the nebular continuum component (automatically turned off if add_neb_emission is False).

  • add_stellar_remnants – (default: True) Switch to add stellar remnants in the stellar mass computation.

  • redshift_colors – (default: False) Flag specifying how to compute magnitudes. This has no effect in python-FSPS. Magnitudes are always computed at a fixed redshift specified by zred or the redshift parameter of get_mags. See get_mags for details.

  • compute_light_ages – (default: False) Flag specifying whether to compute light- and mass-weighted ages. If True then the returned spectra are actually light-weighted ages (in Gyr) at every wavelength, the returned magnitudes are filter transmission weighted averages of these, the log_lbol attribute is the bolometric luminosity weighted age, and the stellar_mass attribute gives the mass-weighted age.

  • nebemlineinspec – (default: True) Flag to include the emission line fluxes in the spectrum. Turning this off is a significant speedup in model calculation time. If not set, the line luminosities are still computed.

  • smooth_velocity – (default: True) Switch to choose smoothing in velocity space (True) or wavelength space.

  • smooth_lsf – (default: False) Switch to apply smoothing of the SSPs by a wavelength dependent line spread function. See the set_lsf() method for details. Only takes effect if smooth_velocity is True.

  • cloudy_dust – (default: False) Switch to include dust in the Cloudy tables.

  • agb_dust – (default: 1.0) Scales the circumstellar AGB dust emission.

  • tpagb_norm_type

    (default: 2) Flag specifying TP-AGB normalization scheme:

    • 0: default Padova 2007 isochrones

    • 1: Conroy & Gunn 2010 normalization

    • 2: Villaume, Conroy, Johnson 2015 normalization

  • dell – (default: 0.0) Shift in \(\log L_\mathrm{bol}\) of the TP-AGB isochrones. Note that the meaning of this parameter and the one below has changed to reflect the updated calibrations presented in Conroy & Gunn (2009). That is, these parameters now refer to a modification about the calibrations presented in that paper. Only has effect if tpagb_norm_type=1.

  • delt – (default: 0.0) Shift in \(\log T_\mathrm{eff}\) of the TP-AGB isochrones. Only has effect if tpagb_norm_type=1.

  • redgb – (default: 1.0) Modify weight given to RGB. Only available with BaSTI isochrone set.

  • agb – (default: 1.0) Modify weight given to TP-AGB. This only has effect for FSPS v3.1 or higher.

  • fcstar – (default: 1.0) Fraction of stars that the Padova isochrones identify as Carbon stars that FSPS assigns to a Carbon star spectrum. Set this to 0.0 if for example the users wishes to turn all Carbon stars into regular M-type stars.

  • sbss – (default: 0.0) Specific frequency of blue straggler stars. See Conroy et al. (2009a) for details and a plausible range.

  • fbhb – (default: 0.0) Fraction of horizontal branch stars that are blue. The blue HB stars are uniformly spread in \(\log T_\mathrm{eff}\) to \(10^4\) K. See Conroy et al. (2009a) for details and a plausible range.

  • pagb – (default: 1.0) Weight given to the post–AGB phase. A value of 0.0 turns off post-AGB stars; a value of 1.0 implies that the Vassiliadis & Wood (1994) tracks are implemented as–is.

  • frac_xrb – (default: 1.0) Scaling factor for the X-ray source spectrum to be added to the SSPs.

  • zred – (default: 0.0) Redshift. If this value is non-zero and if redshift_colors=1, the magnitudes will be computed for the spectrum placed at redshift zred.

  • zmet – (default: 1) The metallicity is specified as an integer ranging between 1 and nz. If zcontinuous > 0 then this parameter is ignored.

  • logzsol – (default: 0.0) Parameter describing the metallicity, given in units of \(\log (Z/Z_\odot)\). Only used if zcontinuous > 0.

  • pmetals – (default: 2.0) The power for the metallicty distribution function. The MDF is given by \((Z \, e^{-Z})^{\mathrm{pmetals}}\) where \(Z = z/(z_\odot \, 10^{\mathrm{logzsol}})\) and z is the metallicity in linear units (i.e., \(z_\odot = 0.019\)). Using a negative value will result in smoothing of the SSPs by a three-point triangular kernel before linear interpolation (in \(\log Z\)) to the requested metallicity. Only used if zcontinuous = 2.

  • imf_type

    (default: 2) Common variable defining the IMF type:

    • 0: Salpeter (1955)

    • 1: Chabrier (2003)

    • 2: Kroupa (2001)

    • 3: van Dokkum (2008)

    • 4: Dave (2008)

    • 5: tabulated piece-wise power law IMF, specified in imf.dat file located in the data directory

  • imf_upper_limit – (default: 120) The upper limit of the IMF, in solar masses. Note that if this is above the maximum mass in the isochrones then those stars will not contribute to the spectrum but will affect the overall IMF normalization.

  • imf_lower_limit – (default: 0.08) The lower limit of the IMF, in solar masses. Note that if this is below the minimum mass in the isochrones then those stars will not contribute to the spectrum but will affect the overall IMF normalization.

  • imf1 – (default: 1.3) Logarithmic slope of the IMF over the range \(0.08 < M < 0.5 M_\odot\). Only used if imf_type=2.

  • imf2 – (default: 2.3) Logarithmic slope of the IMF over the range \(0.5 < M < 1 M_\odot\). Only used if imf_type=2.

  • imf3 – (default: 2.3) Logarithmic slope of the IMF over the range \(1.0 < M < \mathrm{imf\_upper\_limit} M_\odot\). Only used if imf_type=2.

  • vdmc – (default: 0.08) IMF parameter defined in van Dokkum (2008). Only used if imf_type=3.

  • mdave – (default: 0.5) IMF parameter defined in Dave (2008). Only used if imf_type=4.

  • evtype – (default: -1) Compute SSPs for only the given evolutionary type. All phases used when set to -1.

  • use_wr_spectra – (default: 1) Turn on/off the WR spectral library. If off (0), will use the main default library instead

  • logt_wmb_hot – (default: 0.0) Use the Eldridge (2017) WMBasic hot star library above this value of \(\log T_\mathrm{eff}\) or 25,000K, whichever is larger.

  • add_xrb_emission – (default: 0) Turn on/off the x-ray binary population spectra from Garofali et al.

  • masscut – (default: 150.0) Truncate the IMF above this value.

  • sigma_smooth – (default: 0.0) If smooth_velocity is True, this gives the velocity dispersion in km/s. Otherwise, it gives the width of the gaussian wavelength smoothing in Angstroms. These widths are in terms of \(\sigma\), not FWHM.

  • min_wave_smooth – (default: 1e3) Minimum wavelength to consider when smoothing the spectrum.

  • max_wave_smooth – (default: 1e4) Maximum wavelength to consider when smoothing the spectrum.

  • gas_logu – (default: -2) Log of the gas ionization parameter; relevant only for the nebular emission model.

  • gas_logz – (default: 0.0) Log of the gas-phase metallicity; relevant only for the nebular emission model. In units of \(\log (Z/Z_\odot)\).

  • igm_factor – (default: 1.0) Factor used to scale the IGM optical depth.

  • sfh

    (default: 0) Defines the type of star formation history, normalized such that one solar mass of stars is formed over the full SFH. Default value is 0.

    • 0: Compute a simple stellar population (SSP).

    • 1: Tau-model. A six parameter SFH (tau model plus a constant component and a burst) with parameters tau, const, sf_start, sf_trunc, tburst, and fburst (see below).

    • 2: This option is not supported in Python-FSPS.

    • 3: Compute a tabulated SFH, which is supplied through the set_tabular_sfh method. See that method for details.

    • 4: Delayed tau-model. This is the same as option 1 except that the tau-model component takes the form \(t\,e^{−t/\tau}\).

    • 5: Delayed tau-model with a transition at a time sf_trunc to a linearly decreasing SFH with the slope specified by sf_slope. See Simha et al. 2014 for details.

  • tau – (default: 1.0) Defines e-folding time for the SFH, in Gyr. Only used if sfh=1 or sfh=4.

  • const – (default: 0.0) Defines the constant component of the SFH. This quantity is defined as the fraction of mass formed in a constant mode of SF; the range is therefore \(0 \le C \le 1\). Only used if sfh=1 or sfh=4.

  • sf_start – (default: 0.0) Start time of the SFH, in Gyr. Only used if sfh=1 or sfh=4 or sfh=5.

  • sf_trunc – (default: 0.0) Truncation time of the SFH, in Gyr. If set to 0.0, there is no trunction. Only used if sfh=1 or sfh=4 or sfh=5.

  • tage – (default: 0.0) If set to a non-zero value, the fsps.StellarPopulation.compute_csp() method will compute the spectra and magnitudes only at this age, and will therefore only output one age result. The units are Gyr. (The default is to compute and return results from \(t \approx 0\) to the maximum age in the isochrones).

  • fburst – (default: 0.0) Defines the fraction of mass formed in an instantaneous burst of star formation. Only used if sfh=1 or sfh=4.

  • tburst – (default: 11.0) Defines the age of the Universe when the burst occurs. If tburst > tage then there is no burst. Only used if sfh=1 or sfh=4.

  • sf_slope – (default: 0.0) For sfh=5, this is the slope of the SFR after time sf_trunc.

  • dust_type

    (default: 0) Common variable defining the attenuation curve for dust around ‘old’ stars:

    • 0: power law with index dust index set by dust_index.

    • 1: Milky Way extinction law (with the \(R = A_V /E(B - V)\) value given by mwr) parameterized by Cardelli et al. (1989), with variable UV bump strength (see uvb below).

    • 2: Calzetti et al. (2000) attenuation curve. Note that if this value is set then the dust attenuation is applied to all starlight equally (not split by age), and therefore the only relevant parameter is dust2, which sets the overall normalization (you must set dust1=0.0 for this to work correctly).

    • 3: allows the user to access a variety of attenuation curve models from Witt & Gordon (2000) using the parameters wgp1 and wgp2. In this case the parameters dust1 and dust2 have no effect because the WG00 models specify the full attenuation curve.

    • 4: Kriek & Conroy (2013) attenuation curve. In this model the slope of the curve, set by the parameter dust_index, is linked to the strength of the UV bump and is the offset in slope from Calzetti.

    • 5: The SMC bar extinction curve from Gordon et al. (2003)

    • 6: The Reddy et al. (2015) attenuation curve.

  • dust_tesc – (default: 7.0) Stars younger than dust_tesc are attenuated by both dust1 and dust2, while stars older are attenuated by dust2 only. Units are \(\log (\mathrm{yrs})\).

  • dust1 – (default: 0.0) Dust parameter describing the attenuation of young stellar light, i.e. where t <= dust_tesc (for details, see Conroy et al. 2009a).

  • dust2 – (default: 0.0) Dust parameter describing the attenuation of old stellar light, i.e. where t > dust_tesc (for details, see Conroy et al. 2009a).

  • dust3 – (default: 0.0) Dust parameter describing extra attenuation of old stellar light that does _not_ afect the young (t < dust_tesc ) star light. The attenuation curve will be the one specified by dust_type

  • dust_clumps – (default: -99.) Dust parameter describing the dispersion of a Gaussian PDF density distribution for the old dust. Setting this value to -99.0 sets the distribution to a uniform screen. See Conroy et al. (2009b) for details. Values other than -99 are no longer supported.

  • frac_nodust – (default: 0.0) Fraction of starlight that is not attenuated by the diffuse dust component (i.e. that is not affected by dust2).

  • frac_obrun – (default: 0.0) Fraction of the young stars (age < dust_tesc) that are not attenuated by dust1 and that do not contribute to any nebular emission, representing runaway OB stars or escaping ionizing radiation. These stars are still attenuated by dust2.

  • dust_index – (default: -0.7) Power law index of the attenuation curve. Only used when dust_type=0.

  • dust1_index – (default: -1.0) Power law index of the attenuation curve affecting stars younger than dust_tesc corresponding to dust1. Used for all dust types.

  • mwr – (default: 3.1) The ratio of total to selective absorption which characterizes the MW extinction curve: \(R = A_V /E(B - V)\). Only used when dust_type=1.

  • uvb – (default: 1.0) Parameter characterizing the strength of the 2175A extinction feature with respect to the standard Cardelli et al. determination for the MW. Only used when dust_type=1.

  • wgp1 – (default: 1) Integer specifying the optical depth in the Witt & Gordon (2000) models. Values range from 1 − 18, corresponding to optical depths of 0.25, 0.50, 0.75, 1.00, 1.50, 2.00, 2.50, 3.00, 3.50, 4.00, 4.50, 5.00, 5.50, 6.00, 7.00, 8.00, 9.00, 10.0. Note that these optical depths are defined differently from the optical depths defined by the parameters dust1 and dust2. See Witt & Gordon (2000) for details.

  • wgp2 – (default: 1) Integer specifying the type of large-scale geometry and extinction curve. Values range from 1-6, corresponding to MW+dusty, MW+shell, MW+cloudy, SMC+dusty, SMC+shell, SMC+cloudy. Dusty, shell, and cloudy specify the geometry and are described in Witt & Gordon (2000).

  • wgp3 – (default: 1) Integer specifying the local geometry for the Witt & Gordon (2000) dust models. A value of 0 corresponds to a homogeneous distribution, and a value of 1 corresponds to a clumpy distribution. See Witt & Gordon (2000) for details.

  • duste_gamma – (default: 0.01) Parameter of the Draine & Li (2007) dust emission model. Specifies the relative contribution of dust heated at a radiation field strength of \(U_\mathrm{min}\) and dust heated at \(U_\mathrm{min} < U \le U_\mathrm{max}\). Allowable range is 0.0 - 1.0.

  • duste_umin – (default: 1.0) Parameter of the Draine & Li (2007) dust emission model. Specifies the minimum radiation field strength in units of the MW value. Valid range is 0.1 - 25.0.

  • duste_qpah – (default: 3.5) Parameter of the Draine & Li (2007) dust emission model. Specifies the grain size distribution through the fraction of grain mass in PAHs. This parameter has units of % and a valid range of 0.0 - 10.0.

  • fagn – (default: 0.0) The total luminosity of the AGN, expressed as a fraction of the bolometric stellar luminosity (so it can be greater than 1). The shape of the AGN SED is from the Nenkova et al. 2008 templates.

  • agn_tau – (default: 10) Optical depth of the AGN dust torus, which affects the shape of the AGN SED. Outside the range (5, 150) the AGN SED is an extrapolation.

property dust_mass

Dust mass, in solar masses.

property duste_library

The name of the dust emission SED library being used in FSPS.

property emline_luminosity

emission line luminosities, in \(L_\odot\). shape=(ne)

property emline_wavelengths

Emission line wavelengths, in Angstroms

filter_data()

Return effective wavelengths, and vega and solar magnitudes of all filters.

Returns lambda_eff:

Effective wavelength of each filter.

Returns magvega:

The AB magnitude of Vega (used to convert between AB and Vega systems).

Returns magsun:

The AB absolute magnitude of the Sun.

property formed_mass

Integral of the SFH, in solar masses.

get_mags(zmet=None, tage=0.0, redshift=None, bands=None)

Get the magnitude of the CSP.

Parameters:
  • zmet – (default: None) The (integer) index of the metallicity to use. By default, use the current value of self.params["zmet"].

  • tage – (default: 0.0) The age of the stellar population in Gyr. By default, this will compute a grid of ages from \(t \approx 0\) to the maximum age in the isochrones.

  • redshift – (default: None) Optionally redshift the spectrum first. If not supplied, the redshift given by StellarPopulation.params["zred"] is assumed. If supplied, the value of zred is ignored (and IGM attenuation will not work properly).

  • bands – (default: None) The names of the filters that you would like to compute the magnitude for. This should correspond to the result of fsps.find_filter().

Returns mags:

The magnitude grid. If an age was was provided by the tage parameter then the result is a 1D array with NBANDS values. Otherwise, it is a 2D array with shape (NTFULL, NBANDS). If a particular set of bands was requested then this return value will be properly compressed along that axis, ordered according to the bands argument. If redshift is not 0, the units are apparent observed frame magnitude \(m\) assuming \(\Omega_m=0.3, \Omega_\Lambda=0.7\)

get_spectrum(zmet=None, tage=0.0, peraa=False)

Return spectra for the current CSP.

Parameters:
  • zmet – (default: None) The (integer) index of the metallicity to use. By default, use the current value of self.params["zmet"].

  • tage – (default: 0.0) The age of the stellar population in Gyr) for which to obtain a spectrum. By default, this will compute a grid of ages from \(t \approx 0\) to the maximum age in the isochrones.

  • peraa – (default: False) If True, return the spectrum in \(L_\odot/A\). Otherwise, return the spectrum in the FSPS standard \(L_\odot/\mathrm{Hz}\).

Returns wavelengths:

The wavelength grid in Angstroms.

Returns spectrum:

The spectrum in \(L_\odot/\mathrm{Hz}\) or \(L_\odot/A\). If an age was provided by the tage parameter then the result is a 1D array with NSPEC values. Otherwise, it is a 2D array with shape (NTFULL, NSPEC).

property isoc_library

The name of the stellar isochrone library being used in FSPS.

isochrones(outfile='pyfsps_tmp')

Write the isochrone data (age, mass, weights, phases, magnitudes, etc.) to a .cmd file, then read it into a huge numpy array. Only parameters listed in StellarPopulation.params.ssp_params affect the output of this method. This method does not work for the BPASS isochrones.

Parameters:

outfile – (default: ‘pyfsps_tmp’) The file root name of the .cmd file, which will be placed in the $SPS_HOME/OUTPUTS/ directory

Returns dat:

A huge numpy structured array containing information about every isochrone point for the current metallicity. In general the columns may be isochrone specific, but for Padova they are

  • age: log age, yrs

  • log(Z): log metallicity

  • mini: initial stellar mass in solar masses

  • mact: actual stellar mass (accounting for mass loss)

  • logl: log bolometric luminosity, solar units

  • logt: log temperature (K)

  • logg: log gravity

  • phase: (see evtype)

  • log(weight): IMF weight corresponding to a total of 1 Msol formed.

  • log(mdot): mass loss rate (Msol/yr)

property log_age

log10(age/yr).

property log_lbol

log(bolometric luminosity / \(L_\odot\)).

property resolutions

The resolution array, in km/s dispersion. Negative numbers indicate poorly defined, approximate, resolution (based on coarse opacity binning in theoretical spectra)

set_lsf(wave, sigma, wmin=None, wmax=None)

Set a wavelength dependent Gaussian line-spread function that will be applied to the SSPs. Only takes effect if smooth_lsf and smooth_velocity are True.

Parameters:
  • wave – Wavelength in angstroms, sorted ascending. If wmin or wmax are not specified they are taken from the minimum and maximum of this array. ndarray.

  • sigma – The dispersion of the Gaussian LSF at the wavelengths given by wave, in km/s. If 0, no smoothing is applied at that wavelength. ndarray of same shape as wave.

  • wmin – (optional) The minimum wavelength (in AA) for which smoothing will be applied. If not given, it is taken from the minimum of wave.

  • wmax – (optional) The maximum wavelength (in AA) for which smoothing will be applied. If not given, it is taken from the maximum of wave.

set_tabular_sfh(age, sfr, Z=None)

Set a tabular SFH for use with the sfh=3 option. See the FSPS documentation for information about tabular SFHs. This SFH will be piecewise linearly interpolated.

Parameters:
  • age – Time since the beginning of the universe in Gyr. Must be increasing. ndarray of shape (ntab,)

  • sfr – The SFR at each age, in Msun/yr. Must be an ndarray same length as age, and contain at least one non-zero value.

  • Z – (optional) The metallicity at each age, in units of absolute metallicity (e.g. Z=0.019 for solar with the Padova isochrones and MILES stellar library).

property sfr

Star formation rate (\(M_\odot/yr\)).

sfr_avg(times=None, dt=0.1)

The average SFR between time-dt and time, given the SFH parameters, for sfh=1 or sfh=4. Like SFHSTAT in FSPS. Requires scipy, as it uses gamma functions.

Parameters:
  • times – (default, None) The times (in Gyr of cosmic time) at which the average SFR over the last dt is desired. if None, uses the current value of the tage in the parameter set. Scalar or iterable.

  • dt – (default: 0.1) The time interval over which the recent SFR is averaged, in Gyr. defaults to 100 Myr (i.e. sfr8).

Returns sfr_avg:

The SFR at tage averaged over the last dt Gyr, in units of solar masses per year, such that \(1 M_\odot\) formed by tage. Same shape as times. For times < sf_start + dt the value of sfr_avg is nan, for times > tage the value is 0.

smoothspec(wave, spec, sigma, minw=None, maxw=None)

Smooth a spectrum by a gaussian with standard deviation given by sigma. Whether the smoothing is in velocity space or in wavelength space depends on the value of the value of smooth_velocity.

Parameters:
  • wave – The input wavelength grid.

  • spec – The input spectrum.

  • sigma – The standard deviation of the gaussian broadening function.

  • minw – Optionally set the minimum wavelength to consider when broadening.

  • maxw – Optionally set the maximum wavelength to consider when broadening.

Returns outspec:

The smoothed spectrum, on the same wavelength grid as the input.

property solar_metallicity

The definition of solar metallicity, as a mass fraction. E.g. Z_sol sim 0.014-0.02

property spec_library

The name of the stellar spectral library being used in FSPS.

property ssp_ages

The age grid of the SSPs, in log(years), used by FSPS.

property stellar_mass

Surviving stellar mass in solar masses (including remnants if the FSPS parameter add_stellar_remants=1).

property wavelengths

The wavelength scale for the computed spectra, in Angstroms

property zlegend

The available metallicities.