1#include <pybind11/pybind11.h>
2#include <pybind11/stl.h>
3#include <pybind11/stl_bind.h>
14namespace py = pybind11;
20 template <IsDynamicEngine T, IsDynamicEngine BaseT>
21 void registerDynamicEngineDefs(py::class_<T, BaseT> pyClass) {
22 pyClass.def(
"calculateRHSAndEnergy", &T::calculateRHSAndEnergy,
26 "Calculate the right-hand side (dY/dt) and energy generation rate."
28 .def(
"generateJacobianMatrix", py::overload_cast<
const std::vector<double>&,
double,
double>(&T::generateJacobianMatrix, py::const_),
32 "Generate the Jacobian matrix for the current state."
34 .def(
"generateStoichiometryMatrix", &T::generateStoichiometryMatrix)
35 .def(
"calculateMolarReactionFlow",
36 static_cast<double (T::*)(
const gridfire::reaction::Reaction&,
const std::vector<double>&,
const double,
const double) const
>(&T::calculateMolarReactionFlow),
41 "Calculate the molar reaction flow for a given reaction."
43 .def(
"getNetworkSpecies", &T::getNetworkSpecies,
44 "Get the list of species in the network."
46 .def(
"getNetworkReactions", &T::getNetworkReactions,
47 "Get the set of logical reactions in the network."
49 .def (
"setNetworkReactions", &T::setNetworkReactions,
51 "Set the network reactions to a new set of reactions."
53 .def(
"getJacobianMatrixEntry", &T::getJacobianMatrixEntry,
56 "Get an entry from the previously generated Jacobian matrix."
58 .def(
"getStoichiometryMatrixEntry", &T::getStoichiometryMatrixEntry,
59 py::arg(
"speciesIndex"),
60 py::arg(
"reactionIndex"),
61 "Get an entry from the stoichiometry matrix."
63 .def(
"getSpeciesTimescales", &T::getSpeciesTimescales,
67 "Get the timescales for each species in the network."
69 .def(
"getSpeciesDestructionTimescales", &T::getSpeciesDestructionTimescales,
73 "Get the destruction timescales for each species in the network."
75 .def(
"update", &T::update,
77 "Update the engine state based on the provided NetIn object."
79 .def(
"setScreeningModel", &T::setScreeningModel,
80 py::arg(
"screeningModel"),
81 "Set the screening model for the engine."
83 .def(
"getScreeningModel", &T::getScreeningModel,
84 "Get the current screening model of the engine."
86 .def(
"getSpeciesIndex", &T::getSpeciesIndex,
88 "Get the index of a species in the network."
90 .def(
"mapNetInToMolarAbundanceVector", &T::mapNetInToMolarAbundanceVector,
92 "Map a NetIn object to a vector of molar abundances."
94 .def(
"primeEngine", &T::primeEngine,
96 "Prime the engine with a NetIn object to prepare for calculations."
98 .def(
"getDepth", &T::getDepth,
99 "Get the current build depth of the engine."
101 .def(
"rebuild", &T::rebuild,
102 py::arg(
"composition"),
104 "Rebuild the engine with a new composition and build depth."
106 .def(
"isStale", &T::isStale,
108 "Check if the engine is stale based on the provided NetIn object."
119 py::arg(
"composition"),
121 py::arg(
"reverse") =
false,
122 "Build a nuclear network from a composition using ReacLib data."
125 py::enum_<gridfire::PrimingReportStatus>(m,
"PrimingReportStatus")
135 std::stringstream ss;
139 "String representation of the PrimingReport."
142 py::class_<gridfire::PrimingReport>(m,
"PrimingReport")
148 std::stringstream ss;
157 py::class_<gridfire::StepDerivatives<double>>(m,
"StepDerivatives")
161 py::class_<gridfire::SparsityPattern>(m,
"SparsityPattern");
169 py::class_<gridfire::Engine, PyEngine>(m,
"Engine");
173 const auto a = py::class_<gridfire::DynamicEngine, PyDynamicEngine>(m,
"DynamicEngine");
177 py::enum_<gridfire::NetworkBuildDepth>(m,
"NetworkBuildDepth")
186 py::class_<gridfire::BuildDepthType>(m,
"BuildDepthType");
188 auto py_dynamic_engine_bindings = py::class_<gridfire::GraphEngine, gridfire::DynamicEngine>(m,
"GraphEngine");
191 py_dynamic_engine_bindings.def(py::init<const fourdst::composition::Composition &, const gridfire::BuildDepthType>(),
192 py::arg(
"composition"),
194 "Initialize GraphEngine with a composition and build depth."
196 py_dynamic_engine_bindings.def(py::init<const fourdst::composition::Composition &, const gridfire::partition::PartitionFunction &, const gridfire::BuildDepthType>(),
197 py::arg(
"composition"),
198 py::arg(
"partitionFunction"),
200 "Initialize GraphEngine with a composition, partition function and build depth."
202 py_dynamic_engine_bindings.def(py::init<const gridfire::reaction::LogicalReactionSet &>(),
203 py::arg(
"reactions"),
204 "Initialize GraphEngine with a set of reactions."
207 py::arg(
"Y_dynamic"),
210 py::arg(
"sparsityPattern"),
211 "Generate the Jacobian matrix for the current state with a specified sparsity pattern."
215 "Get the net stoichiometry for a given reaction."
219 "Check if a given species is involved in the network."
223 "Export the network to a DOT file for visualization."
227 "Export the network to a CSV file for analysis."
230 py::arg(
"precompute"),
231 "Enable or disable precomputation for the engine."
234 "Check if precomputation is enabled for the engine."
237 "Get the partition function used by the engine."
242 "Calculate the reverse rate for a given reaction at a specific temperature."
247 py::arg(
"forwardRate"),
248 py::arg(
"expFactor"),
249 "Calculate the reverse rate for a two-body reaction at a specific temperature."
254 py::arg(
"reverseRate"),
255 "Calculate the derivative of the reverse rate for a two-body reaction at a specific temperature."
258 "Check if the engine is using reverse reactions."
261 py::arg(
"useReverse"),
262 "Enable or disable the use of reverse reactions in the engine."
267 registerDynamicEngineDefs<gridfire::GraphEngine, gridfire::DynamicEngine>(py_dynamic_engine_bindings);
271 auto py_defined_engine_view_bindings = py::class_<gridfire::DefinedEngineView, gridfire::DynamicEngine>(m,
"DefinedEngineView");
275 py::arg(
"baseEngine"),
276 "Construct a defined engine view with a list of tracked reactions and a base engine.");
278 "Get the base engine associated with this defined engine view.");
280 registerDynamicEngineDefs<gridfire::DefinedEngineView, gridfire::DynamicEngine>(py_defined_engine_view_bindings);
282 auto py_file_defined_engine_view_bindings = py::class_<gridfire::FileDefinedEngineView, gridfire::DefinedEngineView>(m,
"FileDefinedEngineView");
283 py_file_defined_engine_view_bindings.def(py::init<gridfire::DynamicEngine&, const std::string&, const gridfire::io::NetworkFileParser&>(),
284 py::arg(
"baseEngine"),
287 "Construct a defined engine view from a file and a base engine."
290 "Get the network file associated with this defined engine view."
293 "Get the parser used for this defined engine view."
296 "Get the base engine associated with this file defined engine view.");
298 registerDynamicEngineDefs<gridfire::FileDefinedEngineView, gridfire::DefinedEngineView>(py_file_defined_engine_view_bindings);
300 auto py_priming_engine_view_bindings = py::class_<gridfire::NetworkPrimingEngineView, gridfire::DefinedEngineView>(m,
"NetworkPrimingEngineView");
301 py_priming_engine_view_bindings.def(py::init<const std::string&, gridfire::DynamicEngine&>(),
302 py::arg(
"primingSymbol"),
303 py::arg(
"baseEngine"),
304 "Construct a priming engine view with a priming symbol and a base engine.");
305 py_priming_engine_view_bindings.def(py::init<const fourdst::atomic::Species&, gridfire::DynamicEngine&>(),
306 py::arg(
"primingSpecies"),
307 py::arg(
"baseEngine"),
308 "Construct a priming engine view with a priming species and a base engine.");
310 "Get the base engine associated with this priming engine view.");
312 registerDynamicEngineDefs<gridfire::NetworkPrimingEngineView, gridfire::DefinedEngineView>(py_priming_engine_view_bindings);
314 auto py_adaptive_engine_view_bindings = py::class_<gridfire::AdaptiveEngineView, gridfire::DynamicEngine>(m,
"AdaptiveEngineView");
315 py_adaptive_engine_view_bindings.def(py::init<gridfire::DynamicEngine&>(),
316 py::arg(
"baseEngine"),
317 "Construct an adaptive engine view with a base engine.");
319 "Get the base engine associated with this adaptive engine view.");
321 registerDynamicEngineDefs<gridfire::AdaptiveEngineView, gridfire::DynamicEngine>(py_adaptive_engine_view_bindings);
323 auto py_qse_cache_config = py::class_<gridfire::QSECacheConfig>(m,
"QSECacheConfig");
324 auto py_qse_cache_key = py::class_<gridfire::QSECacheKey>(m,
"QSECacheKey");
326 py_qse_cache_key.def(py::init<
double,
double,
const std::vector<double>&>(),
333 "Get the pre-computed hash value of the key");
338 "bin a value based on a tolerance");
339 py_qse_cache_key.def(
"__eq__", &gridfire::QSECacheKey::operator==,
341 "Check if two QSECacheKeys are equal");
343 auto py_multiscale_engine_view_bindings = py::class_<gridfire::MultiscalePartitioningEngineView, gridfire::DynamicEngine>(m,
"MultiscalePartitioningEngineView");
344 py_multiscale_engine_view_bindings.def(py::init<gridfire::GraphEngine&>(),
345 py::arg(
"baseEngine"),
346 "Construct a multiscale partitioning engine view with a base engine.");
348 "Get the base engine associated with this multiscale partitioning engine view.");
350 py::arg(
"timescale_pools"),
354 "Analyze the connectivity of timescale pools in the network.");
359 "Partition the network based on species timescales and connectivity.");
362 "Partition the network based on a NetIn object.");
368 "Export the network to a DOT file for visualization.");
370 "Get the list of fast species in the network.");
372 "Get the list of dynamic species in the network.");
377 "Equilibrate the network based on species abundances and conditions.");
380 "Equilibrate the network based on a NetIn object.");
382 registerDynamicEngineDefs<gridfire::MultiscalePartitioningEngineView, gridfire::DynamicEngine>(py_multiscale_engine_view_bindings);
const DynamicEngine & getBaseEngine() const override
Gets the base engine.
const DynamicEngine & getBaseEngine() const override
Access the underlying engine instance.
Abstract class for engines supporting Jacobian and stoichiometry operations.
std::string getNetworkFile() const
const io::NetworkFileParser & getParser() const
bool isPrecomputationEnabled() const
double calculateReverseRateTwoBody(const reaction::Reaction &reaction, const double T9, const double forwardRate, const double expFactor) const
double calculateReverseRate(const reaction::Reaction &reaction, double T9) const
void setUseReverseReactions(bool useReverse)
void setPrecomputation(bool precompute)
void exportToCSV(const std::string &filename) const
Exports the network to a CSV file for analysis.
static std::unordered_map< fourdst::atomic::Species, int > getNetReactionStoichiometry(const reaction::Reaction &reaction)
Gets the net stoichiometry for a given reaction.
bool involvesSpecies(const fourdst::atomic::Species &species) const
Checks if a given species is involved in the network.
void generateJacobianMatrix(const std::vector< double > &Y_dynamic, const double T9, const double rho) const override
Generates the Jacobian matrix for the current state.
void exportToDot(const std::string &filename) const
Exports the network to a DOT file for visualization.
const partition::PartitionFunction & getPartitionFunction() const
bool isUsingReverseReactions() const
double calculateReverseRateTwoBodyDerivative(const reaction::Reaction &reaction, const double T9, const double reverseRate) const
const std::vector< fourdst::atomic::Species > & getDynamicSpecies() const
Gets the dynamic species in the network.
const DynamicEngine & getBaseEngine() const override
Gets the base engine.
std::vector< fourdst::atomic::Species > getFastSpecies() const
Gets the fast species in the network.
fourdst::composition::Composition equilibrateNetwork(const std::vector< double > &Y, double T9, double rho)
Equilibrates the network by partitioning and solving for QSE abundances.
void partitionNetwork(const std::vector< double > &Y, double T9, double rho)
Partitions the network into dynamic and algebraic (QSE) groups based on timescales.
void exportToDot(const std::string &filename, const std::vector< double > &Y, const double T9, const double rho) const
Exports the network to a DOT file for visualization.
std::vector< std::vector< size_t > > analyzeTimescalePoolConnectivity(const std::vector< std::vector< size_t > > ×cale_pools, const std::vector< double > &Y, double T9, double rho) const
Analyzes the connectivity of timescale pools.
Represents a single nuclear reaction from a specific data source.
void abs_stype_register_engine_bindings(pybind11::module &m)
void abs_stype_register_dynamic_engine_bindings(pybind11::module &m)
void con_stype_register_graph_engine_bindings(pybind11::module &m)
void register_engine_view_bindings(pybind11::module &m)
void register_base_engine_bindings(pybind11::module &m)
void register_engine_bindings(py::module &m)
Core header for the GridFire reaction network engine module.
std::map< PrimingReportStatus, std::string > PrimingReportStatusStrings
Mapping from PrimingReportStatus codes to human-readable strings.
std::vector< std::pair< size_t, size_t > > SparsityPattern
PrimingReportStatus
Enumerates outcome codes for a network priming operation.
@ FAILED_TO_FIND_PRIMING_REACTIONS
@ FAILED_TO_FIND_CREATION_CHANNEL
@ BASE_NETWORK_TOO_SHALLOW
@ FAILED_TO_FINALIZE_COMPOSITION
reaction::LogicalReactionSet build_reaclib_nuclear_network(const fourdst::composition::Composition &composition, BuildDepthType maxLayers=NetworkBuildDepth::Full, bool reverse=false)
Builds a nuclear reaction network from the Reaclib library based on an initial composition.
Captures the result of a network priming operation.
fourdst::composition::Composition primedComposition
std::vector< std::pair< fourdst::atomic::Species, double > > massFractionChanges
PrimingReportStatus status
size_t hash() const
Computes the hash value for this key.
static long bin(double value, double tol)
Converts a value to a discrete bin based on a tolerance.
T nuclearEnergyGenerationRate
Specific energy generation rate (e.g., erg/g/s).
std::vector< T > dydt
Derivatives of abundances (dY/dt for each species).