1#include <pybind11/pybind11.h>
2#include <pybind11/stl.h>
3#include <pybind11/stl_bind.h>
4#include <pybind11/numpy.h>
12namespace py = pybind11;
16 py::enum_<serif::network::NetworkFormat>(network_submodule,
"NetworkFormat")
23 py::class_<serif::network::NetIn>(network_submodule,
"NetIn")
33 ss <<
"NetIn(composition=" << netIn.composition
34 <<
", tMax=" << netIn.tMax
35 <<
", dt0=" << netIn.dt0
36 <<
", temperature=" << netIn.temperature
37 <<
", density=" << netIn.density
38 <<
", energy=" << netIn.energy <<
")";
42 py::class_<serif::network::NetOut>(network_submodule,
"NetOut")
48 ss <<
"NetOut(composition=" << netOut.composition
49 <<
", num_steps=" << netOut.num_steps
50 <<
", energy=" << netOut.energy <<
")";
54 py::class_<serif::network::Network>(network_submodule,
"Network")
55 .def(py::init<serif::network::NetworkFormat>(), py::arg(
"format"))
65 auto approx8Module = network_submodule.def_submodule(
"approx8",
"Approx8 nuclear reaction network module");
73 py::class_<vec7>(network_submodule,
"vec7")
74 .def(
"__getitem__", [](
const vec7 &v,
const size_t i) {
75 if (i >= v.size())
throw py::index_error();
78 .def(
"__setitem__", [](
vec7 &v,
const size_t i,
const double value) {
79 if (i >= v.size())
throw py::index_error();
82 .def(
"__len__", [](
const vec7 &v) {
return v.size(); })
83 .def(
"__repr__", [](
const vec7 &v) {
86 for (
size_t i = 0; i < v.size(); ++i) {
88 if (i < v.size() - 1) ss <<
", ";
94 network_submodule.def(
"sum_product", &
sum_product,
"Find the sum of the element wise product of two vec7 arrays", py::arg(
"a"), py::arg(
"b"));
95 network_submodule.def(
"get_T9_array", &
get_T9_array,
"Return an array of T9 terms for the nuclear reaction rate fit", py::arg(
"T"));
96 network_submodule.def(
"rate_fit", &
rate_fit,
"Evaluate the nuclear reaction rate given the T9 array and coefficients", py::arg(
"T9"), py::arg(
"coef"));
97 network_submodule.def(
"pp_rate", &
pp_rate,
"Calculate the rate for the reaction p + p -> d", py::arg(
"T9"));
98 network_submodule.def(
"dp_rate", &
dp_rate,
"Calculate the rate for the reaction p + d -> he3", py::arg(
"T9"));
99 network_submodule.def(
"he3he3_rate", &
he3he3_rate,
"Calculate the rate for the reaction he3 + he3 -> he4 + 2p", py::arg(
"T9"));
100 network_submodule.def(
"he3he4_rate", &
he3he4_rate,
"Calculate the rate for the reaction he3(he3,2p)he4", py::arg(
"T9"));
101 network_submodule.def(
"triple_alpha_rate", &
triple_alpha_rate,
"Calculate the rate for the reaction he4 + he4 + he4 -> c12", py::arg(
"T9"));
102 network_submodule.def(
"c12p_rate", &
c12p_rate,
"Calculate the rate for the reaction c12 + p -> n13", py::arg(
"T9"));
103 network_submodule.def(
"c12a_rate", &
c12a_rate,
"Calculate the rate for the reaction c12 + he4 -> o16", py::arg(
"T9"));
104 network_submodule.def(
"n14p_rate", &
n14p_rate,
"Calculate the rate for the reaction n14(p,g)o15 - o15 + p -> c12", py::arg(
"T9"));
105 network_submodule.def(
"n14a_rate", &
n14a_rate,
"Calculate the rate for the reaction n14(a,g)f18 assumed to go on to ne20", py::arg(
"T9"));
106 network_submodule.def(
"n15pa_rate", &
n15pa_rate,
"Calculate the rate for the reaction n15(p,a)c12 (CNO I)", py::arg(
"T9"));
107 network_submodule.def(
"n15pg_rate", &
n15pg_rate,
"Calculate the rate for the reaction n15(p,g)o16 (CNO II)", py::arg(
"T9"));
108 network_submodule.def(
"n15pg_frac", &
n15pg_frac,
"Calculate the fraction for the reaction n15(p,g)o16", py::arg(
"T9"));
109 network_submodule.def(
"o16p_rate", &
o16p_rate,
"Calculate the rate for the reaction o16(p, g)f17 then f17 -> o17(p,a)n14", py::arg(
"T9"));
110 network_submodule.def(
"o16a_rate", &
o16a_rate,
"Calculate the rate for the reaction o16(a,g)ne20", py::arg(
"T9"));
111 network_submodule.def(
"ne20a_rate", &
ne20a_rate,
"Calculate the rate for the reaction ne20(a,g)mg24", py::arg(
"T9"));
112 network_submodule.def(
"c12c12_rate", &
c12c12_rate,
"Calculate the rate for the reaction c12(c12,a)ne20", py::arg(
"T9"));
113 network_submodule.def(
"c12o16_rate", &
c12o16_rate,
"Calculate the rate for the reaction c12(o16,a)mg24", py::arg(
"T9"));
115 py::class_<Approx8Network>(network_submodule,
"Approx8Network")
117 .def(
"evaluate", &
Approx8Network::evaluate, py::arg(
"netIn"),
"Evaluate the Approx8 nuclear reaction network with the given input")
121 std::stringstream ss;
122 ss <<
"Approx8Network(stiff=" << (
network.isStiff() ?
"True" :
"False") <<
")";
Header file for the Approx8 nuclear reaction network.
Class for network evaluation.
virtual NetOut evaluate(const NetIn &netIn)
Evaluate the network based on the input parameters.
NetworkFormat setFormat(const NetworkFormat format)
NetworkFormat getFormat() const
Class for the Approx8 nuclear reaction network.
NetOut evaluate(const NetIn &netIn) override
Evaluates the nuclear network.
void setStiff(bool stiff) override
Sets whether the solver should use a stiff method.
bool isStiff() const override
Checks if the solver is using a stiff method.
double ne20a_rate(const vec7 &T9)
Calculates the rate for the reaction ne20(a,g)mg24.
double n14p_rate(const vec7 &T9)
Calculates the rate for the reaction n14(p,g)o15 - o15 + p -> c12 + he4.
double dp_rate(const vec7 &T9)
Calculates the rate for the reaction p + d -> he3.
double rate_fit(const vec7 &T9, const vec7 &coef)
double pp_rate(const vec7 &T9)
Calculates the rate for the reaction p + p -> d.
double o16a_rate(const vec7 &T9)
Calculates the rate for the reaction o16(a,g)ne20.
double o16p_rate(const vec7 &T9)
Calculates the rate for the reaction o16(p,g)f17 then f17 -> o17(p,a)n14.
double n15pg_rate(const vec7 &T9)
Calculates the rate for the reaction n15(p,g)o16 (CNO II).
double he3he4_rate(const vec7 &T9)
Calculates the rate for the reaction he3(he3,2p)he4.
double c12c12_rate(const vec7 &T9)
Calculates the rate for the reaction c12(c12,a)ne20.
double c12o16_rate(const vec7 &T9)
Calculates the rate for the reaction c12(o16,a)mg24.
double n15pa_rate(const vec7 &T9)
Calculates the rate for the reaction n15(p,a)c12 (CNO I).
double c12a_rate(const vec7 &T9)
Calculates the rate for the reaction c12 + he4 -> o16.
vec7 get_T9_array(const double &T)
double sum_product(const vec7 &a, const vec7 &b)
double triple_alpha_rate(const vec7 &T9)
Calculates the rate for the reaction he4 + he4 + he4 -> c12.
double n14a_rate(const vec7 &T9)
Calculates the rate for the reaction n14(a,g)f18 assumed to go on to ne20.
double he3he3_rate(const vec7 &T9)
Calculates the rate for the reaction he3 + he3 -> he4 + 2p.
double c12p_rate(const vec7 &T9)
Calculates the rate for the reaction c12 + p -> n13.
double n15pg_frac(const vec7 &T9)
Calculates the fraction for the reaction n15(p,g)o16.
std::array< double, 7 > vec7
Alias for a std::array of 7 doubles.
@ APPROX8
Approx8 nuclear reaction network format.
static std::unordered_map< NetworkFormat, std::string > FormatStringLookup
void register_network_bindings(pybind11::module &network_submodule)
void register_approx8_bindings(pybind11::module &network_submodule)
Input structure for the network evaluation.
double temperature
Temperature in Kelvin.
double dt0
Initial time step.
double energy
Energy in ergs.
serif::composition::Composition composition
Composition of the network.
double density
Density in g/cm^3.
Output structure for the network evaluation.
double energy
Energy in ergs after evaluation.
int num_steps
Number of steps taken in the evaluation.
serif::composition::Composition composition
Composition of the network after evaluation.