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 "helm.h"
8#include "bindings.h"
9#include "EOSio.h"
10#include "helm.h"
11
12namespace serif::eos {
13 class EOSio;
14}
15
16namespace py = pybind11;
17
18
19void register_eos_bindings(pybind11::module &eos_submodule) {
20 py::class_<serif::eos::EOSio>(eos_submodule, "EOSio")
21 .def(py::init<std::string>(), py::arg("filename"))
22 // .def("load", &EOSio::load)
23 .def("getFormat", &serif::eos::EOSio::getFormat, "Get the format of the EOS table.")
24 .def("getTable", [](serif::eos::EOSio &self) -> serif::eos::helmholtz::HELMTable* {
25 auto& table_variant = self.getTable();
26
27 // Use std::get_if to safely access the contents of the variant.
28 // This returns a pointer to the value if the variant holds that type, otherwise nullptr.
29 if (auto* ptr_to_unique_ptr = std::get_if<std::unique_ptr<serif::eos::helmholtz::HELMTable>>(&table_variant)) {
30 return (*ptr_to_unique_ptr).get();
31 }
32
33 return nullptr;
34 }, py::return_value_policy::reference_internal, // IMPORTANT: Keep this policy!
35 "Get the EOS table data.")
36 .def("__repr__", [](const serif::eos::EOSio &eos) {
37 return "<EOSio(filename='" + eos.getFilename() + "', format='" + eos.getFormatName() + "')>";
38 });
39
40 py::class_<serif::eos::EOSTable>(eos_submodule, "EOSTable");
41
42 py::class_<serif::eos::helmholtz::HELMTable>(eos_submodule, "HELMTable")
43 .def_readonly("loaded", &serif::eos::helmholtz::HELMTable::loaded)
44 .def_readonly("imax", &serif::eos::helmholtz::HELMTable::imax)
45 .def_readonly("jmax", &serif::eos::helmholtz::HELMTable::jmax)
46 .def_readonly("t", &serif::eos::helmholtz::HELMTable::t)
47 .def_readonly("d", &serif::eos::helmholtz::HELMTable::d)
48 .def("__repr__", [](const serif::eos::helmholtz::HELMTable &table) {
49 return "<HELMTable(loaded=" + std::to_string(table.loaded) + ", imax=" + std::to_string(table.imax) +
50 ", jmax=" + std::to_string(table.jmax) + ")>";
51 })
52 .def_property_readonly("f", [](serif::eos::helmholtz::HELMTable &table) -> py::array_t<double> {
53 // --- Check Preconditions ---
54 // 1. Check if dimensions are valid
55 if (table.imax <= 0 || table.jmax <= 0) {
56 // Return empty array or throw error for invalid dimensions
57 throw std::runtime_error("HELMTable dimensions (imax, jmax) are non-positive.");
58 // Alternatively: return py::array_t<double>();
59 }
60 // 2. Check if pointer 'f' and the data block 'f[0]' are non-null
61 // (Essential check assuming f could be null if not loaded/initialized)
62 if (!table.f || !table.f[0]) {
63 throw std::runtime_error("HELMTable data buffer 'f' is null or not initialized.");
64 // Alternatively: return py::array_t<double>();
65 }
66
67 // --- Get necessary info ---
68 py::ssize_t rows = static_cast<py::ssize_t>(table.imax);
69 py::ssize_t cols = static_cast<py::ssize_t>(table.jmax);
70 double* data_ptr = table.f[0]; // Pointer to the start of contiguous data block
71
72 // --- Define NumPy array shape and strides ---
73 std::vector<py::ssize_t> shape = {rows, cols};
74 std::vector<py::ssize_t> strides = {
75 static_cast<py::ssize_t>(cols * sizeof(double)), // Stride to next row
76 static_cast<py::ssize_t>( sizeof(double)) // Stride to next element in row
77 };
78
79 // --- Create and return the py::array_t ---
80 // py::cast(table) creates a py::object that acts as the 'base'.
81 // This tells NumPy not to manage the memory of 'data_ptr' and
82 // ensures the 'table' object stays alive as long as the NumPy array view exists.
83 return py::array_t<double>(
84 shape, // The dimensions of the array
85 strides, // How many bytes to step in each dimension
86 data_ptr, // Pointer to the actual data
87 py::cast(table) // Owner object (keeps C++ object alive)
88 );
89 }, py::return_value_policy::reference_internal); // Keep parent 'table' alive
90
91 py::class_<serif::eos::helmholtz::HELMEOSOutput>(eos_submodule, "EOS")
92 .def(py::init<>())
94 .def_readonly("etaele", &serif::eos::helmholtz::HELMEOSOutput::etaele)
95 .def_readonly("xnefer", &serif::eos::helmholtz::HELMEOSOutput::xnefer)
96
97 .def_readonly("ptot", &serif::eos::helmholtz::HELMEOSOutput::ptot)
98 .def_readonly("pgas", &serif::eos::helmholtz::HELMEOSOutput::pgas)
99 .def_readonly("prad", &serif::eos::helmholtz::HELMEOSOutput::prad)
100
101 .def_readonly("etot", &serif::eos::helmholtz::HELMEOSOutput::etot)
102 .def_readonly("egas", &serif::eos::helmholtz::HELMEOSOutput::egas)
103 .def_readonly("erad", &serif::eos::helmholtz::HELMEOSOutput::erad)
104
105 .def_readonly("stot", &serif::eos::helmholtz::HELMEOSOutput::stot)
106 .def_readonly("sgas", &serif::eos::helmholtz::HELMEOSOutput::sgas)
107 .def_readonly("srad", &serif::eos::helmholtz::HELMEOSOutput::srad)
108
109 .def_readonly("dpresdd", &serif::eos::helmholtz::HELMEOSOutput::dpresdd)
110 .def_readonly("dpresdt", &serif::eos::helmholtz::HELMEOSOutput::dpresdt)
111 .def_readonly("dpresda", &serif::eos::helmholtz::HELMEOSOutput::dpresda)
112 .def_readonly("dpresdz", &serif::eos::helmholtz::HELMEOSOutput::dpresdz)
113
114 .def_readonly("dentrdd", &serif::eos::helmholtz::HELMEOSOutput::dentrdd)
115 .def_readonly("dentrdt", &serif::eos::helmholtz::HELMEOSOutput::dentrdt)
116 .def_readonly("dentrda", &serif::eos::helmholtz::HELMEOSOutput::dentrda)
117 .def_readonly("dentrdz", &serif::eos::helmholtz::HELMEOSOutput::dentrdz)
118
119 .def_readonly("denerdd", &serif::eos::helmholtz::HELMEOSOutput::denerdd)
120 .def_readonly("denerdt", &serif::eos::helmholtz::HELMEOSOutput::denerdt)
121 .def_readonly("denerda", &serif::eos::helmholtz::HELMEOSOutput::denerda)
122 .def_readonly("denerdz", &serif::eos::helmholtz::HELMEOSOutput::denerdz)
123
124 .def_readonly("chiT", &serif::eos::helmholtz::HELMEOSOutput::chiT)
125 .def_readonly("chiRho", &serif::eos::helmholtz::HELMEOSOutput::chiRho)
126 .def_readonly("csound", &serif::eos::helmholtz::HELMEOSOutput::csound)
127 .def_readonly("grad_ad", &serif::eos::helmholtz::HELMEOSOutput::grad_ad)
128 .def_readonly("gamma1", &serif::eos::helmholtz::HELMEOSOutput::gamma1)
129 .def_readonly("gamma2", &serif::eos::helmholtz::HELMEOSOutput::gamma2)
130 .def_readonly("gamma3", &serif::eos::helmholtz::HELMEOSOutput::gamma3)
131 .def_readonly("cV", &serif::eos::helmholtz::HELMEOSOutput::cV)
132 .def_readonly("cP", &serif::eos::helmholtz::HELMEOSOutput::cP)
133 .def_readonly("dse", &serif::eos::helmholtz::HELMEOSOutput::dse)
134 .def_readonly("dpe", &serif::eos::helmholtz::HELMEOSOutput::dpe)
135 .def_readonly("dsp", &serif::eos::helmholtz::HELMEOSOutput::dsp)
136
137 .def("__repr__", [](const serif::eos::helmholtz::HELMEOSOutput &eos) {
138 return "<EOS (output from helmholtz eos)>";
139 });
140
141 py::class_<serif::eos::helmholtz::HELMEOSInput>(eos_submodule, "HELMEOSInput")
142 .def(py::init<>())
143 .def_readwrite("T", &serif::eos::helmholtz::HELMEOSInput::T)
144 .def_readwrite("rho", &serif::eos::helmholtz::HELMEOSInput::rho)
145 .def_readwrite("abar", &serif::eos::helmholtz::HELMEOSInput::abar)
146 .def_readwrite("zbar", &serif::eos::helmholtz::HELMEOSInput::zbar)
147 .def("__repr__", [](const serif::eos::helmholtz::HELMEOSInput &input) {
148 return "<HELMEOSInput(T=" + std::to_string(input.T) +
149 ", rho=" + std::to_string(input.rho) +
150 ", abar=" + std::to_string(input.abar) +
151 ", zbar=" + std::to_string(input.zbar) + ")>";
152 });
153
154 eos_submodule.def("get_helm_eos",
156 py::arg("q"), py::arg("table"),
157 "Calculate the Helmholtz EOS components based on input parameters and table data.");
158
159}
Handles the input/output operations for EOS tables.
Definition EOSio.h:57
std::string getFormatName() const
Gets the format name (as a string) of the EOS table.
Definition EOSio.cpp:46
EOSFormat getFormat() const
Definition EOSio.cpp:41
EOSTable & getTable()
Gets the EOS table.
Definition EOSio.cpp:51
std::string getFilename() const
Definition EOSio.h:106
void register_eos_bindings(pybind11::module &eos_submodule)
Definition bindings.cpp:19
serif::eos::helmholtz::HELMEOSOutput get_helm_EOS(serif::eos::helmholtz::HELMEOSInput &q, const serif::eos::helmholtz::HELMTable &table)
Calculate the Helmholtz EOS components.
Definition helm.cpp:243
Structure to hold the input parameters for the EOS calculation.
Definition helm.h:176
double abar
Mean atomic mass.
Definition helm.h:179
double zbar
Mean atomic number.
Definition helm.h:180
Structure to hold the Helmholtz EOS table data.
Definition helm.h:63