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