GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
engine_graph.cpp
Go to the documentation of this file.
3#include "gridfire/network.h"
4
5#include "fourdst/composition/species.h"
6#include "fourdst/composition/atomicSpecies.h"
7
8#include "quill/LogMacros.h"
9
10#include <cstdint>
11#include <iostream>
12#include <set>
13#include <stdexcept>
14#include <string>
15#include <string_view>
16#include <unordered_map>
17#include <utility>
18#include <vector>
19#include <fstream>
20
21#include <boost/numeric/odeint.hpp>
22
23
24namespace gridfire {
26 const fourdst::composition::Composition &composition
27 ):
28 m_reactions(build_reaclib_nuclear_network(composition, false)) {
30 }
31
36
38 const std::vector<double> &Y,
39 const double T9,
40 const double rho
41 ) const {
42 return calculateAllDerivatives<double>(Y, T9, rho);
43 }
44
45
54
55 // --- Network Graph Construction Methods ---
57 m_networkSpecies.clear();
58 m_networkSpeciesMap.clear();
59
60 std::set<std::string_view> uniqueSpeciesNames;
61
62 for (const auto& reaction: m_reactions) {
63 for (const auto& reactant: reaction.reactants()) {
64 uniqueSpeciesNames.insert(reactant.name());
65 }
66 for (const auto& product: reaction.products()) {
67 uniqueSpeciesNames.insert(product.name());
68 }
69 }
70
71 for (const auto& name: uniqueSpeciesNames) {
72 auto it = fourdst::atomic::species.find(std::string(name));
73 if (it != fourdst::atomic::species.end()) {
74 m_networkSpecies.push_back(it->second);
75 m_networkSpeciesMap.insert({name, it->second});
76 } else {
77 LOG_ERROR(m_logger, "Species '{}' not found in global atomic species database.", name);
78 m_logger->flush_log();
79 throw std::runtime_error("Species not found in global atomic species database: " + std::string(name));
80 }
81 }
82
83 }
84
86 LOG_TRACE_L1(m_logger, "Populating reaction ID map for REACLIB graph network (serif::network::GraphNetwork)...");
87 m_reactionIDMap.clear();
88 for (auto& reaction: m_reactions) {
89 m_reactionIDMap.emplace(reaction.id(), &reaction);
90 }
91 LOG_TRACE_L1(m_logger, "Populated {} reactions in the reaction ID map.", m_reactionIDMap.size());
92 }
93
95 m_speciesToIndexMap.clear();
96 for (size_t i = 0; i < m_networkSpecies.size(); ++i) {
98 }
99 }
100
102 // The implementation of this function (and others) constrains this nuclear network to a constant temperature and density during
103 // each evaluation.
104 size_t numSpecies = m_networkSpecies.size();
105 m_jacobianMatrix.clear();
106 m_jacobianMatrix.resize(numSpecies, numSpecies, false); // Sparse matrix, no initial values
107 LOG_TRACE_L2(m_logger, "Jacobian matrix resized to {} rows and {} columns.",
108 m_jacobianMatrix.size1(), m_jacobianMatrix.size2());
109 }
110
111 // --- Basic Accessors and Queries ---
112 const std::vector<fourdst::atomic::Species>& GraphEngine::getNetworkSpecies() const {
113 // Returns a constant reference to the vector of unique species in the network.
114 LOG_TRACE_L3(m_logger, "Providing access to network species vector. Size: {}.", m_networkSpecies.size());
115 return m_networkSpecies;
116 }
117
119 // Returns a constant reference to the set of reactions in the network.
120 LOG_TRACE_L3(m_logger, "Providing access to network reactions set. Size: {}.", m_reactions.size());
121 return m_reactions;
122 }
123
124 bool GraphEngine::involvesSpecies(const fourdst::atomic::Species& species) const {
125 // Checks if a given species is present in the network's species map for efficient lookup.
126 const bool found = m_networkSpeciesMap.contains(species.name());
127 LOG_DEBUG(m_logger, "Checking if species '{}' is involved in the network: {}.", species.name(), found ? "Yes" : "No");
128 return found;
129 }
130
131 // --- Validation Methods ---
133 LOG_TRACE_L1(m_logger, "Validating mass (A) and charge (Z) conservation across all reactions in the network.");
134
135 for (const auto& reaction : m_reactions) {
136 uint64_t totalReactantA = 0;
137 uint64_t totalReactantZ = 0;
138 uint64_t totalProductA = 0;
139 uint64_t totalProductZ = 0;
140
141 // Calculate total A and Z for reactants
142 for (const auto& reactant : reaction.reactants()) {
143 auto it = m_networkSpeciesMap.find(reactant.name());
144 if (it != m_networkSpeciesMap.end()) {
145 totalReactantA += it->second.a();
146 totalReactantZ += it->second.z();
147 } else {
148 // This scenario indicates a severe data integrity issue:
149 // a reactant is part of a reaction but not in the network's species map.
150 LOG_ERROR(m_logger, "CRITICAL ERROR: Reactant species '{}' in reaction '{}' not found in network species map during conservation validation.",
151 reactant.name(), reaction.id());
152 return false;
153 }
154 }
155
156 // Calculate total A and Z for products
157 for (const auto& product : reaction.products()) {
158 auto it = m_networkSpeciesMap.find(product.name());
159 if (it != m_networkSpeciesMap.end()) {
160 totalProductA += it->second.a();
161 totalProductZ += it->second.z();
162 } else {
163 // Similar critical error for product species
164 LOG_ERROR(m_logger, "CRITICAL ERROR: Product species '{}' in reaction '{}' not found in network species map during conservation validation.",
165 product.name(), reaction.id());
166 return false;
167 }
168 }
169
170 // Compare totals for conservation
171 if (totalReactantA != totalProductA) {
172 LOG_ERROR(m_logger, "Mass number (A) not conserved for reaction '{}': Reactants A={} vs Products A={}.",
173 reaction.id(), totalReactantA, totalProductA);
174 return false;
175 }
176 if (totalReactantZ != totalProductZ) {
177 LOG_ERROR(m_logger, "Atomic number (Z) not conserved for reaction '{}': Reactants Z={} vs Products Z={}.",
178 reaction.id(), totalReactantZ, totalProductZ);
179 return false;
180 }
181 }
182
183 LOG_TRACE_L1(m_logger, "Mass (A) and charge (Z) conservation validated successfully for all reactions.");
184 return true; // All reactions passed the conservation check
185 }
186
187 void GraphEngine::validateComposition(const fourdst::composition::Composition &composition, double culling, double T9) {
188 // Check if the requested network has already been cached.
189 // PERF: Rebuilding this should be pretty fast but it may be a good point of optimization in the future.
190 const reaction::LogicalReactionSet validationReactionSet = build_reaclib_nuclear_network(composition, false);
191 // TODO: need some more robust method here to
192 // A. Build the basic network from the composition's species with non zero mass fractions.
193 // B. rebuild a new composition from the reaction set's reactants + products (with the mass fractions from the things that are only products set to 0)
194 // C. Rebuild the reaction set from the new composition
195 // D. Cull reactions where all reactants have mass fractions below the culling threshold.
196 // E. Be careful about maintaining caching through all of this
197
198
199 // This allows for dynamic network modification while retaining caching for networks which are very similar.
200 if (validationReactionSet != m_reactions) {
201 LOG_DEBUG(m_logger, "Reaction set not cached. Rebuilding the reaction set for T9={} and culling={}.", T9, culling);
202 m_reactions = std::move(validationReactionSet);
203 syncInternalMaps(); // Re-sync internal maps after updating reactions. Note this will also retrace the AD tape.
204 }
205 }
206
207 // --- Generate Stoichiometry Matrix ---
209 LOG_TRACE_L1(m_logger, "Generating stoichiometry matrix...");
210
211 // Task 1: Set dimensions and initialize the matrix
212 size_t numSpecies = m_networkSpecies.size();
213 size_t numReactions = m_reactions.size();
214 m_stoichiometryMatrix.resize(numSpecies, numReactions, false);
215
216 LOG_TRACE_L1(m_logger, "Stoichiometry matrix initialized with dimensions: {} rows (species) x {} columns (reactions).",
217 numSpecies, numReactions);
218
219 // Task 2: Populate the stoichiometry matrix
220 // Iterate through all reactions, assign them a column index, and fill in their stoichiometric coefficients.
221 size_t reactionColumnIndex = 0;
222 for (const auto& reaction : m_reactions) {
223 // Get the net stoichiometry for the current reaction
224 std::unordered_map<fourdst::atomic::Species, int> netStoichiometry = reaction.stoichiometry();
225
226 // Iterate through the species and their coefficients in the stoichiometry map
227 for (const auto& [species, coefficient] : netStoichiometry) {
228 // Find the row index for this species
229 auto it = m_speciesToIndexMap.find(species);
230 if (it != m_speciesToIndexMap.end()) {
231 const size_t speciesRowIndex = it->second;
232 // Set the matrix element. Boost.uBLAS handles sparse insertion.
233 m_stoichiometryMatrix(speciesRowIndex, reactionColumnIndex) = coefficient;
234 } else {
235 // This scenario should ideally not happen if m_networkSpeciesMap and m_speciesToIndexMap are correctly synced
236 LOG_ERROR(m_logger, "CRITICAL ERROR: Species '{}' from reaction '{}' stoichiometry not found in species to index map.",
237 species.name(), reaction.id());
238 m_logger -> flush_log();
239 throw std::runtime_error("Species not found in species to index map: " + std::string(species.name()));
240 }
241 }
242 reactionColumnIndex++; // Move to the next column for the next reaction
243 }
244
245 LOG_TRACE_L1(m_logger, "Stoichiometry matrix population complete. Number of non-zero elements: {}.",
246 m_stoichiometryMatrix.nnz()); // Assuming nnz() exists for compressed_matrix
247 }
248
250 const std::vector<double> &Y_in,
251 const double T9,
252 const double rho
253 ) const {
254 return calculateAllDerivatives<double>(Y_in, T9, rho);
255 }
256
258 const std::vector<ADDouble> &Y_in,
259 const ADDouble &T9,
260 const ADDouble &rho
261 ) const {
262 return calculateAllDerivatives<ADDouble>(Y_in, T9, rho);
263 }
264
267 const std::vector<double> &Y,
268 const double T9,
269 const double rho
270 ) const {
272 }
273
275 const std::vector<double> &Y,
276 const double T9,
277 const double rho
278 ) {
279
280 LOG_TRACE_L1(m_logger, "Generating jacobian matrix for T9={}, rho={}..", T9, rho);
281 const size_t numSpecies = m_networkSpecies.size();
282
283 // 1. Pack the input variables into a vector for CppAD
284 std::vector<double> adInput(numSpecies + 2, 0.0); // +2 for T9 and rho
285 for (size_t i = 0; i < numSpecies; ++i) {
286 adInput[i] = Y[i];
287 }
288 adInput[numSpecies] = T9; // T9
289 adInput[numSpecies + 1] = rho; // rho
290
291 // 2. Calculate the full jacobian
292 const std::vector<double> dotY = m_rhsADFun.Jacobian(adInput);
293
294 // 3. Pack jacobian vector into sparse matrix
295 m_jacobianMatrix.clear();
296 for (size_t i = 0; i < numSpecies; ++i) {
297 for (size_t j = 0; j < numSpecies; ++j) {
298 const double value = dotY[i * (numSpecies + 2) + j];
299 if (std::abs(value) > MIN_JACOBIAN_THRESHOLD) {
300 m_jacobianMatrix(i, j) = value;
301 }
302 }
303 }
304 LOG_TRACE_L1(m_logger, "Jacobian matrix generated with dimensions: {} rows x {} columns.", m_jacobianMatrix.size1(), m_jacobianMatrix.size2());
305 }
306
307 double GraphEngine::getJacobianMatrixEntry(const int i, const int j) const {
308 return m_jacobianMatrix(i, j);
309 }
310
311 std::unordered_map<fourdst::atomic::Species, int> GraphEngine::getNetReactionStoichiometry(
313 ) {
314 return reaction.stoichiometry();
315 }
316
318 const int speciesIndex,
319 const int reactionIndex
320 ) const {
321 return m_stoichiometryMatrix(speciesIndex, reactionIndex);
322 }
323
324 void GraphEngine::exportToDot(const std::string &filename) const {
325 LOG_TRACE_L1(m_logger, "Exporting network graph to DOT file: {}", filename);
326
327 std::ofstream dotFile(filename);
328 if (!dotFile.is_open()) {
329 LOG_ERROR(m_logger, "Failed to open file for writing: {}", filename);
330 m_logger->flush_log();
331 throw std::runtime_error("Failed to open file for writing: " + filename);
332 }
333
334 dotFile << "digraph NuclearReactionNetwork {\n";
335 dotFile << " graph [rankdir=LR, splines=true, overlap=false, bgcolor=\"#f0f0f0\"];\n";
336 dotFile << " node [shape=circle, style=filled, fillcolor=\"#a7c7e7\", fontname=\"Helvetica\"];\n";
337 dotFile << " edge [fontname=\"Helvetica\", fontsize=10];\n\n";
338
339 // 1. Define all species as nodes
340 dotFile << " // --- Species Nodes ---\n";
341 for (const auto& species : m_networkSpecies) {
342 dotFile << " \"" << species.name() << "\" [label=\"" << species.name() << "\"];\n";
343 }
344 dotFile << "\n";
345
346 // 2. Define all reactions as intermediate nodes and connect them
347 dotFile << " // --- Reaction Edges ---\n";
348 for (const auto& reaction : m_reactions) {
349 // Create a unique ID for the reaction node
350 std::string reactionNodeId = "reaction_" + std::string(reaction.id());
351
352 // Define the reaction node (small, black dot)
353 dotFile << " \"" << reactionNodeId << "\" [shape=point, fillcolor=black, width=0.1, height=0.1, label=\"\"];\n";
354
355 // Draw edges from reactants to the reaction node
356 for (const auto& reactant : reaction.reactants()) {
357 dotFile << " \"" << reactant.name() << "\" -> \"" << reactionNodeId << "\";\n";
358 }
359
360 // Draw edges from the reaction node to products
361 for (const auto& product : reaction.products()) {
362 dotFile << " \"" << reactionNodeId << "\" -> \"" << product.name() << "\" [label=\"" << reaction.qValue() << " MeV\"];\n";
363 }
364 dotFile << "\n";
365 }
366
367 dotFile << "}\n";
368 dotFile.close();
369 LOG_TRACE_L1(m_logger, "Successfully exported network to {}", filename);
370 }
371
372 void GraphEngine::exportToCSV(const std::string &filename) const {
373 LOG_TRACE_L1(m_logger, "Exporting network graph to CSV file: {}", filename);
374
375 std::ofstream csvFile(filename, std::ios::out | std::ios::trunc);
376 if (!csvFile.is_open()) {
377 LOG_ERROR(m_logger, "Failed to open file for writing: {}", filename);
378 m_logger->flush_log();
379 throw std::runtime_error("Failed to open file for writing: " + filename);
380 }
381 csvFile << "Reaction;Reactants;Products;Q-value;sources;rates\n";
382 for (const auto& reaction : m_reactions) {
383 // Dynamic cast to REACLIBReaction to access specific properties
384 csvFile << reaction.id() << ";";
385 // Reactants
386 int count = 0;
387 for (const auto& reactant : reaction.reactants()) {
388 csvFile << reactant.name();
389 if (++count < reaction.reactants().size()) {
390 csvFile << ",";
391 }
392 }
393 csvFile << ";";
394 count = 0;
395 for (const auto& product : reaction.products()) {
396 csvFile << product.name();
397 if (++count < reaction.products().size()) {
398 csvFile << ",";
399 }
400 }
401 csvFile << ";" << reaction.qValue() << ";";
402 // Reaction coefficients
403 auto sources = reaction.sources();
404 count = 0;
405 for (const auto& source : sources) {
406 csvFile << source;
407 if (++count < sources.size()) {
408 csvFile << ",";
409 }
410 }
411 csvFile << ";";
412 // Reaction coefficients
413 count = 0;
414 for (const auto& rates : reaction) {
415 csvFile << rates;
416 if (++count < reaction.size()) {
417 csvFile << ",";
418 }
419 }
420 csvFile << "\n";
421 }
422 csvFile.close();
423 LOG_TRACE_L1(m_logger, "Successfully exported network graph to {}", filename);
424 }
425
426 std::unordered_map<fourdst::atomic::Species, double> GraphEngine::getSpeciesTimescales(const std::vector<double> &Y, const double T9,
427 const double rho) const {
428 auto [dydt, _] = calculateAllDerivatives<double>(Y, T9, rho);
429 std::unordered_map<fourdst::atomic::Species, double> speciesTimescales;
430 speciesTimescales.reserve(m_networkSpecies.size());
431 for (size_t i = 0; i < m_networkSpecies.size(); ++i) {
432 double timescale = std::numeric_limits<double>::infinity();
433 const auto species = m_networkSpecies[i];
434 if (std::abs(dydt[i]) > 0.0) {
435 timescale = std::abs(Y[i] / dydt[i]);
436 }
437 speciesTimescales.emplace(species, timescale);
438 }
439 return speciesTimescales;
440 }
441
443 LOG_TRACE_L1(m_logger, "Recording AD tape for the RHS calculation...");
444
445 // Task 1: Set dimensions and initialize the matrix
446 const size_t numSpecies = m_networkSpecies.size();
447 if (numSpecies == 0) {
448 LOG_ERROR(m_logger, "Cannot record AD tape: No species in the network.");
449 m_logger->flush_log();
450 throw std::runtime_error("Cannot record AD tape: No species in the network.");
451 }
452 const size_t numADInputs = numSpecies + 2; // Note here that by not letting T9 and rho be independent variables, we are constraining the network to a constant temperature and density during each evaluation.
453
454 // --- CppAD Tape Recording ---
455 // 1. Declare independent variable (adY)
456 // We also initialize the dummy variable for tape recording (these tell CppAD what the derivative chain looks like).
457 // Their numeric values are irrelevant except for in so far as they avoid numerical instabilities.
458
459 // Distribute total mass fraction uniformly between species in the dummy variable space
460 const auto uniformMassFraction = static_cast<CppAD::AD<double>>(1.0 / numSpecies);
461 std::vector<CppAD::AD<double>> adInput(numADInputs, uniformMassFraction);
462 adInput[numSpecies] = 1.0; // Dummy T9
463 adInput[numSpecies + 1] = 1.0; // Dummy rho
464
465 // 3. Declare independent variables (what CppAD will differentiate wrt.)
466 // This also beings the tape recording process.
467 CppAD::Independent(adInput);
468
469 std::vector<CppAD::AD<double>> adY(numSpecies);
470 for(size_t i = 0; i < numSpecies; ++i) {
471 adY[i] = adInput[i];
472 }
473 const CppAD::AD<double> adT9 = adInput[numSpecies];
474 const CppAD::AD<double> adRho = adInput[numSpecies + 1];
475
476
477 // 5. Call the actual templated function
478 // We let T9 and rho be constant, so we pass them as fixed values.
479 auto [dydt, nuclearEnergyGenerationRate] = calculateAllDerivatives<CppAD::AD<double>>(adY, adT9, adRho);
480
481 m_rhsADFun.Dependent(adInput, dydt);
482
483 LOG_TRACE_L1(m_logger, "AD tape recorded successfully for the RHS calculation. Number of independent variables: {}.",
484 adInput.size());
485 }
486}
const std::vector< fourdst::atomic::Species > & getNetworkSpecies() const override
Gets the list of species in the network.
std::unordered_map< fourdst::atomic::Species, double > getSpeciesTimescales(const std::vector< double > &Y, double T9, double rho) const override
Computes timescales for all species in the network.
void populateReactionIDMap()
Populates the reaction ID map.
CppAD::ADFun< double > m_rhsADFun
CppAD function for the right-hand side of the ODE.
boost::numeric::ublas::compressed_matrix< double > m_jacobianMatrix
Jacobian matrix (species x species).
double getJacobianMatrixEntry(const int i, const int j) const override
Gets an entry from the previously generated Jacobian matrix.
std::unordered_map< std::string_view, fourdst::atomic::Species > m_networkSpeciesMap
Map from species name to Species object.
void populateSpeciesToIndexMap()
Populates the species-to-index map.
quill::Logger * m_logger
void reserveJacobianMatrix()
Reserves space for the Jacobian matrix.
std::unordered_map< std::string_view, reaction::Reaction * > m_reactionIDMap
Map from reaction ID to REACLIBReaction. //PERF: This makes copies of REACLIBReaction and could be a ...
int getStoichiometryMatrixEntry(const int speciesIndex, const int reactionIndex) const override
Gets an entry from the stoichiometry matrix.
void exportToCSV(const std::string &filename) const
Exports the network to a CSV file for analysis.
StepDerivatives< double > calculateRHSAndEnergy(const std::vector< double > &Y, const double T9, const double rho) const override
Calculates the right-hand side (dY/dt) and energy generation rate.
static std::unordered_map< fourdst::atomic::Species, int > getNetReactionStoichiometry(const reaction::Reaction &reaction)
Gets the net stoichiometry for a given reaction.
double calculateMolarReactionFlow(const reaction::Reaction &reaction, const std::vector< double > &Y, const double T9, const double rho) const override
Calculates the molar reaction flow for a given reaction.
std::vector< fourdst::atomic::Species > m_networkSpecies
Vector of unique species in the network.
void recordADTape()
Records the AD tape for the right-hand side of the ODE.
GraphEngine(const fourdst::composition::Composition &composition)
Constructs a GraphEngine from a composition.
bool involvesSpecies(const fourdst::atomic::Species &species) const
Checks if a given species is involved in the network.
reaction::LogicalReactionSet m_reactions
Set of REACLIB reactions in the network.
void syncInternalMaps()
Synchronizes the internal maps.
bool validateConservation() const
Validates mass and charge conservation across all reactions.
boost::numeric::ublas::compressed_matrix< int > m_stoichiometryMatrix
Stoichiometry matrix (species x reactions).
const reaction::LogicalReactionSet & getNetworkReactions() const override
Gets the set of logical reactions in the network.
std::unordered_map< fourdst::atomic::Species, size_t > m_speciesToIndexMap
Map from species to their index in the stoichiometry matrix.
void exportToDot(const std::string &filename) const
Exports the network to a DOT file for visualization.
void generateJacobianMatrix(const std::vector< double > &Y, const double T9, const double rho) override
Generates the Jacobian matrix for the current state.
void generateStoichiometryMatrix() override
Generates the stoichiometry matrix for the network.
void collectNetworkSpecies()
Collects the unique species in the network.
void validateComposition(const fourdst::composition::Composition &composition, double culling, double T9)
Validates the composition against the current reaction set.
StepDerivatives< T > calculateAllDerivatives(const std::vector< T > &Y_in, T T9, T rho) const
Calculates all derivatives (dY/dt) and the energy generation rate.
A collection of LogicalReaction objects.
Definition reaction.h:554
Represents a single nuclear reaction from a specific data source.
Definition reaction.h:71
CppAD::AD< double > ADDouble
Alias for CppAD AD type for double precision.
reaction::LogicalReactionSet build_reaclib_nuclear_network(const fourdst::composition::Composition &composition, bool reverse)
Definition network.cpp:64
static constexpr double MIN_JACOBIAN_THRESHOLD
Minimum value for Jacobian matrix entries.
STL namespace.
Defines classes for representing and managing nuclear reactions.
Structure holding derivatives and energy generation for a network step.