pynucastro.networks package

Contents

pynucastro.networks package#

The pynucastro modules that support the creation of networks. There are several main submodules here:

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.

base_cxx_network: the support routines to generate a standalone integrable network in pure C++.

fortran_network: the support routines for creating Fortran wrappers around the simple C++ network.

nse_network: a network specialized in solving for nuclear statistical equilibrium.

numpy_network: a network that caches rates to allow for more efficient computation.

python_network: the support routines to generate a full, integrable network in python.

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.

simple_cxx_network: the support routines to generate a simple pure C++ network for interfacing with simulation codes.

sympy_network_support: common set of sympy routines that are used for code generation from networks.

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

Base class for a C++ network. This takes the same arguments as RateCollection and interprets the collection of rates and nuclei to produce the C++ code needed to integrate the network.

compose_jacobian()[source]#

Create the Jacobian matrix, df/dY, where f is a dY/dt and Y is a molar fraction

The Jacobian is stored as a list with each entry representing a Jacobian element. We also store whether the entry is null.

compose_ydot()[source]#

Create the expressions for dY/dt for each nucleus, where Y is the molar fraction.

This stores the result in a dict where the key is a Nucleus, and the value is a list of tuples, with the forward-reverse pairs of a rate

get_indent_amt(l, k)[source]#

Determine the amount of spaces to indent a line. This looks for a string of the form “<key>(#)”, where # is the a number that is the amount of indent.

Parameters:
  • l (str) – a line from a template file to check for an indent

  • k (str) – a keyword of the form “<key>” that appears in the line

Return type:

int

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]#

Generate a table of NSE properties. For every combination of density, temperature, and Ye, we solve for the NSE state and output composition properties to a file.

Parameters:
  • rho_values (numpy.ndarray) – values of density to use in the tabulation

  • T_values (numpy.ndarray) – values of temperature to use in the tabulation

  • Ye_values (numpy.ndarray) – values of electron fraction to use in the tabulation

  • comp_reduction_func (Callable) – a function that takes the NSE composition and return a reduced composition

  • verbose (bool) – output progress on creating the table as we go along

  • outfile (str) – filename for the NSE table output

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 (float) – NSE state density

  • T (float) – NSE state Temperature

  • ye (float) – prescribed electron fraction

  • init_guess ((tuple, list)) – initial guess of chemical potential of proton and neutron, [mu_p, mu_n]

  • tol (float) – tolerance of scipy.fsolve

  • use_coulomb_corr (bool) – whether to include coulomb correction terms

  • return_sol (bool) – whether to return the solution of the proton and neutron chemical potential.

Returns:

  • comp (Composition) – the NSE composition

  • u (ndarray) – the chemical potentials found by solving the nonlinear system.

class pynucastro.networks.nse_network.NSETableEntry(rho, T, Ye, *, comp=None, ydots=None, enu=None, comp_reduction_func=None)[source]#

Bases: object

A simple container to hold a single entry in the NSE table.

Parameters:
  • rho (float) – the density of the NSE state

  • T (float) – the temperature of the NSE state

  • Ye (float) – the electron fraction

  • comp (Composition) – the NSE composition

  • ydots (dict) – a dictionary of dY/dt keyed by Nucleus. This is meant to be the weak nuclear rates only, since those affect the NSE state.

  • enu (float) – the weak rate neutrino energy loss

  • comp_reduction_function (Callable) – a function that converts the NSE composition into a smaller set of nuclei. It takes a Composition object and returns a dictionary with the nucleus name (like “Ni56”) as the key and the corresponding mass fraction as the value. It should be ordered in the way you want the nuclei output into the NSE table file.

value()[source]#

a simple integer used for sorting. This has the format (logrho)(logT)(1-Ye)

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.

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.

clear_arrays()[source]#

Clear all cached arrays stored by the update_yfac_arr() and update_prefac_arr() member functions, freeing up memory.

property coef_arr#

Reaclib rate coefficient array, with shape (number_of_rates, number_of_sets, 7)

property coef_mask#

Boolean mask array determining how many sets to include in the final rate evaluation, with shape (number_of_rates, number_of_sets).

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() and update_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. See clear_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() and update_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. See clear_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() and update_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. See clear_arrays() for freeing memory post calculation.

property nuc_cons_count#

Array storing the count of each nucleus in rates consuming that nucleus, with shape (number_of_species, number_of_rates).

property nuc_prod_count#

Array storing the count of each nucleus in rates producing that nucleus, with shape (number_of_species, number_of_rates)

property nuc_used#

A boolean matrix of whether the nucleus is involved in the reaction or not, with shape (number_of_species, number_of_rates)

update_prefac_arr(rho, composition)[source]#

