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, do_screening=True, verbose=False)[source]#

Bases: RateCollection

A pure python reaction network. This can create a python module as a file that contains everything needed to evaluate the reaction rates and construct the righthand side and Jacobian functions.

full_jacobian_element_string(ydot_i_nucleus, y_j_nucleus, indent='')[source]#

Construct a string containing the python code for a single element of the Jacobian, dYdot(ydot_i_nucleus)/dY(y_j_nucleus)

Parameters:
  • ydot_i_nucleus (Nucleus) – The nucleus representing the dY/dt term we are differentiating. This is the row of the Jacobian.

  • ydot_j_nucleus (Nucleus) – The nucleus we are differentiating with respect to. This is the column of the Jacobian.

  • indent (str) – A string that will be prepended to each line of the output, typically consisting of just spaces representing the amount of indent desired.

Return type:

str

full_ydot_string(nucleus, indent='')[source]#

Construct a string containing the python code for dY(nucleus)/dt by considering every reaction that involves nucleus, adding terms that result in its creation and subtracting terms representing its destruction.

Parameters:
  • nucleus (Nucleus) – The nucleus we are constructing the time derivative of.

  • indent (str) – A string that will be prepended to each line of the output, typically consisting of just spaces representing the amount of indent desired.

Return type:

str

integrate_network(tmax, rho, T, Y0=None, screen_method=None, initial_comp='uniform', rtol=1e-08, atol=1e-08)[source]#

Integrate the network to tmax given (rho, T, Y0) using SciPy’s solve_ivp() with BDF method.

Parameters:
  • tmax (float) – final integration time.

  • rho (float) – density used to integrate the network

  • T (float) – temperature used to integrate the network

  • Y0 (numpy.ndarray) – initial molar abundance of the nuclei. If not provided, the initial composition is initialized according to initial_comp

  • screen_method (str) – name of the screening function used to evaluate rates when integrating the network. Valid choices are: screen5, chugunov_2007, chugunov_2009, potekhin_1998, and debye_huckel. If None, no screening is applied.

  • initial_comp (str) – different modes to use to set up the initial composition if Y0 is None. Valid choices are: uniform, random, and solar.

  • rtol (float) – relative tolerance for SciPy’s solve_ivp()

  • atol (float) – absolute tolerance for SciPy’s solve_ivp()

Returns:

sol – Solution returned by scipy.integrate.solve_ivp().

Return type:

object

plot_evolution(sol, tmin=None, tmax=None, ymin=None, ymax=None, size=(800, 600), dpi=100, X_cutoff_value=None, label_size=14, legend_size=10, three_level_style=False, outfile=None)[source]#

Plot the time evolution of nuclei mass fractions using the solution returned by SciPy’s solve_ivp().

Parameters:
  • sol (object) – Solution object returned by scipy.integrate.solve_ivp(). The array sol.y is assumed to contain the molar abundances, Y_i, ordered consistently with unique_nuclei.

  • tmin (float) – Minimum time shown on the x-axis. If None, the first value of sol.t is used.

  • tmax (float) – Maximum time shown on the x-axis. If None, the last value of sol.t is used.

  • ymin (float) – Minimum mass fraction shown on the y-axis. If None, use the Matplotlib autoscaled value.

  • ymax (float) – Maximum mass fraction shown on the y-axis. If None, use the Matplotlib autoscaled value. The autoscaled value is capped at 1.2

  • dpi (int) – dots per inch used with size to set output image size

  • size ((tuple, list)) – (width, height) of the plot in pixels

  • X_cutoff_value (float) – Minimum peak mass fraction required for a nucleus to be plotted.

  • label_size (int) – Font size for axis labels.

  • legend_size (int) – Font size for the legend.

  • three_level_style (bool) – If True, use three-level linestyle and linewidth based on the peak mass fraction to help distinguish different curves. If False, all curves use the same line style and linewidth.

  • outfile (str) – output name of the plot (extension determines the type)

Return type:

matplotlib.figure.Figure

rates_string(indent='')[source]#

Create the python code that calls the evaluation function for each rate. Typically this is of the form name(rate_eval, ...), where rate_eval is a container class in the output network that stores the rate values. This also calls screening_string() after the main rates are evaluated but before the approximate rates are constructed.

Parameters:

indent (str) – A string that will be prepended to each line of the output, typically consisting of just spaces representing the amount of indent desired.

Return type:

str

screening_string(indent='')[source]#

Create a string containing the python code that sets up the screening (PlasmaState) and calls the screening function on every set of reactants in our network. This computes log(screening) terms and store them in local variables.

Parameters:

indent (str) – A string that will be prepended to each line of the output, typically consisting of just spaces representing the amount of indent desired.

Return type:

str