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"
13
14#include <string>
15#include <unordered_map>
16#include <vector>
17#include <memory>
18
19#include <boost/numeric/ublas/matrix_sparse.hpp>
20
21#include "cppad/cppad.hpp"
22
23// PERF: The function getNetReactionStoichiometry returns a map of species to their stoichiometric coefficients for a given reaction.
24// this makes extra copies of the species, which is not ideal and could be optimized further.
25// Even more relevant is the member m_reactionIDMap which makes copies of a REACLIBReaction for each reaction ID.
26// REACLIBReactions are quite large data structures, so this could be a performance bottleneck.
27
28namespace gridfire {
34 typedef CppAD::AD<double> ADDouble;
35
36 using fourdst::config::Config;
37 using fourdst::logging::LogManager;
38 using fourdst::constant::Constants;
39
47 static constexpr double MIN_DENSITY_THRESHOLD = 1e-18;
48
56 static constexpr double MIN_ABUNDANCE_THRESHOLD = 1e-18;
57
64 static constexpr double MIN_JACOBIAN_THRESHOLD = 1e-24;
65
89 class GraphEngine final : public DynamicEngine{
90 public:
101 explicit GraphEngine(const fourdst::composition::Composition &composition);
102
111 explicit GraphEngine(reaction::LogicalReactionSet reactions);
112
127 const std::vector<double>& Y,
128 const double T9,
129 const double rho
130 ) const override;
131
146 const std::vector<double>& Y,
147 const double T9,
148 const double rho
149 ) override;
150
157 void generateStoichiometryMatrix() override;
158
171 [[nodiscard]] double calculateMolarReactionFlow(
173 const std::vector<double>&Y,
174 const double T9,
175 const double rho
176 ) const override;
177
182 [[nodiscard]] const std::vector<fourdst::atomic::Species>& getNetworkSpecies() const override;
183
188 [[nodiscard]] const reaction::LogicalReactionSet& getNetworkReactions() const override;
189
201 [[nodiscard]] double getJacobianMatrixEntry(
202 const int i,
203 const int j
204 ) const override;
205
212 [[nodiscard]] static std::unordered_map<fourdst::atomic::Species, int> getNetReactionStoichiometry(
214 );
215
227 [[nodiscard]] int getStoichiometryMatrixEntry(
228 const int speciesIndex,
229 const int reactionIndex
230 ) const override;
231
243 [[nodiscard]] std::unordered_map<fourdst::atomic::Species, double> getSpeciesTimescales(
244 const std::vector<double>& Y,
245 double T9,
246 double rho
247 ) const override;
248
249 void update(const NetIn& netIn) override;
250
257 [[nodiscard]] bool involvesSpecies(
258 const fourdst::atomic::Species& species
259 ) const;
260
277 void exportToDot(
278 const std::string& filename
279 ) const;
280
297 void exportToCSV(
298 const std::string& filename
299 ) const;
300
302
303 [[nodiscard]] screening::ScreeningType getScreeningModel() const override;
304
305
306 private:
308 std::unordered_map<std::string_view, reaction::Reaction*> m_reactionIDMap;
309
310 std::vector<fourdst::atomic::Species> m_networkSpecies;
311 std::unordered_map<std::string_view, fourdst::atomic::Species> m_networkSpeciesMap;
312 std::unordered_map<fourdst::atomic::Species, size_t> m_speciesToIndexMap;
313
314 boost::numeric::ublas::compressed_matrix<int> m_stoichiometryMatrix;
315 boost::numeric::ublas::compressed_matrix<double> m_jacobianMatrix;
316
317 CppAD::ADFun<double> m_rhsADFun;
318
320 std::unique_ptr<screening::ScreeningModel> m_screeningModel = screening::selectScreeningModel(m_screeningType);
321
322 Config& m_config = Config::getInstance();
323 Constants& m_constants = Constants::getInstance();
324 quill::Logger* m_logger = LogManager::getInstance().getLogger("log");
325
326 private:
334 void syncInternalMaps();
335
343
351
359
368
378 void recordADTape();
379
389 [[nodiscard]] bool validateConservation() const;
390
403 const fourdst::composition::Composition &composition,
404 double culling,
405 double T9
406 );
407
421 template <IsArithmeticOrAD T>
424 const std::vector<T> &Y,
425 const T T9,
426 const T rho
427 ) const;
428
441 template<IsArithmeticOrAD T>
443 const std::vector<T> &Y_in,
444 T T9,
445 T rho
446 ) const;
447
461 const std::vector<double>& Y_in,
462 const double T9,
463 const double rho
464 ) const;
465
479 const std::vector<ADDouble>& Y_in,
480 const ADDouble &T9,
481 const ADDouble &rho
482 ) const;
483 };
484
485
486 template<IsArithmeticOrAD T>
488 const std::vector<T> &Y_in, T T9, T rho) const {
489 std::vector<T> screeningFactors = m_screeningModel->calculateScreeningFactors(
492 Y_in,
493 T9,
494 rho
495 );
496
497 // --- Setup output derivatives structure ---
498 StepDerivatives<T> result;
499 result.dydt.resize(m_networkSpecies.size(), static_cast<T>(0.0));
500
501 // --- AD Pre-setup (flags to control conditionals in an AD safe / branch aware manner) ---
502 // ----- Constants for AD safe calculations ---
503 const T zero = static_cast<T>(0.0);
504 const T one = static_cast<T>(1.0);
505
506 // ----- Initialize variables for molar concentration product and thresholds ---
507 // Note: the logic here is that we use CppAD::CondExprLt to test thresholds and if they are less we set the flag
508 // to zero so that the final returned reaction flow is 0. This is as opposed to standard if statements
509 // which create branches that break the AD tape.
510 const T rho_threshold = static_cast<T>(MIN_DENSITY_THRESHOLD);
511
512 // --- Check if the density is below the threshold where we ignore reactions ---
513 T threshold_flag = CppAD::CondExpLt(rho, rho_threshold, zero, one); // If rho < threshold, set flag to 0
514
515 std::vector<T> Y = Y_in;
516 for (size_t i = 0; i < m_networkSpecies.size(); ++i) {
517 // We use CppAD::CondExpLt to handle AD taping and prevent branching
518 // Note that while this is syntactically more complex this is equivalent to
519 // if (Y[i] < 0) {Y[i] = 0;}
520 // The issue is that this would introduce a branch which would require the auto diff tape to be re-recorded
521 // each timestep, which is very inefficient.
522 Y[i] = CppAD::CondExpLt(Y[i], zero, zero, Y[i]); // Ensure no negative abundances
523 }
524
525 const T u = static_cast<T>(m_constants.get("u").value); // Atomic mass unit in grams
526 const T N_A = static_cast<T>(m_constants.get("N_a").value); // Avogadro's number in mol^-1
527 const T c = static_cast<T>(m_constants.get("c").value); // Speed of light in cm/s
528
529 // --- SINGLE LOOP OVER ALL REACTIONS ---
530 for (size_t reactionIndex = 0; reactionIndex < m_reactions.size(); ++reactionIndex) {
531 const auto& reaction = m_reactions[reactionIndex];
532
533 // 1. Calculate reaction rate
534 const T molarReactionFlow = screeningFactors[reactionIndex] * calculateMolarReactionFlow<T>(reaction, Y, T9, rho);
535
536 // 2. Use the rate to update all relevant species derivatives (dY/dt)
537 for (size_t speciesIndex = 0; speciesIndex < m_networkSpecies.size(); ++speciesIndex) {
538 const T nu_ij = static_cast<T>(m_stoichiometryMatrix(speciesIndex, reactionIndex));
539 result.dydt[speciesIndex] += threshold_flag * nu_ij * molarReactionFlow / rho;
540 }
541 }
542
543 T massProductionRate = static_cast<T>(0.0); // [mol][s^-1]
544 for (const auto& [species, index] : m_speciesToIndexMap) {
545 massProductionRate += result.dydt[index] * species.mass() * u;
546 }
547
548 result.nuclearEnergyGenerationRate = -massProductionRate * N_A * c * c; // [cm^2][s^-3] = [erg][s^-1][g^-1]
549
550 return result;
551 }
552
553
554 template <IsArithmeticOrAD T>
557 const std::vector<T> &Y,
558 const T T9,
559 const T rho
560 ) const {
561
562 // --- Pre-setup (flags to control conditionals in an AD safe / branch aware manner) ---
563 // ----- Constants for AD safe calculations ---
564 const T zero = static_cast<T>(0.0);
565 const T one = static_cast<T>(1.0);
566
567 // ----- Initialize variables for molar concentration product and thresholds ---
568 // Note: the logic here is that we use CppAD::CondExprLt to test thresholds and if they are less we set the flag
569 // to zero so that the final returned reaction flow is 0. This is as opposed to standard if statements
570 // which create branches that break the AD tape.
571 const T Y_threshold = static_cast<T>(MIN_ABUNDANCE_THRESHOLD);
572 T threshold_flag = one;
573
574 // --- Calculate the molar reaction rate (in units of [s^-1][cm^3(N-1)][mol^(1-N)] for N reactants) ---
575 const T k_reaction = reaction.calculate_rate(T9);
576
577 // --- Cound the number of each reactant species to account for species multiplicity ---
578 std::unordered_map<std::string, int> reactant_counts;
579 reactant_counts.reserve(reaction.reactants().size());
580 for (const auto& reactant : reaction.reactants()) {
581 reactant_counts[std::string(reactant.name())]++;
582 }
583
584 // --- Accumulator for the molar concentration ---
585 auto molar_concentration_product = static_cast<T>(1.0);
586
587 // --- Loop through each unique reactant species and calculate the molar concentration for that species then multiply that into the accumulator ---
588 for (const auto& [species_name, count] : reactant_counts) {
589 // --- Resolve species to molar abundance ---
590 // PERF: Could probably optimize out this lookup
591 const auto species_it = m_speciesToIndexMap.find(m_networkSpeciesMap.at(species_name));
592 const size_t species_index = species_it->second;
593 const T Yi = Y[species_index];
594
595 // --- Check if the species abundance is below the threshold where we ignore reactions ---
596 threshold_flag *= CppAD::CondExpLt(Yi, Y_threshold, zero, one);
597
598 // --- Convert from molar abundance to molar concentration ---
599 T molar_concentration = Yi * rho;
600
601 // --- 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 ---
602 molar_concentration_product *= CppAD::pow(molar_concentration, static_cast<T>(count)); // ni^count
603
604 // --- Apply factorial correction for identical reactions ---
605 if (count > 1) {
606 molar_concentration_product /= static_cast<T>(std::tgamma(static_cast<double>(count + 1))); // Gamma function for factorial
607 }
608 }
609 // --- Final reaction flow calculation [mol][s^-1][cm^-3] ---
610 // Note: If the threshold flag ever gets set to zero this will return zero.
611 // This will result basically in multiple branches being written to the AD tape, which will make
612 // the tape more expensive to record, but it will also mean that we only need to record it once for
613 // the entire network.
614 return molar_concentration_product * k_reaction * threshold_flag;
615 }
616};
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 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
Abstract interfaces for reaction network engines in GridFire.
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.
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).