Calculate and store rate prefactors, which include both statistical prefactors and mass density raised to the corresponding density exponents. The results are stored in the prefac array.

update_yfac_arr(composition)[source]#

Calculate and store molar fraction component of each rate (Y of each reactant raised to the appropriate power). The results are stored in the yfac array.

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. 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

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, and updating the reaction rate values stored in the network.

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

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

Parameters:
  • nuclei (list, tuple) – an iterable of Nucleus objects

  • small (float) – a floor for nuclei mass fractions, used as the default value

property A#

Nucleus molar masses

Returns:

A – {Nucleus : A} pairs

Return type:

dict

property X#

backwards-compatible getter for self.X

property Z#

Nucleus charge

Returns:

Z – {Nucleus : Z} pairs

Return type:

dict

property abar#

Return the mean molecular weight

Return type:

float

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.

Parameters:
  • nuclei (list) – Input nuclei (either as string names or Nucleus objects) defining the new composition.

  • verbose (bool) – Output more information

  • exclude (bool) – List of nuclei in nuclei that only exact matches from the original composition can map into

Returns:

new_composition – The new binned composition

Return type:

Composition

get_molar()[source]#

Return a dictionary of molar fractions, Y = X/A

Returns:

molar – {Nucleus : Y}

Return type:

dict

get_nuclei()[source]#

Return a list of Nuclei objects that make up this composition.

Return type:

list

get_sum_X()[source]#

return the sum of the mass fractions

Return type:

float

normalize()[source]#

Normalize the mass fractions to sum to 1

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 (float) – the threshold to consider a component to be trace.

  • hard_limit (float) – limit below which an abundance will not be included in the trace nuclei wedget of the plot.

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

Return type:

matplotlib.figure.Figure

set_all(xval)[source]#

Set all species to the same scalar value.

Parameters:

xval (float) – mass fraction value for all species

set_array(arr)[source]#

Set the mass fractions of all species to the values in arr, get_nuclei()

Parameters:

arr (list, tuple, numpy.ndarray) – input values of mass fractions

set_equal()[source]#

Set all species to be equal

set_nuc(name, xval)[source]#

Set nuclei name to the mass fraction xval

Parameters:
set_random(alpha=None, seed=None)[source]#

Set all species using a Dirichlet distribution with parameters alpha and specified rng seed.

Parameters:
  • alpha (list, tuple, numpy.ndarray) – distribution length for the Dirichlet distribution

  • seed (float) – seed for the random number generator

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

Parameters:

Z (float) – The desired metalicity

property ye#

Return the electron fraction of the composition

Return type:

float

property zbar#

Return the mean charge, Zbar

Return type:

float

class pynucastro.networks.rate_collection.Explorer(rc, comp, **kwargs)[source]#

Bases: object

A simple class that enables interactive exploration a RateCollection, presenting density and temperature sliders to update the reaction rate values.

Parameters:
  • rc (RateCollection) – The RateCollection we will visualize.

  • comp (Composition) – A composition that will be used for evaluating the rates

  • kwargs (dict) – Additional parameters that will be passed through to the RateCollection plot() function. Note that “T” and “rho” will be ignored.

explore(logrho=(2, 6, 0.1), logT=(7, 9, 0.1))[source]#

Create the interactive visualization. This uses ipywidgets.interact to create an interactive visualization.

Parameters:
  • logrho (list, tuple) – a tuple of (starting log(rho), ending log(rho), dlogrho) that defines the range of densities to explore with an interactive slider.

  • logT (list, tuple) – a tuple of (starting log(T), ending log(T), dlogT) that defines the range of temperatures to explore with an interactive slider.

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. There are several arguments to the constructor – any combination may be supplied.

Parameters:
  • rate_files (str, list, tuple) – a string or iterable of strings of file names that define valid rates. This can include Reaclib library files storing multiple rates.

  • libraries (Library, list, tuple) – a Library or iterable of Library objects

  • rates (Rate, list, tuple) – a Rate or iterable of Rate objects

  • inert_nuclei (list, tuple) – an iterable of Nuclei that should be part of the collection but are not linked via reactions to the other Nuclei in the network.

  • symmetric_screening (bool) – symmetric screening means that we screen the reverse rates using the same factor as the forward rates, for rates computed via detailed balance.

  • do_screening (bool) – should we consider screening at all – this mainly affects whether we build the screening map

add_rates(rates)[source]#

Add new rates to the network. If the rate already exists, it will not be added. The network is then regenerated using the updated rates

Parameters:

rates (Rate, list) – a single Rate object or a list of Rate objects specifying the rates to be added to the network.

dedupe_partition_function_temperatures()[source]#

Return a list of unique temperature arrays needed by partition function tables, along with a dictionary mapping each Nucleus to the corresponding index into that list

