GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
engine_graph.h
Go to the documentation of this file.
1#pragma once
2
3#include "fourdst/composition/atomicSpecies.h"
4#include "fourdst/composition/composition.h"
5#include "fourdst/logging/logging.h"
6#include "fourdst/config/config.h"
7
8#include "gridfire/network.h"
11
12#include <string>
13#include <unordered_map>
14#include <vector>
15
16#include <boost/numeric/ublas/matrix_sparse.hpp>
17
18#include "cppad/cppad.hpp"
19
20// PERF: The function getNetReactionStoichiometry returns a map of species to their stoichiometric coefficients for a given reaction.
21// this makes extra copies of the species, which is not ideal and could be optimized further.
22// Even more relevant is the member m_reactionIDMap which makes copies of a REACLIBReaction for each reaction ID.
23// REACLIBReactions are quite large data structures, so this could be a performance bottleneck.
24
25namespace gridfire {
31 typedef CppAD::AD<double> ADDouble;
32
33 using fourdst::config::Config;
34 using fourdst::logging::LogManager;
35 using fourdst::constant::Constants;
36
44 static constexpr double MIN_DENSITY_THRESHOLD = 1e-18;
45
53 static constexpr double MIN_ABUNDANCE_THRESHOLD = 1e-18;
54
61 static constexpr double MIN_JACOBIAN_THRESHOLD = 1e-24;
62
86 class GraphEngine final : public DynamicEngine{
87 public:
98 explicit GraphEngine(const fourdst::composition::Composition &composition);
99
108 explicit GraphEngine(reaction::LogicalReactionSet reactions);
109
124 const std::vector<double>& Y,
125 const double T9,
126 const double rho
127 ) const override;
128
143 const std::vector<double>& Y,
144 const double T9,
145 const double rho
146 ) override;
147
154 void generateStoichiometryMatrix() override;
155
170 const std::vector<double>&Y,
171 const double T9,
172 const double rho
173 ) const override;
174
179 [[nodiscard]] const std::vector<fourdst::atomic::Species>& getNetworkSpecies() const override;
180
185 [[nodiscard]] const reaction::LogicalReactionSet& getNetworkReactions() const override;
186
198 [[nodiscard]] double getJacobianMatrixEntry(
199 const int i,
200 const int j
201 ) const override;
202
209 [[nodiscard]] static std::unordered_map<fourdst::atomic::Species, int> getNetReactionStoichiometry(
211 );
212
224 [[nodiscard]] int getStoichiometryMatrixEntry(
225 const int speciesIndex,
226 const int reactionIndex
227 ) const override;
228
240 [[nodiscard]] std::unordered_map<fourdst::atomic::Species, double> getSpeciesTimescales(
241 const std::vector<double>& Y,
242 double T9,
243 double rho
244 ) const override;
245
252 [[nodiscard]] bool involvesSpecies(
253 const fourdst::atomic::Species& species
254 ) const;
255
272 void exportToDot(
273 const std::string& filename
274 ) const;
275
292 void exportToCSV(
293 const std::string& filename
294 ) const;
295
296
297 private:
299 std::unordered_map<std::string_view, reaction::Reaction*> m_reactionIDMap;
300
301 std::vector<fourdst::atomic::Species> m_networkSpecies;
302 std::unordered_map<std::string_view, fourdst::atomic::Species> m_networkSpeciesMap;
303 std::unordered_map<fourdst::atomic::Species, size_t> m_speciesToIndexMap;
304
305 boost::numeric::ublas::compressed_matrix<int> m_stoichiometryMatrix;
306 boost::numeric::ublas::compressed_matrix<double> m_jacobianMatrix;
307
308 CppAD::ADFun<double> m_rhsADFun;
309
310 Config& m_config = Config::getInstance();
311 Constants& m_constants = Constants::getInstance();
312 quill::Logger* m_logger = LogManager::getInstance().getLogger("log");
313
314 private:
322 void syncInternalMaps();
323
331
339
347
356
366 void recordADTape();
367
377 [[nodiscard]] bool validateConservation() const;
378
391 const fourdst::composition::Composition &composition,
392 double culling,
393 double T9
394 );
395
409 template <IsArithmeticOrAD T>
412 const std::vector<T> &Y,
413 const T T9,
414 const T rho
415 ) const;
416
429 template<IsArithmeticOrAD T>
431 const std::vector<T> &Y_in,
432 T T9,
433 T rho
434 ) const;
435
449 const std::vector<double>& Y_in,
450 const double T9,
451 const double rho
452 ) const;
453
467 const std::vector<ADDouble>& Y_in,
468 const ADDouble &T9,
469 const ADDouble &rho
470 ) const;
471 };
472
473
474 template<IsArithmeticOrAD T>
476 const std::vector<T> &Y_in, T T9, T rho) const {
477
478 // --- Setup output derivatives structure ---
479 StepDerivatives<T> result;
480 result.dydt.resize(m_networkSpecies.size(), static_cast<T>(0.0));
481
482 // --- AD Pre-setup (flags to control conditionals in an AD safe / branch aware manner) ---
483 // ----- Constants for AD safe calculations ---
484 const T zero = static_cast<T>(0.0);
485 const T one = static_cast<T>(1.0);
486
487 // ----- Initialize variables for molar concentration product and thresholds ---
488 // Note: the logic here is that we use CppAD::CondExprLt to test thresholds and if they are less we set the flag
489 // to zero so that the final returned reaction flow is 0. This is as opposed to standard if statements
490 // which create branches that break the AD tape.
491 const T rho_threshold = static_cast<T>(MIN_DENSITY_THRESHOLD);
492
493 // --- Check if the density is below the threshold where we ignore reactions ---
494 T threshold_flag = CppAD::CondExpLt(rho, rho_threshold, zero, one); // If rho < threshold, set flag to 0
495
496 std::vector<T> Y = Y_in;
497 for (size_t i = 0; i < m_networkSpecies.size(); ++i) {
498 // We use CppAD::CondExpLt to handle AD taping and prevent branching
499 // Note that while this is syntactically more complex this is equivalent to
500 // if (Y[i] < 0) {Y[i] = 0;}
501 // The issue is that this would introduce a branch which would require the auto diff tape to be re-recorded
502 // each timestep, which is very inefficient.
503 Y[i] = CppAD::CondExpLt(Y[i], zero, zero, Y[i]); // Ensure no negative abundances
504 }
505
506 const T u = static_cast<T>(m_constants.get("u").value); // Atomic mass unit in grams
507 const T N_A = static_cast<T>(m_constants.get("N_a").value); // Avogadro's number in mol^-1
508 const T c = static_cast<T>(m_constants.get("c").value); // Speed of light in cm/s
509
510 // --- SINGLE LOOP OVER ALL REACTIONS ---
511 for (size_t reactionIndex = 0; reactionIndex < m_reactions.size(); ++reactionIndex) {
512 const auto& reaction = m_reactions[reactionIndex];
513
514 // 1. Calculate reaction rate
515 const T molarReactionFlow = calculateMolarReactionFlow<T>(reaction, Y, T9, rho);
516
517 // 2. Use the rate to update all relevant species derivatives (dY/dt)
518 for (size_t speciesIndex = 0; speciesIndex < m_networkSpecies.size(); ++speciesIndex) {
519 const T nu_ij = static_cast<T>(m_stoichiometryMatrix(speciesIndex, reactionIndex));
520 result.dydt[speciesIndex] += threshold_flag * nu_ij * molarReactionFlow / rho;
521 }
522 }
523
524 T massProductionRate = static_cast<T>(0.0); // [mol][s^-1]
525 for (const auto& [species, index] : m_speciesToIndexMap) {
526 massProductionRate += result.dydt[index] * species.mass() * u;
527 }
528
529 result.nuclearEnergyGenerationRate = -massProductionRate * N_A * c * c; // [cm^2][s^-3] = [erg][s^-1][g^-1]
530
531 return result;
532 }
533
534
535 template <IsArithmeticOrAD T>
538 const std::vector<T> &Y,
539 const T T9,
540 const T rho
541 ) const {
542
543 // --- Pre-setup (flags to control conditionals in an AD safe / branch aware manner) ---
544 // ----- Constants for AD safe calculations ---
545 const T zero = static_cast<T>(0.0);
546 const T one = static_cast<T>(1.0);
547
548 // ----- Initialize variables for molar concentration product and thresholds ---
549 // Note: the logic here is that we use CppAD::CondExprLt to test thresholds and if they are less we set the flag
550 // to zero so that the final returned reaction flow is 0. This is as opposed to standard if statements
551 // which create branches that break the AD tape.
552 const T Y_threshold = static_cast<T>(MIN_ABUNDANCE_THRESHOLD);
553 T threshold_flag = one;
554
555 // --- Calculate the molar reaction rate (in units of [s^-1][cm^3(N-1)][mol^(1-N)] for N reactants) ---
556 const T k_reaction = reaction.calculate_rate(T9);
557
558 // --- Cound the number of each reactant species to account for species multiplicity ---
559 std::unordered_map<std::string, int> reactant_counts;
560 reactant_counts.reserve(reaction.reactants().size());
561 for (const auto& reactant : reaction.reactants()) {
562 reactant_counts[std::string(reactant.name())]++;
563 }
564
565 // --- Accumulator for the molar concentration ---
566 auto molar_concentration_product = static_cast<T>(1.0);
567
568 // --- Loop through each unique reactant species and calculate the molar concentration for that species then multiply that into the accumulator ---
569 for (const auto& [species_name, count] : reactant_counts) {
570 // --- Resolve species to molar abundance ---
571 // PERF: Could probably optimize out this lookup
572 const auto species_it = m_speciesToIndexMap.find(m_networkSpeciesMap.at(species_name));
573 const size_t species_index = species_it->second;
574 const T Yi = Y[species_index];
575
576 // --- Check if the species abundance is below the threshold where we ignore reactions ---
577 threshold_flag *= CppAD::CondExpLt(Yi, Y_threshold, zero, one);
578
579 // --- Convert from molar abundance to molar concentration ---
580 T molar_concentration = Yi * rho;
581
582 // --- If count is > 1 , we need to raise the molar concentration to the power of count since there are really count bodies in that reaction ---
583 molar_concentration_product *= CppAD::pow(molar_concentration, static_cast<T>(count)); // ni^count
584
585 // --- Apply factorial correction for identical reactions ---
586 if (count > 1) {
587 molar_concentration_product /= static_cast<T>(std::tgamma(static_cast<double>(count + 1))); // Gamma function for factorial
588 }
589 }
590 // --- Final reaction flow calculation [mol][s^-1][cm^-3] ---
591 // Note: If the threshold flag ever gets set to zero this will return zero.
592 // This will result basically in multiple branches being written to the AD tape, which will make
593 // the tape more expensive to record, but it will also mean that we only need to record it once for
594 // the entire network.
595 return molar_concentration_product * k_reaction * threshold_flag;
596 }
597};
Abstract class for engines supporting Jacobian and stoichiometry operations.
Constants & m_constants
Access to physical constants.
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
Abstract interfaces for reaction network engines in GridFire.
CppAD::AD< double > ADDouble
Alias for CppAD AD type for double precision.
static constexpr double MIN_ABUNDANCE_THRESHOLD
Minimum abundance threshold below which species are ignored.
static constexpr double MIN_DENSITY_THRESHOLD
Minimum density threshold below which reactions are ignored.
static constexpr double MIN_JACOBIAN_THRESHOLD
Minimum value for Jacobian matrix entries.
Defines classes for representing and managing nuclear reactions.
Structure holding derivatives and energy generation for a network step.
T nuclearEnergyGenerationRate
Specific energy generation rate (e.g., erg/g/s).
std::vector< T > dydt
Derivatives of abundances (dY/dt for each species).