pynucastro.networks package#
The pynucastro modules that support the creation of networks. There are several main submodules here:
rate_collection
: this is
simply a collection of rates that knows about the links connecting
nuclei. This is used as the base for the different classes the write
code to output networks for integration.
python_network
: the
support routines to generate a full, integrable network in python.
base_cxx_network
:
the support routines to generate a standalone integrable network in
pure C++.
simple_cxx_network
:
the support routines to generate a simple pure C++ network for
interfacing with simulation codes.
amrexastro_cxx_network
:
the support routines to generate a C++ network that can be incorporated
into the AMReX-Astro Microphysics routines supported by astrophysical
hydrodynamics codes.
Submodules#
pynucastro.networks.amrexastro_cxx_network module#
A C++ reaction network for integration into the AMReX Astro Microphysics set of reaction networks used by astrophysical hydrodynamics codes
- class pynucastro.networks.amrexastro_cxx_network.AmrexAstroCxxNetwork(*args, **kwargs)[source]#
Bases:
BaseCxxNetwork
pynucastro.networks.base_cxx_network module#
Support for a pure C++ reaction network. These functions will write the C++ code necessary to integrate a reaction network comprised of the rates that are passed in.
- class pynucastro.networks.base_cxx_network.BaseCxxNetwork(*args, **kwargs)[source]#
Bases:
ABC
,RateCollection
Interpret the collection of rates and nuclei and produce the C++ code needed to integrate the network.
pynucastro.networks.fortran_network module#
A Fortran wrapper around the SimpleCxxNetwork
- class pynucastro.networks.fortran_network.FortranNetwork(*args, **kwargs)[source]#
Bases:
SimpleCxxNetwork
pynucastro.networks.nse_network module#
- class pynucastro.networks.nse_network.NSENetwork(rate_files=None, libraries=None, rates=None, inert_nuclei=None, symmetric_screening=False, do_screening=True)[source]#
Bases:
RateCollection
a network for solving for the NSE composition and outputting tabulated NSE quantities
- generate_table(rho_values=None, T_values=None, Ye_values=None, comp_reduction_func=None, verbose=False, outfile='nse.tbl')[source]#
- get_comp_nse(rho, T, ye, init_guess=(-3.5, -15), tol=1e-11, use_coulomb_corr=False, return_sol=False)[source]#
Returns the NSE composition given density, temperature and prescribed electron fraction using scipy.fsolve.
Parameters:#
rho: NSE state density
T: NSE state Temperature
ye: prescribed electron fraction
init_guess: optional, initial guess of chemical potential of proton and neutron, [mu_p, mu_n]
tol: optional, sets the tolerance of scipy.fsolve
use_coulomb_corr: Whether to include coulomb correction terms
return_sol: Whether to return the solution of the proton and neutron chemical potential.
pynucastro.networks.numpy_network module#
- class pynucastro.networks.numpy_network.NumpyNetwork(rate_files=None, libraries=None, rates=None, inert_nuclei=None, symmetric_screening=False, do_screening=True)[source]#
Bases:
RateCollection
A network that uses numpy arrays to evaluate rates more efficiently.
Cached arrays#
- coef_arr#
Reaclib rate coefficient array, with shape
(number_of_rates, number_of_sets, 7)
.
- coef_mask#
Boolean mask array determining how many sets to include in the final rate evaluation, with shape
(number_of_rates, number_of_sets)
.
- nuc_prod_count#
Array storing the count of each nucleus in rates producing that nucleus, with shape
(number_of_species, number_of_rates)
.
- nuc_cons_count#
Array storing the count of each nucleus in rates consuming that nucleus, with shape
(number_of_species, number_of_rates)
.
- nuc_used#
A boolean matrix of whether the nucleus is involved in the reaction or not, with shape
(number_of_species, number_of_rates)
.
- yfac#
Array storing the molar fraction component of each rate (Y of each reactant raised to the appropriate power).
Depends on composition only.
- prefac#
Array storing the prefactor for each rate, which includes both the statistical prefactor and mass density raised to the corresponding density exponent.
Depends on composition and density.
Methods#
- clear_arrays()[source]#
Clear all cached arrays stored by the
update_yfac_arr()
andupdate_prefac_arr()
member functions, freeing up memory.
- evaluate_activity_arr(T)[source]#
Sum over all of the terms contributing to dY/dt for a specific temperature, neglecting sign, assuming necessary precalculations have been carried out (calling the methods
update_yfac_arr()
andupdate_prefac_arr()
to set the composition and density).This performs a vectorized calculation, and returns an array ordered by the nuclei in the
unique_nuclei
member variable. This does not support tabular rates.See
evaluate_activity()
for the non-vectorized version. Relative performance between the two varies based on the setup. Seeclear_arrays()
for freeing memory post calculation.
- evaluate_rates_arr(T)[source]#
Evaluate the rates in the network for a specific temperature, assuming necessary precalculations have been carried out (calling the methods
update_yfac_arr()
andupdate_prefac_arr()
to set the composition and density).This performs a vectorized calculation, and returns an array ordered by the rates in the
rates
member variable. This does not support tabular rates.See
evaluate_rates()
for the non-vectorized version. Relative performance between the two varies based on the setup. Seeclear_arrays()
for freeing memory post calculation.
- evaluate_ydots_arr(T)[source]#
Evaluate net rate of change of molar abundance for each nucleus in the network for a specific temperature, assuming necessary precalculations have been carried out (calling the methods
update_yfac_arr()
andupdate_prefac_arr()
to set the composition and density).This performs a vectorized calculation, and returns an array ordered by the nuclei in the
unique_nuclei
member variable. This does not support tabular rates.See
evaluate_ydots()
for the non-vectorized version. Relative performance between the two varies based on the setup. Seeclear_arrays()
for freeing memory post calculation.
pynucastro.networks.python_network module#
Support modules to write a pure python reaction network ODE source
- class pynucastro.networks.python_network.PythonNetwork(rate_files=None, libraries=None, rates=None, inert_nuclei=None, symmetric_screening=False, do_screening=True)[source]#
Bases:
RateCollection
A pure python reaction network.
pynucastro.networks.rate_collection module#
A collection of classes and methods to deal with collections of rates that together make up a network.
- class pynucastro.networks.rate_collection.Composition(nuclei, small=1e-16)[source]#
Bases:
UserDict
a composition holds the mass fractions of the nuclei in a network – useful for evaluating the rates
- property A#
return nuclei: molar mass pairs for elements in composition
- property X#
backwards-compatible getter for self.X
- property Z#
return nuclei: charge pairs for elements in composition
- property abar#
return the mean molecular weight
- bin_as(nuclei, *, verbose=False, exclude=None)[source]#
given a list of nuclei, return a new Composition object with the current composition mass fractions binned into the new nuclei.
if a list of nuclei is provided by exclude, then only exact matches will be binned into the nuclei in that list
- plot(trace_threshold=0.1, hard_limit=None, size=(9, 5))[source]#
Make a pie chart of Composition. group trace nuclei together and explode into bar chart
parameters#
trace_threshold : the threshold to consider a component to be trace.
hard_limit : hard limit for nuclei to be labeled in the plot. Default is None which will set the limit to 5% of total trace abundance
size: tuple giving width x height of the plot in inches
- set_array(arr)[source]#
set all species from a sequence of mass fractions, in the same order as returned by get_nuclei()
- set_random(alpha=None, seed=None)[source]#
set all species using a Dirichlet distribution with parameters alpha and specified rng seed
- set_solar_like(Z=0.02)[source]#
approximate a solar abundance, setting p to 0.7, He4 to 0.3 - Z and the remainder evenly distributed with Z
- property ye#
return the electron fraction
- property zbar#
return the mean charge, Zbar
- class pynucastro.networks.rate_collection.Explorer(rc, comp, **kwargs)[source]#
Bases:
object
interactively explore a rate collection
- class pynucastro.networks.rate_collection.RateCollection(rate_files=None, libraries=None, rates=None, inert_nuclei=None, symmetric_screening=False, do_screening=True)[source]#
Bases:
object
a collection of rates that together define a network
- dedupe_partition_function_temperatures()[source]#
return a list of unique temperature arrays, along with a dictionary mapping each Nucleus to the corresponding index into that list
- evaluate_activity(rho, T, composition, screen_func=None)[source]#
sum over all of the terms contributing to ydot, neglecting sign
- evaluate_energy_generation(rho, T, composition, screen_func=None, return_enu=False)[source]#
evaluate the specific energy generation rate of the network for a specific density, temperature and composition
screen_func: (optional) a function object to call to apply screening
return_enu: (optional) return both enuc and enu – the energy loss from neutrinos from weak reactions
- evaluate_jacobian(rho, T, comp, screen_func=None)[source]#
return an array of the form J_ij = dYdot_i/dY_j for the network
- evaluate_rates(rho, T, composition, screen_func=None)[source]#
evaluate the rates for a specific density, temperature, and composition, with optional screening
Note: this returns that rate as dY/dt, where Y is the molar fraction. For a 2 body reaction, a + b, this will be of the form:
rho Y_a Y_b N_A <sigma v> / (1 + delta_{ab}
where delta is the Kronecker delta that accounts for a = b.
If you want dn/dt, where n is the number density (so you get n_a n_b <sigma v>), then you need to multiply the results here by rho N_A (where N_A is Avogadro’s number).
- evaluate_screening(rho, T, composition, screen_func)[source]#
Evaluate the screening factors for each rate, using one of the methods in
pynucastro.screening
- evaluate_ydots(rho, T, composition, screen_func=None, rate_filter=None)[source]#
evaluate net rate of change of molar abundance for each nucleus for a specific density, temperature, and composition
rho: the density (g/cm^3)
T: the temperature (K)
composition: a Composition object holding the mass fractions of the nuclei
screen_func: (optional) the name of a screening function to call to add electron screening to the rates
rate_filter_funcion: (optional) a function that takes a Rate object and returns true or false if it is to be shown as an edge.
- find_duplicate_links()[source]#
report on an rates where another rate exists that has the same reactants and products. These may not be the same Rate object (e.g., one could be tabular the other a simple decay), but they will present themselves in the network as the same link.
We return a list, where each entry is a list of all the rates that share the same link
- find_reverse(forward_rate, reverse_rates=None)[source]#
given a forward rate, locate the rate that is its reverse
- find_unimportant_rates(states, cutoff_ratio, screen_func=None)[source]#
evaluate the rates at multiple thermodynamic states, and find the rates that are always less than cutoff_ratio times the fastest rate for each state
Here, states is a list of tuple of the form (density, temperature, composition), where composition is of type Composition.
- get_nuclei_needing_partition_functions()[source]#
return a list of Nuclei that require partition functions for one or more DerivedRates in the collection
- get_rate(rid)[source]#
Return a rate matching the id provided. Here rid should be the string return by Rate.fname
- get_rate_by_nuclei(reactants, products)[source]#
given a list of reactants and products, return any matching rates
- get_rate_pairs()[source]#
return a list of RatePair objects, grouping the rates together by forward and reverse
- gridplot(comp=None, color_field='X', rho=None, T=None, **kwargs)[source]#
Plot nuclides as cells on a grid of Z vs. N, colored by color_field. If called without a composition, the function will just plot the grid with no color field.
- Parameters:
comp – Composition of the environment.
color_field – Field to color by. Must be one of ‘X’ (mass fraction), ‘Y’ (molar abundance), ‘Xdot’ (time derivative of X), ‘Ydot’ (time derivative of Y), or ‘activity’ (sum of contributions to Ydot of all rates, ignoring sign).
rho – Density to evaluate rates at. Needed for fields involving time derivatives.
T – Temperature to evaluate rates at. Needed for fields involving time derivatives.
- Keyword Arguments:
scale – One of ‘linear’, ‘log’, and ‘symlog’. Linear by default.
small – If using logarithmic scaling, zeros will be replaced with this value. 1e-30 by default.
linthresh – Linearity threshold for symlog scaling.
linscale – The number of decades to use for each half of the linear range. Stretches linear range relative to the logarithmic range.
filter_function – A callable to filter Nucleus objects with. Should return True if the nuclide should be plotted.
outfile – Output file to save the plot to. The plot will be shown if not specified.
dpi – DPI to save the image file at.
cmap – Name of the matplotlib colormap to use. Default is ‘magma’.
edgecolor – Color of grid cell edges.
area – Area of the figure without the colorbar, in square inches. 64 by default.
no_axes – Set to True to omit axis spines.
no_ticks – Set to True to omit tickmarks.
no_cbar – Set to True to omit colorbar.
cbar_label – Colorbar label.
cbar_bounds – Explicit colorbar bounds.
cbar_format – Format string or Formatter object for the colorbar ticks.
- linking_nuclei(nuclei, return_type=None, **kwargs)[source]#
Return a new RateCollection/Network object containing only rates linking the given nuclei (parameter nuclei). Nuclei can be provided as an iterable of Nucleus objects or a list of abbreviations. The return_type parameter allows the caller to specify a different constructor (e.g. superclass constructor) if the current class does not take a ‘libraries’ keyword. See method of same name in Library class for valid keyword arguments.
- make_ap_pg_approx(intermediate_nuclei=None)[source]#
combine the rates A(a,g)B and A(a,p)X(p,g)B (and the reverse) into a single effective approximate rate.
- plot(rho=None, T=None, comp=None, *, outfile=None, size=(800, 600), dpi=100, title=None, ydot_cutoff_value=None, show_small_ydot=False, node_size=1000, node_font_size=13, node_color='#A0CBE2', node_shape='o', curved_edges=False, N_range=None, Z_range=None, rotated=False, always_show_p=False, always_show_alpha=False, hide_xp=False, hide_xalpha=False, edge_labels=None, highlight_filter_function=None, nucleus_filter_function=None, rate_filter_function=None)[source]#
Make a plot of the network structure showing the links between nuclei. If a full set of thermodymamic conditions are provided (rho, T, comp), then the links are colored by rate strength.
parameters#
outfile: output name of the plot – extension determines the type
rho: density to evaluate rates with
T: temperature to evaluate rates with
comp: composition to evaluate rates with
size: tuple giving width x height of the plot in pixels
dpi: pixels per inch used by matplotlib in rendering bitmap
title: title to display on the plot
ydot_cutoff_value: rate threshold below which we do not show a line corresponding to a rate
show_small_ydot: if true, then show visible lines for rates below ydot_cutoff_value
node_size: size of a node
node_font_size: size of the font used to write the isotope in the node
node_color: color to make the nodes
node_shape: shape of the node (using matplotlib marker names)
curved_edges: do we use arcs to connect the nodes?
N_range: range of neutron number to zoom in on
Z_range: range of proton number to zoom in on
rotate: if True, we plot A - 2Z vs. Z instead of the default Z vs. N
always_show_p: include p as a node on the plot even if we don’t have p+p reactions
always_show_alpha: include He4 as a node on the plot even if we don’t have 3-alpha
hide_xalpha=False: dont connect the links to alpha for heavy nuclei reactions of the form A(alpha,X)B or A(X,alpha)B, except if alpha is the heaviest product.
hide_xp=False: dont connect the links to p for heavy nuclei reactions of the form A(p,X)B or A(X,p)B.
edge_labels: a dictionary of the form {(n1, n2): “label”} that gives labels for the edges in the network connecting nucleus n1 to n2.
highlight_filter_function: name of a custom function that takes a Rate object and returns true or false if we want to highlight the rate edge.
nucleus_filter_funcion: name of a custom function that takes a Nucleus object and returns true or false if it is to be shown as a node.
rate_filter_funcion: a function that takes a Rate object and returns true or false if it is to be shown as an edge.
- plot_jacobian(rho, T, comp, *, outfile=None, screen_func=None, rate_scaling=10000000000.0, size=(800, 800), dpi=100)[source]#
plot the Jacobian matrix of the system. Here, rate_scaling is used to set the cutoff of values that we show, relative to the peak. Any Jacobian element smaller than this will not be shown.
- plot_network_chart(rho=None, T=None, comp=None, *, outfile=None, size=(800, 800), dpi=100, force_one_column=False)[source]#
- pynucastro_dir = PosixPath('/opt/hostedtoolcache/Python/3.11.10/x64/lib/python3.11/site-packages/pynucastro')#
- remove_nuclei(nuc_list)[source]#
remove the nuclei in nuc_list from the network along with any rates that directly involve them (this doesn’t affect approximate rates that may have these nuclei as hidden intermediate links)
- remove_rates(rates)[source]#
remove the Rate objects in rates from the network. Note, if rate list is a dict, then the keys are assumed to be the rates to remove
- validate(other_library, *, forward_only=True)[source]#
perform various checks on the library, comparing to other_library, to ensure that we are not missing important rates. The idea is that self should be a reduced library where we filtered out a few rates and then we want to compare to the larger other_library to see if we missed something important.
pynucastro.networks.simple_cxx_network module#
A simple C++ reaction network for integrating into other C++ codes
- class pynucastro.networks.simple_cxx_network.SimpleCxxNetwork(*args, **kwargs)[source]#
Bases:
BaseCxxNetwork
pynucastro.networks.sympy_network_support module#
This is a module that interprets the rates, ydots, and Jacobian through sympy
- class pynucastro.networks.sympy_network_support.SympyRates[source]#
Bases:
object
- cxxify(s)[source]#
Given string s, will replace the symbols appearing as keys in self.symbol_ludict with their corresponding entries.
- jacobian_term_symbol(rate, ydot_j, y_i)[source]#
return a sympy expression containing the term in a jacobian matrix in a reaction network corresponding to this rate
Returns the derivative of the j-th YDOT wrt. the i-th Y If the derivative is zero, returns 0.
ydot_j and y_i are objects of the class ‘Nucleus’