Returns:

  • temp_arrays (list) – a list of NumPy ndarray specifying the temperature values for a particular partition function tabulation.

  • temp_indices (dict) – a dictionary that keyed on Nucleus that maps a nucleus to the index in temp_arrays containing the temperature array for its partition function data.

evaluate_activity(rho, T, composition, screen_func=None)[source]#

Compute the activity for each nucleus–the sum of abs(creation rate) + abs(destruction rate), i.e., this neglects the sign of the terms.

Parameters:
  • rho (float) – density used to evaluate rates

  • T (float) – temperature used to evaluate rates

  • composition (Composition) – composition used to evaluate rates

  • screen_func (Callable) – one of the screening functions from pynucastro.screening – if provided, then the evaluated rates will include the screening correction.

Return type:

dict

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

Parameters:
  • rho (float) – density to evaluate the rates with

  • T (float) – temperature to evaluate the rates with

  • composition (Composition) – composition to evaluate the rates with

  • screen_func (Callable) – a function from pynucastro.screening to call to compute the screening factor

  • return_enu (bool) – return both enuc and enu – the energy loss from neutrinos from weak reactions

Returns:

  • enuc (float) – the energy generation rate

  • enu (float) – the neutrino loss rate 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

Parameters:
  • rho (float) – density used to evaluate Jacobian terms

  • T (float) – temperature used to evaluate Jacobian terms

  • comp (Composition) – composition used to evaluate Jacobian terms

  • screen_func (Callable) – one of the screening functions from pynucastro.screening – if provided, then the evaluated rates will include the screening correction.

Return type:

numpy.ndarray

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).

Parameters:
  • rho (float) – density used to evaluate rates

  • T (float) – temperature used to evaluate rates

  • composition (Composition) – composition used to evaluate rates

  • screen_func (Callable) – one of the screening functions from pynucastro.screening – if provided, then the evaluated rates will include the screening correction.

Return type:

dict

evaluate_screening(rho, T, composition, screen_func)[source]#

Evaluate the screening factors for each rate.

Parameters:
  • rho (float) – density used to evaluate screening

  • T (float) – temperature used to evaluate screening

  • composition (Composition) – composition used to evaluate screening

  • screen_func (Callable) – one of the screening functions from pynucastro.screening

Return type:

dict

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

Parameters:
  • rho (float) – density used to evaluate rates

  • T (float) – temperature used to evaluate rates

  • composition (Composition) – composition used to evaluate rates

  • screen_func (Callable) – a function from pynucastro.screening used to compute the screening enhancement for the rates.

  • rate_filter (Callable) – a function that takes a Rate object and returns True or False if it is to be shown as an edge.

Return type:

dict

Check the network to see if there are multiple rates that share 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.

Return type:

list

find_reverse(forward_rate, reverse_rates=None)[source]#

Given a forward rate, locate the rate that is its reverse.

Return type:

Rate

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. This returns a dict keyed by Rate giving the ratio of the rate to the largest rate.

Parameters:
  • states (list, tuple) – A tuple of the form (density, temperature, composition), where composition is a Composition object

  • cutoff_ratio (float) – The ratio of a rate to the fastest rate, below which we consider this rate to be unimportant.

  • screen_func (Callable) – one of the screening functions from pynucastro.screening – if provided, then the evaluated rates will include the screening correction.

Return type:

dict

get_forward_rates()[source]#

Return a list of the forward (exothermic) rates in the network

Return type:

list

get_nuclei()[source]#

Get all the nuclei that are part of the network.

Return type:

list

get_nuclei_latex_string()[source]#

return a string listing the nuclei in latex format

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_name(name)[source]#

given a rate in the form ‘A(x,y)B’ return the Rate

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

get_rates()[source]#

get a list of the reaction rates in this network

get_rates_latex_table_string()[source]#
get_reverse_rates()[source]#

Return a list of the reverse (endothermic) rates). Note these may not be the same as the reverse rates identified by ReacLib.

Return type:

list

gridplot(rho=None, T=None, comp=None, color_field='X', **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:
  • rho (float) – density used to evaluate color_field

  • T (float) – temperature used to evaluate color_field

  • comp (Composition) – composition used to evaluate color_field

  • color_field (str) – 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).

  • kwargs (dict) –

    • “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.

Return type:

matplotlib.figure.Figure

linking_nuclei(nuclei, return_type=None, **kwargs)[source]#

Return a new network containing only rates linking the given nuclei.

Parameters:
  • nuclei (list, tuple) – An iterable of Nucleus objects or string names of nuclei.

  • return_type (Callable) – A different constructor (e.g., a superclass constructor) to use if the current class does not take a libraries keyword.

  • kwargs (dict) – Additional arguments to pass onto the library linking_nuclei method. See pynucastro.rates.library.Library.linking_nuclei

