SERiF 0.0.1a
3+1D Stellar Structure and Evolution
Loading...
Searching...
No Matches
bindings.cpp
Go to the documentation of this file.
1#include <pybind11/pybind11.h>
2#include <pybind11/stl.h> // Needed for vectors, maps, sets, strings
3#include <pybind11/stl_bind.h> // Needed for binding std::vector, std::map etc if needed directly
4#include <pybind11/numpy.h>
5
6#include <string>
7#include "bindings.h"
8
9#include "approx8.h"
10#include "network.h"
11
12namespace py = pybind11;
13
14
15void register_network_bindings(pybind11::module &network_submodule) {
16 py::enum_<serif::network::NetworkFormat>(network_submodule, "NetworkFormat")
19 .def("__str__", [](const serif::network::NetworkFormat format) {
21 });
22
23 py::class_<serif::network::NetIn>(network_submodule, "NetIn")
24 .def(py::init<>())
25 .def_readwrite("composition", &serif::network::NetIn::composition)
26 .def_readwrite("tMax", &serif::network::NetIn::tMax)
27 .def_readwrite("dt0", &serif::network::NetIn::dt0)
28 .def_readwrite("temperature", &serif::network::NetIn::temperature)
29 .def_readwrite("density", &serif::network::NetIn::density)
30 .def_readwrite("energy", &serif::network::NetIn::energy)
31 .def("__repr__", [](const serif::network::NetIn &netIn) {
32 std::stringstream ss;
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 << ")";
39 return ss.str();
40 });
41
42 py::class_<serif::network::NetOut>(network_submodule, "NetOut")
43 .def_readonly("composition", &serif::network::NetOut::composition)
44 .def_readonly("num_steps", &serif::network::NetOut::num_steps)
45 .def_readonly("energy", &serif::network::NetOut::energy)
46 .def("__repr__", [](const serif::network::NetOut &netOut) {
47 std::stringstream ss;
48 ss << "NetOut(composition=" << netOut.composition
49 << ", num_steps=" << netOut.num_steps
50 << ", energy=" << netOut.energy << ")";
51 return ss.str();
52 });
53
54 py::class_<serif::network::Network>(network_submodule, "Network")
55 .def(py::init<serif::network::NetworkFormat>(), py::arg("format"))
56 .def("evaluate", &serif::network::Network::evaluate, py::arg("netIn"))
57 .def("getFormat", &serif::network::Network::getFormat)
58 .def("setFormat", &serif::network::Network::setFormat, py::arg("format"))
59 .def("__repr__", [](const serif::network::Network &network) {
60 std::stringstream ss;
61 ss << "Network(format=" << serif::network::FormatStringLookup[network.getFormat()] << ")";
62 return ss.str();
63 });
64
65 auto approx8Module = network_submodule.def_submodule("approx8", "Approx8 nuclear reaction network module");
66 register_approx8_bindings(approx8Module);
67
68}
69
70void register_approx8_bindings(pybind11::module &network_submodule) {
71 using namespace serif::network::approx8;
72
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();
76 return v[i];
77 })
78 .def("__setitem__", [](vec7 &v, const size_t i, const double value) {
79 if (i >= v.size()) throw py::index_error();
80 v[i] = value;
81 })
82 .def("__len__", [](const vec7 &v) { return v.size(); })
83 .def("__repr__", [](const vec7 &v) {
84 std::stringstream ss;
85 ss << "vec7(";
86 for (size_t i = 0; i < v.size(); ++i) {
87 ss << v[i];
88 if (i < v.size() - 1) ss << ", ";
89 }
90 ss << ")";
91 return ss.str();
92 });
93
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"));
114
115 py::class_<Approx8Network>(network_submodule, "Approx8Network")
116 .def(py::init<>())
117 .def("evaluate", &Approx8Network::evaluate, py::arg("netIn"), "Evaluate the Approx8 nuclear reaction network with the given input")
118 .def("setStiff", &Approx8Network::setStiff, py::arg("stiff"), "Set whether to use a stiff solver or not")
119 .def("isStiff", &Approx8Network::isStiff, "Get whether the network is set to use a stiff solver or not")
120 .def("__repr__", [](const Approx8Network &network) {
121 std::stringstream ss;
122 ss << "Approx8Network(stiff=" << (network.isStiff() ? "True" : "False") << ")";
123 return ss.str();
124 });
125}
Header file for the Approx8 nuclear reaction network.
Class for network evaluation.
Definition network.h:114
virtual NetOut evaluate(const NetIn &netIn)
Evaluate the network based on the input parameters.
Definition network.cpp:54
NetworkFormat setFormat(const NetworkFormat format)
Definition network.cpp:48
NetworkFormat getFormat() const
Definition network.cpp:44
Class for the Approx8 nuclear reaction network.
Definition approx8.h:294
NetOut evaluate(const NetIn &netIn) override
Evaluates the nuclear network.
Definition approx8.cpp:449
void setStiff(bool stiff) override
Sets whether the solver should use a stiff method.
Definition approx8.cpp:509
bool isStiff() const override
Checks if the solver is using a stiff method.
Definition approx8.h:315
double ne20a_rate(const vec7 &T9)
Calculates the rate for the reaction ne20(a,g)mg24.
Definition approx8.cpp:222
double n14p_rate(const vec7 &T9)
Calculates the rate for the reaction n14(p,g)o15 - o15 + p -> c12 + he4.
Definition approx8.cpp:168
double dp_rate(const vec7 &T9)
Calculates the rate for the reaction p + d -> he3.
Definition approx8.cpp:126
double rate_fit(const vec7 &T9, const vec7 &coef)
Definition approx8.cpp:114
double pp_rate(const vec7 &T9)
Calculates the rate for the reaction p + p -> d.
Definition approx8.cpp:119
double o16a_rate(const vec7 &T9)
Calculates the rate for the reaction o16(a,g)ne20.
Definition approx8.cpp:214
double o16p_rate(const vec7 &T9)
Calculates the rate for the reaction o16(p,g)f17 then f17 -> o17(p,a)n14.
Definition approx8.cpp:208
double n15pg_rate(const vec7 &T9)
Calculates the rate for the reaction n15(p,g)o16 (CNO II).
Definition approx8.cpp:194
double he3he4_rate(const vec7 &T9)
Calculates the rate for the reaction he3(he3,2p)he4.
Definition approx8.cpp:139
double c12c12_rate(const vec7 &T9)
Calculates the rate for the reaction c12(c12,a)ne20.
Definition approx8.cpp:231
double c12o16_rate(const vec7 &T9)
Calculates the rate for the reaction c12(o16,a)mg24.
Definition approx8.cpp:237
double n15pa_rate(const vec7 &T9)
Calculates the rate for the reaction n15(p,a)c12 (CNO I).
Definition approx8.cpp:185
double c12a_rate(const vec7 &T9)
Calculates the rate for the reaction c12 + he4 -> o16.
Definition approx8.cpp:161
vec7 get_T9_array(const double &T)
Definition approx8.cpp:97
double sum_product(const vec7 &a, const vec7 &b)
Definition approx8.cpp:82
double triple_alpha_rate(const vec7 &T9)
Calculates the rate for the reaction he4 + he4 + he4 -> c12.
Definition approx8.cpp:146
double n14a_rate(const vec7 &T9)
Calculates the rate for the reaction n14(a,g)f18 assumed to go on to ne20.
Definition approx8.cpp:177
double he3he3_rate(const vec7 &T9)
Calculates the rate for the reaction he3 + he3 -> he4 + 2p.
Definition approx8.cpp:133
double c12p_rate(const vec7 &T9)
Calculates the rate for the reaction c12 + p -> n13.
Definition approx8.cpp:154
double n15pg_frac(const vec7 &T9)
Calculates the fraction for the reaction n15(p,g)o16.
Definition approx8.cpp:201
std::array< double, 7 > vec7
Alias for a std::array of 7 doubles.
Definition approx8.h:57
@ APPROX8
Approx8 nuclear reaction network format.
Definition network.h:37
static std::unordered_map< NetworkFormat, std::string > FormatStringLookup
Definition network.h:42
void register_network_bindings(pybind11::module &network_submodule)
Definition bindings.cpp:15
void register_approx8_bindings(pybind11::module &network_submodule)
Definition bindings.cpp:70
Input structure for the network evaluation.
Definition network.h:65
double temperature
Temperature in Kelvin.
Definition network.h:69
double dt0
Initial time step.
Definition network.h:68
double tMax
Maximum time.
Definition network.h:67
double energy
Energy in ergs.
Definition network.h:71
serif::composition::Composition composition
Composition of the network.
Definition network.h:66
double density
Density in g/cm^3.
Definition network.h:70
Output structure for the network evaluation.
Definition network.h:89
double energy
Energy in ergs after evaluation.
Definition network.h:92
int num_steps
Number of steps taken in the evaluation.
Definition network.h:91
serif::composition::Composition composition
Composition of the network after evaluation.
Definition network.h:90