SYGMA Module¶
GCE SYGMA (Stellar Yields for Galaxy Modelling Applications) module
Functionality¶
This tool allows the modeling of simple stellar populations. Creating a SYGMA instance runs the simulation while extensive analysis can be done with the plot_* functions found in the chem_evol_plot module. See the DOC directory for a detailed documentation.
Made by¶
v0.1 NOV2013: C. Fryer, C. Ritter
v0.2 JAN2014: C. Ritter
- v0.3 APR2014: C. Ritter, J. F. Navarro, F. Herwig, C. Fryer, E. Starkenburg,
- M. Pignatari, S. Jones, K. Venn1, P. A. Denissenkov & the NuGrid collaboration
v0.4 FEB2015: C. Ritter, B. Cote v0.5 MAY2016: C. Ritter; Note: See from now the official releases via github.
Usage¶
Simple run example:
>>> import sygma as s
Get help with
>>> help s
Now start a calculation by providing the initial metal fraction iniZ, the final evolution time tend, and the total mass of the SSP:
>>> s1=s.sygma(iniZ=0.0001, tend=5e9, mgal=1e5)
For more information regarding the input parameter see below or try
>>> s.sygma?
For plotting utilize the plotting functions (plot_*) as described further below, for example:
>>> s1.plot_totmasses()
You can write out the information about the composition of the total ejecta of a SSP via
>>> s.write_evol_table(elements=['H','C','O'])
Yield tables are available in the NUPYCEE subdirectory yield/textunderscore tables. Add your yield tables to this directory and SYGMA will be able to read the table if you have specified the $table$ variable. Only for tables with Z=0 is the variable $pop3/textunderscore table$ used. Both tables need yields specified in the SYGMA (and OMEGA) yield input format. See the default table for the structure. It is important to provide an initial abundance file which must match the number of species provided in the yield tables. Provide the file in the iniAbu directory inside the directory yield/extunderscore tables. The input variable with which the table file can be specified is $iniabu/textunderscore table$. For the necessary structure see again the default choice of that variable.
For example with artificial yields of only H-1, you can try
>>> s2 = s.sygma(iniZ=0.0001,dt=1e8,tend=1.5e10, mgal=1e11,table='yield_tables/isotope_yield_table_h1.txt',
sn1a_table='yield_tables/sn1a_h1.txt',iniabu_table='yield_tables/iniab1.0E-04GN93_alpha_h1.ppn.txt')
-
class
sygma.
sygma
(sfr='input', imf_type='kroupa', alphaimf=2.35, imf_bdys=[0.1, 100], sn1a_rate='power_law', iniZ=0.0, dt=1000000.0, special_timesteps=30, nsmerger_bdys=[8, 100], tend=13000000000.0, mgal=10000.0, transitionmass=8, iolevel=0, ini_alpha=True, table='yield_tables/isotope_yield_table_MESA_only_fryer12_delay.txt', hardsetZ=-1, sn1a_on=True, sn1a_table='yield_tables/sn1a_t86.txt', sn1a_energy=1e+51, ns_merger_on=True, bhns_merger_on=False, f_binary=1.0, f_merger=0.0008, t_merger_max=10000000000.0, m_ej_nsm=0.025, nsm_dtd_power=[], m_ej_bhnsm=0.025, bhnsmerger_table='yield_tables/r_process_rosswog_2014.txt', nsmerger_table='yield_tables/r_process_rosswog_2014.txt', iniabu_table='', extra_source_on=False, nb_nsm_per_m=-1.0, t_nsm_coal=-1.0, extra_source_table='yield_tables/extra_source.txt', f_extra_source=1.0, pop3_table='yield_tables/popIII_heger10.txt', imf_bdys_pop3=[0.1, 100], imf_yields_range_pop3=[10, 30], starbursts=[], beta_pow=-1.0, gauss_dtd=[1000000000.0, 660000000.0], exp_dtd=2000000000.0, nb_1a_per_m=0.001, direct_norm_1a=-1, Z_trans=0.0, f_arfo=1.0, imf_yields_range=[1, 30], exclude_masses=[], netyields_on=False, wiersmamod=False, yield_interp='lin', stellar_param_on=False, t_dtd_poly_split=-1.0, stellar_param_table='yield_tables/isotope_yield_table_MESA_only_param.txt', tau_ferrini=False, dt_in=array([], dtype=float64), nsmerger_dtd_array=array([], dtype=float64), bhnsmerger_dtd_array=array([], dtype=float64), ytables_in=array([], dtype=float64), zm_lifetime_grid_nugrid_in=array([], dtype=float64), isotopes_in=array([], dtype=float64), ytables_pop3_in=array([], dtype=float64), zm_lifetime_grid_pop3_in=array([], dtype=float64), ytables_1a_in=array([], dtype=float64), mass_sampled=array([], dtype=float64), scale_cor=array([], dtype=float64), poly_fit_dtd_5th=array([], dtype=float64), poly_fit_range=array([], dtype=float64), ytables_nsmerger_in=array([], dtype=float64))[source]¶ Bases:
chem_evol.chem_evol
- sfr : string
Description of the star formation, usually an instantaneous burst.
Choices :
- ‘input’ - read and use the sfr_input file to set the percentage
- of gas that is converted into stars at each timestep.
‘schmidt’ - use an adapted Schmidt law (see Timmes95)
Default value : ‘input’
- special_timesteps : integer
Number of special timesteps. This option (already activated by default) is activated when special_timesteps > 0. It uses a logarithm timestep scheme which increases the duration of timesteps throughout the simulation.
Default value : 30
- dt : float
Duration of the first timestep [yr] if special_timesteps is activated. Duration of each timestep if special_timesteps is desactivated.
Default value : 1.0e6
- tend : float
Total duration of the simulation [yr].
Default value : 13.0e9
- imf_bdys : list
Upper and lower mass limits of the initial mass function (IMF) [Mo].
Default value : [0.1,100]
- imf_yields_range : list
Initial mass of stars that contribute to stellar ejecta [Mo].
Default value : [1,30]
- imf_type : string
Choices : ‘salpeter’, ‘chabrier’, ‘kroupa’, ‘alphaimf’ ‘alphaimf’ creates a custom IMF with a single power-law covering imf_bdys.
Default value : ‘kroupa’
- alphaimf : float
Aplha index of the custom IMF, dN/dM = Constant * M^-alphaimf
Default value : 2.35
- imf_bdys_pop3 : list
Upper and lower mass limits of the IMF of PopIII stars [Mo].
Default value : [0.1,100]
- imf_yields_range_pop3 : list
Initial mass of stars that contribute to PopIII stellar ejecta [Mo]. PopIII stars ejecta taken from Heger et al. (2010)
Default value : [10,30]
- iniZ : float
Initial metallicity of the gas in mass fraction (e.g. Solar Z = 0.02).
- Choices : 0.0, 0.0001, 0.001, 0.006, 0.01, 0.02
- (-1.0 to use non-default yield tables)
Default value : 0.0
- Z_trans : float
Variable used when interpolating stellar yields as a function of Z. Transition Z below which PopIII yields are used, and above which default yields are used.
Default value : -1 (not active)
- mgal : float
Initial mass of gas in the simulation [Mo].
Default value : 1.6e11
- sn1a_on : boolean
True or False to include or exclude the contribution of SNe Ia.
Default value : True
- sn1a_rate : string
SN Ia delay-time distribution function used to calculate the SN Ia rate.
- Choices :
‘power_law’ - custom power law, set parameter with beta_pow (similar to Maoz & Mannucci 2012)
‘gauss’ - gaussian DTD, set parameter with gauss_dtd
‘exp’ - exponential DTD, set parameter with exp_dtd
‘maoz’ - specific power law from Maoz & Mannucci (2012)
Default value : ‘power_law’
- sn1a_energy : float
Energy ejected by single SNIa event. Units in erg.
Default value : 1e51
- ns_merger_on : boolean
True or False to include or exclude the contribution of neutron star mergers.
Default value: True
- beta_pow : float
Slope of the power law for custom SN Ia rate, R = Constant * t^-beta_pow.
Default value : -1.0
- gauss_dtd : list
Contains parameter for the gaussian DTD: first the characteristic time [yrs] (gaussian center) and then the width of the distribution [yrs].
Default value : [3.3e9,6.6e8]
- exp_dtd : float
- Characteristic delay time [yrs] for the e-folding DTD.
- nb_1a_per_m : float
Number of SNe Ia per stellar mass formed in a simple stellar population.
Default value : 1.0e-03
- direct_norm_1a : float
Normalization coefficient for SNIa rate integral.
Default: deactived but replaces the usage of teh nb_1a_per_m when its value is larger than zero.
- transitionmass : float
Initial mass which marks the transition from AGB to massive stars [Mo].
Default value : 8.0
- exclude_masses : list
Contains initial masses in yield tables to be excluded from the simulation;
Default value : []
- table : string
Path pointing toward the stellar yield tables for massive and AGB stars.
Default value : ‘yield_tables/isotope_yield_table_MESA_only_fryer12_delay.txt’ (NuGrid)
- sn1a_table : string
Path pointing toward the stellar yield table for SNe Ia.
Default value : ‘yield_tables/sn1a_t86.txt’ (Tielemann et al. 1986)
- nsmerger_table : string
Path pointing toward the r-process yield tables for neutron star mergers
Default value : ‘yield_tables/r_process_rosswog_2014.txt’ (Rosswog et al. 2013)
- iniabu_table : string
Path pointing toward the table of initial abuncances in mass fraction.
Default value : ‘yield_tables/iniabu/iniab2.0E-02GN93.ppn’
- yield_interp : string
if ‘None’ : no yield interpolation, no interpolation of total ejecta
if ‘lin’ - Simple linear yield interpolation.
if ‘wiersma’ - Interpolation method which makes use of net yields as used e.g. in Wiersma+ (2009); Does not require net yields.
- if netyields_on is true than makes use of given net yields
- else calculates net yields from given X0 in yield table.
Default : ‘lin’
- netyields_on : boolean
if true assumes that yields (input from table parameter) are net yields.
Default : false.
- total_ejecta_interp : boolean
- if true then interpolates total ejecta given in yield tables
- over initial mass range.
Default : True
- stellar_param_on : boolean
if true reads in additional stellar parameter given in table stellar_param_table.
Default : true in sygma and false in omega
- stellar_param_table: string
Path pointoing toward the table hosting the evolution of stellar parameter derived from stellar evolution calculations.
Default table : ‘yield_tables/isotope_yield_table_MESA_only_param.txt’
- iolevel : int
Specifies the amount of output for testing purposes (up to 3).
Default value : 0
- poly_fit_dtd : list
Array of polynomial coefficients of a customized delay-time distribution function for SNe Ia. The polynome can be of any order. Example : [0.2, 0.3, 0.1] for rate_snIa(t) = 0.2*t**2 + 0.3*t + 0.1 Note : Must be used with the poly_fit_range parameter (see below)
Default value : np.array([]) –> Desactivated
- poly_fit_range : list –> [t_min,t_max]
Time range where the customized delay-time distribution function for SNe Ia will be applied for a simple stellar population.
Default value : np.array([]) –> Desactivated
- mass_sampled : list
Stellar masses that are sampled to eject yields in a stellar population. Warning : The use of this parameter bypasses the IMF calculation and do not ensure a correlation with the star formation rate. Each sampled mass will eject the exact amount of mass give in the stellar yields.
Default value : np.array([]) –> Desactivated
- scale_cor : 2D list
- Determine the fraction of yields ejected for any given stellar mass bin. Example : [ [1.0,8], [0.5,100] ] means that stars with initial mass between 0 and 8 Msu will eject 100% of their yields, and stars with initial mass between 8 and 100 will eject 50% of their yields. There is no limit for the number of [%,M_upper_limit] arrays used. Default value : np.array([]) –> Desactivated
Methods
-
plot_iso_ratio
(fig=18, source='all', marker='', shape='', color='', label='', fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14, return_x_y=False, xaxis='age', yaxis='C-12/C-13', markevery=500, solar_ab='yield_tables/iniabu/iniab2.0E-02GN93.ppn', solar_iso='solar_normalization/Asplund_et_al_2009_iso.txt')[source]¶ Plots ratios of isotopes.
Parameters: yaxis: string
Ratio of isotopes, e.g. C-12/C-13
Examples
>>> s.plot_iso_ratio(yaxis='C-12/C-13')
-
plot_mass
(fig=0, specie='C', source='all', norm='no', label='', shape='', marker='', color='', markevery=20, multiplot=False, return_x_y=False, fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14, linewidth=2)[source]¶ mass evolution (in Msun) of an element or isotope vs time. Note: Used in WENDI.
Parameters: specie : string
- isotope or element name, in the form ‘C’ or ‘C-12’
- ratio of element or isotope, e.g. ‘C/O’, ‘C-12/O-12’
source : string
Specifies if yields come from all sources (‘all’), including AGB+SN1a, massive stars. Or from distinctive sources: only agb stars (‘agb’), only SN1a (‘SN1a’), or only massive stars (‘massive’)
norm : string
if ‘no’, no normalization If ‘current’, normalize to current total amount of specie
label : string
figure label
marker : string
figure marker
shape : string
line style
color : string
color of line
fig : string,float
to name the plot figure
Examples
>>> s.plot_mass(specie='C')
-
plot_mass_range_contributions
(fig=7, specie='C', prodfac=False, rebin=1, time=-1, label='', shape='-', marker='o', color='r', markevery=20, extralabel=False, log=False, fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14)[source]¶ Plots yield contribution (Msun) of a certain mass range versus initial mass. Each stellar ejecta in one mass range is represented by the same yields, yields from certain stellar simulation. Be aware that a larger mass range means also a larger amount of yield for that range. Note: Used in WENDI.
Parameters: specie : string
Element or isotope name in the form ‘C’ or ‘C-12’
prodfac : boolean
if true, calculate stellar production factor for each mass range. It is the final stellar yield divided by the initial (IMF-weighted) total mass (ISM+stars) in that region.
rebin : change the bin size to uniform mass intervals of size ‘rebin’
default 0, bin size is defined by ranges of (stellar) yield inputs
log : boolean
if true, use log yaxis
label : string
figure label
marker : string
figure marker
shape : string
line style
color : string
color of line
fig : string,float
to name the plot figure
Examples
>>> s.plot_mass_range_contribution(element='C')
-
plot_mass_ratio
(fig=0, species_ratio='C/N', source='all', norm=False, label='', shape='', marker='', color='', markevery=20, multiplot=False, return_x_y=False, fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14, logy=True)[source]¶ Mass ratio of two species indicated by species_ratio over time. Choice can either be elemental ratio or isotopic ratios. Masses of species are in solar masses. Note: Similar to plot_mass but with ratios of masses.
Parameters: specie : string
ratio of element or isotope, e.g. ‘C/O’, ‘C-12/O-12’
source : string
Specifies if yields come from all sources (‘all’), including AGB+SN1a, massive stars. Or from distinctive sources: only agb stars (‘agb’), only SN1a (‘SN1a’), or only massive stars (‘massive’)
norm : boolean
If True, normalize to current total ISM mass
label : string
figure label
marker : string
figure marker
shape : string
line style
color : string
color of line
fig : string,float
to name the plot figure
logy : bool
if yes, choose yaxis in log scale
Examples
>>> s.plot_mass_ratio('C-12')
-
plot_massfrac
(fig=2, xaxis='age', yaxis='O-16', source='all', norm='no', label='', shape='', marker='', color='', markevery=20, fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14)[source]¶ Plots mass fraction of isotope or element vs either time or other isotope or element.
Parameters: xaxis : string
either ‘age’ for time or isotope name, in the form e.g. ‘C-12’
yaxis : string
isotope name, in the same form as for xaxis
source : string
Specifies if yields come from all sources (‘all’), including AGB+SN1a, massive stars. Or from distinctive sources: only agb stars (‘agb’), only SN1a (‘SN1a’), or only massive stars (‘massive’)
norm : string
if ‘no’, no normalization of mass fraction if ‘ini’, normalization ot initial abundance
label : string
figure label
marker : string
figure marker
shape : string
line style
color : string
color of line
fig : string,float
to name the plot figure
Examples
>>> s.plot_massfrac('age','C-12')
-
plot_metallicity
(source='all', label='', marker='o', color='r', shape='-', fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14)[source]¶ Plots the metal fraction defined as the sum of all isotopes except H1, H2, H3, He3, He4, Li6, Li7.
Parameters: source : string
Specifies if yields come from
all sources (‘all’), including AGB+SN1a, massive stars. Or from
distinctive sources:
only agb stars (‘agb’),
only SN1a (‘SN1a’)
only massive stars (‘massive’)
Examples
>>> s.plot_metallicity(source='all')
-
plot_net_yields
(fig=91, species='[C-12/Fe-56]', netyields_iniabu='yield_tables/iniabu/iniab_solar_Wiersma.ppn')[source]¶ Plots net yields as calculated in the code when using netyields_on=True.
To be used only with net yield input.
Parameters: species : string
Isotope ratio in spectroscopic notation.
Examples
>>> s.plot_net_yields(species='[C-12/Fe-56]')
-
plot_sn_distr
(fig=5, rate=True, rate_only='', xaxis='time', fraction=False, label1='SNIa', label2='SN2', shape1=':', shape2='--', marker1='o', marker2='s', color1='k', color2='b', markevery=20, fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14)[source]¶ Plots the SN1a distribution:
The evolution of the number of SN1a and SN2
Parameters: rate : boolean
if true, calculate rate [1/century] else calculate numbers
fraction ; boolean
if true, ignorate rate and calculate number fraction of SNIa per WD
rate_only : string
if empty string, plot both rates (default)
if ‘sn1a’, plot only SN1a rate
if ‘sn2’, plot only SN2 rate
xaxis: string
if ‘time’ : time evolution if ‘redshift’: experimental! use with caution; redshift evolution
label : string
figure label
marker : string
figure marker
shape : string
line style
color : string
color of line
fig : string,float
to name the plot figure
Examples
>>> s.plot_sn_distr()
-
plot_spectro
(fig=3, xaxis='age', yaxis='[Fe/H]', source='all', label='', shape='-', marker='o', color='k', markevery=100, show_data=False, show_sculptor=False, show_legend=True, return_x_y=False, sub_plot=False, linewidth=3, sub=1, fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14, solar_ab='')[source]¶ Plots elements in spectroscopic notation. Note: Used in WENDI.
Parameters: xaxis : string
Spectroscopic notation of elements e.g. [Fe/H] if ‘age’: time evolution in years
yaxis : string
Elements in spectroscopic notation, e.g. [C/Fe]
source : string
If yields come from all sources use ‘all’ (include AGB+SN1a, massive stars.)
If yields come from distinctive source: only agb stars use ‘agb’, only SN1a (‘SN1a’), or only massive stars (‘massive’)
label : string
figure label
marker : string
figure marker
shape : string
line style
color : string
color of line
fig : string,float
to name the plot figure
Examples
>>> plt.plot_spectro('[Fe/H]','[C/Fe]')
-
plot_star_formation_rate
(fig=6, fraction=True, source='all', marker='', shape='', color='', label='', abs_unit=True, fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14)[source]¶ Plots the star formation rate over time. Shows fraction of ISM mass which transforms into stars.
Parameters: fraction : boolean
if true: fraction of ISM which transforms into stars; else: mass of ISM which goes into stars
source : string
either ‘all’ for total star formation rate; ‘agb’ for AGB and ‘massive’ for massive stars
marker : string
marker type
shape : string
line shape type
fig: figure id
Examples
>>> s.star_formation_rate()
-
plot_stellar_param
(fig=8, quantity='Ekindot_wind', label='', marker='o', color='r', shape='-', fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14, markevery=1)[source]¶ Plots the evolution of stellar parameter as provided as input in stellar parameter table (stellar_param_table variable).
Parameters: quantity: string
Name of stellar parameter of interest. Check for available parameter via <sygma instance>.stellar_param_attrs
Examples
>>> s.plot_stellar_param(quantity='Ekin_wind')
-
plot_table_param
(fig=8, ax='', xaxis='mini', quantity='Lifetime', iniZ=0.02, masses=[], label='', marker='o', color='r', shape='-', table='yield_tables/isotope_yield_table.txt', fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14)[source]¶ Plots the yield table quantities such as lifetimes versus initial mass as given in yield input tables.
Parameters: xaxis : string
if ‘mini’: use initial mass if ‘time’: use lifetime
iniZ : float
Metallicity of interest
masses: list
List of initial masses to be plotted
table: string
Yield table
Examples
>>> s.plot_table_param(quantity='Lifetiem')
-
plot_table_remnant
(fig=8, xaxis='mini', iniZ=0.02, masses=[], label='', marker='o', color='r', shape='-', table='yield_tables/isotope_yield_table.txt', fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14)[source]¶ Plots the remnant masses versus initial mass given in yield tables.
Parameters: xaxis : string
if ‘mini’: use initial mass; if of the form [specie1/specie2] use spec. notation of
yaxis : string
iniZ : float
Metallicity to choose.
Examples
>>> s.plot_table_remnant(iniZ=0.02)
-
plot_table_yield
(fig=8, xaxis='mini', yaxis='C-12', iniZ=0.0001, netyields=False, masses=[], label='', marker='o', color='r', shape='-', table='yield_tables/isotope_yield_table.txt', fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14, solar_ab='', netyields_iniabu='')[source]¶ Plots the yields of the yield input grid versus initial mass or Z or [Fe/H]. Yield can be plotted in solar masses or in spectroscopic notation.
Parameters: xaxis : string
if ‘mini’: use initial mass; if of the form [specie1/specie2] use spec. notation
yaxis : string
specifies isotopes or elements with ‘C-12’ or ‘C’: plot yield of isotope; if chosen spectros: use form [specie3/specie4]
iniZ : float
specifies the metallicity to be plotted
masses : list
if empty plot all available masses for metallicity iniZ; else choose only masses in list masses
netyields : bool
if true assume net yields in table and add corresponding initial contribution to get total yields
table : string
table to plot data from; default sygma input table
Examples
>>> s.plot_iso_ratio(yaxis='C-12')
-
plot_table_yield_mass
(fig=8, xaxis='mini', yaxis='C-12', iniZ=0.0001, netyields=False, masses=[], label='', marker='o', color='r', shape='-', table='yield_tables/isotope_yield_table.txt', fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14)[source]¶ Plots yields for isotopes given in yield tables.
Parameters: xaxis : string
if ‘mini’ use initial mass on x axis
yaxis : string
isotope to plot.
Examples
>>> s.plot_iso_ratio(yaxis='C-12')
-
plot_totmasses
(fig=4, mass='gas', source='all', norm='no', label='', shape='', marker='', color='', markevery=20, log=True, fsize=[10, 4.5], fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14)[source]¶ Plots either gas or star mass as fraction of total mass vs time. Note: Used in WENDI.
Parameters: mass : string
either ‘gas’ for ISM gas mass or ‘stars’ for gas locked away in stars (total gas - ISM gas)
norm : string
normalization, either ‘no’ for no normalization (total gas mass in solar masses),
for normalization to the initial gas mass (mgal) with ‘ini’,
for normalization to the current total gas mass ‘current’. The latter case makes sense when comparing different sources (see below)
source : string
specifies if yields come from all sources (‘all’), including AGB+SN1a, massive stars. Or from distinctive sources: only agb stars (‘agb’), only SN1a (‘sn1a’), or only massive stars (‘massive’)
log : boolean
if true plot logarithmic y axis
label : string
figure label
marker : string
figure marker
shape : string
line style
color : string
color of line
fig : string,float
to name the plot figure
Examples
>>> s.plot_totmasses()
-
plot_yield_mtot
(fig=8, plot_imf_mass_ranges=True, fontsize=14, rspace=0.6, bspace=0.15, labelsize=15, legend_fontsize=14)[source]¶ Plots total mass ejected by stars (!). To distinguish between the total mass of yields from the table and fitted total mass
Parameters: plot_imf_mass_ranges : boolean
If true plots the initial mass ranges for which yields are consider.
Examples
>>> s.plot_yield_mtot()
-
write_evol_table
(elements=['H'], isotopes=['H-1'], table_name='gce_table.txt', path='', interact=False)[source]¶ Writes out evolution of the SSP accummulated ejecta in fraction of total mass of the SSP in the following format (each timestep in a line):
&Age &H-1 &H-2 ....
&0.000E+00 &1.000E+00 &7.600E+08 & ...
&1.000E+07 &1.000E+00 &7.600E+08 & ...
This method has a notebook and normal version. If you are using a python notebook the function will open a link to a page containing the table.
Parameters: table_name : string,optional
Name of table. If you use a notebook version, setting a name is not necessary.
elements : array
Containing the elements with the name scheme ‘H’,’C’
isotopes : array
If elements list empty, ignore elements input and use isotopes input; Containing the isotopes with the name scheme ‘H-1’, ‘C-12’
interact: bool
If true, saves file in current directory (notebook dir) and creates HTML link useful in ipython notebook environment
Examples
>>> s.write_evol_table(elements=['H','C','O'],table_name='testoutput.txt')
-
write_stellar_param_table
(table_name='gce_stellar_param_table.txt', path='evol_tables', interact=False)[source]¶ Writes out evolution of stellar parameter such as luminosity and kinetic energy. Stellar parameter quantities are available via <sygma instance>.stellar_param_attrs.
Table structure:
&Age &Quantity1 &Quantity2 ...
&0.000E+00 &0.000E+00 &0.000E+00
&0.000E+00 &0.000E+00 &0.000E+00
Parameters: table_name : string,optional
Name of table. If you use a notebook version, setting a name is not necessary.
path : string
directory where to save table.
interact: bool
If true, saves file in current directory (notebook dir) and creates HTML link useful in ipython notebook environment
Examples
>>> s.write_evol_table(table_name='testoutput.txt')