Return type:

RateCollection

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. The new approximate rates will be added to the network and the original rates will be removed (although they are still carried by the ApproximateRate object.

Parameters:

intermediate_nuclei (list, tuple) – an iterable of Nucleus objects or string names representing the intermediate nucleus we wish to approximate out.

make_nn_g_approx(intermediate_nuclei=None)[source]#

Combine the rates A(n,g)X(n,g)B into a single effective rate. The new approximate rates will be added to the network and the original rates will be removed (although they are still carried by the ApproximateRate object.

Parameters:

intermediate_nuclei (list, tuple) – an iterable of Nucleus objects or string names representing the intermediate nucleus we wish to approximate out.

make_nse_protons(A)[source]#

for rates involving nuclei with mass number >= A, swap any protons for NSE protons. This will decouple these rates from the proton captures at lower mass number, simplifying the linear algebra.

Parameters:

A (int) – mass number above which to swap regular protons for NSE protons.

network_overview()[source]#

Return a verbose network overview showing for each nucleus which rates consume it and which produce it.

Return type:

str

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=12, node_color='#444444', 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, legend_coord=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:
  • rho (float) – density to evaluate rates with

  • T (float) – temperature to evaluate rates with

  • comp (Composition) – composition to evaluate rates with

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

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

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

  • title (str) – title to display on the plot

  • ydot_cutoff_value (float) – rate threshold below which we do not show a line corresponding to a rate

  • show_small_ydot (bool) – show visible dashed lines for rates below ydot_cutoff_value

  • node_size (float) – size of a node (in networkx units)

  • node_font_size (float) – size of the font used to write the isotope in the node

  • node_color (str) – color to make the nodes

  • node_shape (str) – shape of the node (using matplotlib marker names)

  • curved_edges (bool) – do we use arcs to connect the nodes?

  • N_range ((tuple, list)) – range of neutron number to zoom in on

  • Z_range ((tuple, list)) – range of proton number to zoom in on

  • rotate (bool) – plot A - 2Z vs. Z instead of the default Z vs. N

  • always_show_p (bool) – include p as a node on the plot even if we don’t have p+p reactions

  • always_show_alpha (bool) – include He4 as a node on the plot even if we don’t have 3-alpha

  • hide_xalpha (bool) – 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 (bool) – dont connect the links to p for heavy nuclei reactions of the form A(p,X)B or A(X,p)B.

  • edge_labels (dict) – 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 (Callable) – a function that takes a Rate object and returns True or False if we want to highlight the rate edge.

  • nucleus_filter_function (Callable) – a function that takes a Nucleus object and returns True or False if it is to be shown as a node.

  • rate_filter_function (Callable) – a function that takes a Rate object and returns True or False if it is to be shown as an edge.

Return type:

matplotlib.figure.Figure

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.

Parameters:
  • rho (float) – density used to evaluate terms

  • T (float) – temperature used to evaluate terms

  • comp (Composition) – composition used to evaluate terms

  • outfile (str) – output file for plot (extension is used to specify file type)

  • rate_scaling (float) – the cutoff of values that we show, relative to the peak. Any Jacobian element smaller than this will not be shown.

  • size ((tuple, list)) – size in pixels for the output plot

  • dpi (float) – dots per inch for the output plot

Return type:

matplotlib.figure.Figure

plot_network_chart(rho=None, T=None, comp=None, *, outfile=None, size=(800, 800), dpi=100, force_one_column=False)[source]#

Plot a heatmap showing which rates are affected by which nuclei.

Parameters:
  • rho (float) – density used to evaluate rates

  • T (float) – temperature used to evaluate rates

  • comp (Composition) – composition used to evaluate rates

  • outfile

Return type:

matplotlib.figure.Figure

pynucastro_dir = PosixPath('/opt/hostedtoolcache/Python/3.11.11/x64/lib/python3.11/site-packages/pynucastro')#
rate_pair_overview()[source]#

Return a verbose network overview in terms of forward-reverse pairs

Return type:

str

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.

write_network(*args, **kwargs)[source]#

Before writing the network, check to make sure the rates are distinguishable by name.

exception pynucastro.networks.rate_collection.RateDuplicationError[source]#

Bases: Exception

An error of multiple rates linking the same nuclei occurred

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’

specific_rate_symbol(rate)[source]#

return a sympy expression containing the term in a dY/dt equation in a reaction network corresponding to this rate.

Also enter the symbol and substitution in the lookup table.

ydot_term_symbol(rate, y_i)[source]#

return a sympy expression containing this rate’s contribution to the ydot term for nuclide y_i.