GridFire v0.7.0-alpha
General Purpose Nuclear Network
Loading...
Searching...
No Matches
gridfire::GraphEngine Class Referencefinal

A reaction network engine that uses a graph-based representation. More...

#include <engine_graph.h>

Inheritance diagram for gridfire::GraphEngine:
gridfire::DynamicEngine gridfire::Engine

Classes

class  AtomicReverseRate
 
struct  constants
 
struct  PrecomputedReaction
 

Public Member Functions

 GraphEngine (const fourdst::composition::Composition &composition, BuildDepthType=NetworkBuildDepth::Full)
 Constructs a GraphEngine from a composition.
 
 GraphEngine (const fourdst::composition::Composition &composition, const partition::PartitionFunction &partitionFunction, BuildDepthType buildDepth=NetworkBuildDepth::Full)
 
 GraphEngine (const fourdst::composition::Composition &composition, const partition::PartitionFunction &partitionFunction, BuildDepthType buildDepth, NetworkConstructionFlags reactionTypes)
 
 GraphEngine (const reaction::ReactionSet &reactions)
 Constructs a GraphEngine from a set of reactions.
 
std::expected< StepDerivatives< double >, expectations::StaleEngineErrorcalculateRHSAndEnergy (const fourdst::composition::Composition &comp, double T9, double rho) const override
 Calculates the right-hand side (dY/dt) and energy generation rate.
 
std::expected< StepDerivatives< double >, expectations::StaleEngineErrorcalculateRHSAndEnergy (const fourdst::composition::Composition &comp, double T9, double rho, const reaction::ReactionSet &activeReactions) const
 Calculates the right-hand side (dY/dt) and energy generation rate for a subset of reactions.
 
EnergyDerivatives calculateEpsDerivatives (const fourdst::composition::Composition &comp, double T9, double rho) const override
 Calculates the derivatives of the energy generation rate with respect to temperature and density.
 
EnergyDerivatives calculateEpsDerivatives (const fourdst::composition::Composition &comp, double T9, double rho, const reaction::ReactionSet &activeReactions) const
 Calculates the derivatives of the energy generation rate with respect to temperature and density for a subset of reactions.
 
void generateJacobianMatrix (const fourdst::composition::Composition &comp, double T9, double rho) const override
 Generates the Jacobian matrix for the current state.
 
void generateJacobianMatrix (const fourdst::composition::Composition &comp, double T9, double rho, const std::vector< fourdst::atomic::Species > &activeSpecies) const override
 Generates the Jacobian matrix for the current state with a specified set of active species. generally this will be much faster than the full matrix generation. Here we use forward mode to generate the Jacobian only for the active species.
 
void generateJacobianMatrix (const fourdst::composition::Composition &comp, double T9, double rho, const SparsityPattern &sparsityPattern) const override
 Generates the Jacobian matrix for the current state with a specified sparsity pattern.
 
void generateStoichiometryMatrix () override
 Generates the stoichiometry matrix for the network.
 
double calculateMolarReactionFlow (const reaction::Reaction &reaction, const fourdst::composition::Composition &comp, double T9, double rho) const override
 Calculates the molar reaction flow for a given reaction.
 
const std::vector< fourdst::atomic::Species > & getNetworkSpecies () const override
 Gets the list of species in the network.
 
const reaction::ReactionSetgetNetworkReactions () const override
 Gets the set of logical reactions in the network.
 
void setNetworkReactions (const reaction::ReactionSet &reactions) override
 Sets the reactions for the network.
 
double getJacobianMatrixEntry (const fourdst::atomic::Species &rowSpecies, const fourdst::atomic::Species &colSpecies) const override
 Gets an entry from the previously generated Jacobian matrix.
 
int getStoichiometryMatrixEntry (const fourdst::atomic::Species &species, const reaction::Reaction &reaction) const override
 Gets an entry from the stoichiometry matrix.
 
std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineErrorgetSpeciesTimescales (const fourdst::composition::Composition &comp, double T9, double rho) const override
 Computes timescales for all species in the network.
 
std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineErrorgetSpeciesTimescales (const fourdst::composition::Composition &comp, double T9, double rho, const reaction::ReactionSet &activeReactions) const
 Computes timescales for all species in the network considering a subset of reactions.
 
std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineErrorgetSpeciesDestructionTimescales (const fourdst::composition::Composition &comp, double T9, double rho) const override
 Computes destruction timescales for all species in the network.
 
std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineErrorgetSpeciesDestructionTimescales (const fourdst::composition::Composition &comp, double T9, double rho, const reaction::ReactionSet &activeReactions) const
 Computes destruction timescales for all species in the network considering a subset of reactions.
 
fourdst::composition::Composition update (const NetIn &netIn) override
 Updates the state of the network and the composition to be usable for the current network.
 
bool isStale (const NetIn &netIn) override
 Checks if the engine view is stale and needs to be updated.
 
bool involvesSpecies (const fourdst::atomic::Species &species) const
 Checks if a given species is involved in the network.
 
void exportToDot (const std::string &filename) const
 Exports the network to a DOT file for visualization.
 
void exportToCSV (const std::string &filename) const
 Exports the network to a CSV file for analysis.
 
void setScreeningModel (screening::ScreeningType model) override
 Sets the electron screening model for reaction rate calculations.
 
screening::ScreeningType getScreeningModel () const override
 Gets the current electron screening model.
 
void setPrecomputation (bool precompute)
 Sets whether to precompute reaction rates.
 
bool isPrecomputationEnabled () const
 Checks if precomputation of reaction rates is enabled.
 
const partition::PartitionFunctiongetPartitionFunction () const
 Gets the partition function used for reaction rate calculations.
 
double calculateReverseRate (const reaction::Reaction &reaction, double T9, double rho, const fourdst::composition::Composition &comp) const
 Calculates the reverse rate for a given reaction.
 
double calculateReverseRateTwoBody (const reaction::Reaction &reaction, double T9, double forwardRate, double expFactor) const
 Calculates the reverse rate for a two-body reaction.
 
double calculateReverseRateTwoBodyDerivative (const reaction::Reaction &reaction, double T9, double rho, const fourdst::composition::Composition &comp, double reverseRate) const
 Calculates the derivative of the reverse rate for a two-body reaction with respect to temperature.
 
bool isUsingReverseReactions () const
 Checks if reverse reactions are enabled.
 
void setUseReverseReactions (bool useReverse)
 Sets whether to use reverse reactions in the engine.
 
size_t getSpeciesIndex (const fourdst::atomic::Species &species) const override
 Gets the index of a species in the network.
 
std::vector< double > mapNetInToMolarAbundanceVector (const NetIn &netIn) const override
 Maps the NetIn object to a vector of molar abundances.
 
PrimingReport primeEngine (const NetIn &netIn) override
 Prepares the engine for calculations with initial conditions.
 
BuildDepthType getDepth () const override
 Gets the depth of the network.
 
void rebuild (const fourdst::composition::Composition &comp, BuildDepthType depth) override
 Rebuilds the reaction network based on a new composition.
 
fourdst::composition::Composition collectComposition (fourdst::composition::Composition &comp) const override
 This will return the input comp with the molar abundances of any species not registered in that but registered in the engine active species set to 0.0.
 
template<IsArithmeticOrAD T>
StepDerivatives< T > calculateAllDerivatives (const std::vector< T > &Y_in, const T T9, const T rho, const T Ye, const T mue, const std::function< std::optional< size_t >(const fourdst::atomic::Species &)> speciesLookup, const std::function< bool(const reaction::Reaction &)> &reactionLookup) const
 
- Public Member Functions inherited from gridfire::Engine
virtual ~Engine ()=default
 Virtual destructor.
 

Static Public Member Functions

static std::unordered_map< fourdst::atomic::Species, int > getNetReactionStoichiometry (const reaction::Reaction &reaction)
 Gets the net stoichiometry for a given reaction.
 

Private Types

enum class  JacobianMatrixState { UNINITIALIZED , STALE , READY_DENSE , READY_SPARSE }
 

Private Member Functions

void syncInternalMaps ()
 Synchronizes the internal maps.
 
void collectNetworkSpecies ()
 Collects the unique species in the network.
 
void populateReactionIDMap ()
 Populates the reaction ID map.
 
void populateSpeciesToIndexMap ()
 Populates the species-to-index map.
 
void reserveJacobianMatrix () const
 Reserves space for the Jacobian matrix.
 
void recordADTape () const
 Records the AD tape for the right-hand side of the ODE.
 
void collectAtomicReverseRateAtomicBases ()
 
void precomputeNetwork ()
 
bool validateConservation () const
 Validates mass and charge conservation across all reactions.
 
StepDerivatives< double > calculateAllDerivativesUsingPrecomputation (const fourdst::composition::Composition &comp, const std::vector< double > &bare_rates, const std::vector< double > &bare_reverse_rates, double T9, double rho, const reaction::ReactionSet &activeReactions) const
 
template<IsArithmeticOrAD T>
calculateMolarReactionFlow (const reaction::Reaction &reaction, const std::vector< T > &Y, T T9, T rho, T Ye, T mue, const std::function< std::optional< size_t >(const fourdst::atomic::Species &)> &speciesIDLookup) const
 Calculates the molar reaction flow for a given reaction.
 
template<IsArithmeticOrAD T>
calculateReverseMolarReactionFlow (T T9, T rho, std::vector< T > screeningFactors, const std::vector< T > &Y, size_t reactionIndex, const reaction::Reaction &reaction) const
 
template<IsArithmeticOrAD T>
StepDerivatives< T > calculateAllDerivatives (const std::vector< T > &Y_in, T T9, T rho, T Ye, T mue, std::function< std::optional< size_t >(const fourdst::atomic::Species &)> speciesLookup, const std::function< bool(const reaction::Reaction &)> &reactionLookup) const
 Calculates all derivatives (dY/dt) and the energy generation rate.
 

Private Attributes

std::unordered_map< JacobianMatrixState, std::string > m_jacobianMatrixStateNameMap
 
Config & m_config = Config::getInstance()
 
quill::Logger * m_logger = LogManager::getInstance().getLogger("log")
 
constants m_constants
 
rates::weak::WeakRateInterpolator m_weakRateInterpolator
 Interpolator for weak reaction rates.
 
reaction::ReactionSet m_reactions
 Set of REACLIB reactions in the network.
 
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 performance bottleneck.
 
std::vector< fourdst::atomic::Species > m_networkSpecies
 Vector of unique species in the network.
 
std::unordered_map< std::string_view, fourdst::atomic::Species > m_networkSpeciesMap
 Map from species name to Species object.
 
std::unordered_map< fourdst::atomic::Species, size_t > m_speciesToIndexMap
 Map from species to their index in the stoichiometry matrix.
 
std::unordered_map< size_t, fourdst::atomic::Species > m_indexToSpeciesMap
 Map from index to species in the stoichiometry matrix.
 
boost::numeric::ublas::compressed_matrix< int > m_stoichiometryMatrix
 Stoichiometry matrix (species x reactions).
 
boost::numeric::ublas::compressed_matrix< double > m_jacobianMatrix
 Jacobian matrix (species x species).
 
JacobianMatrixState m_jacobianMatrixState = JacobianMatrixState::UNINITIALIZED
 
CppAD::ADFun< double > m_rhsADFun
 CppAD function for the right-hand side of the ODE.
 
CppAD::ADFun< double > m_epsADFun
 CppAD function for the energy generation rate.
 
CppAD::sparse_jac_work m_jac_work
 Work object for sparse Jacobian calculations.
 
CppAD::sparse_rc< std::vector< size_t > > m_full_jacobian_sparsity_pattern
 Full sparsity pattern for the Jacobian matrix.
 
std::set< std::pair< size_t, size_t > > m_full_sparsity_set
 For quick lookups of the base sparsity pattern.
 
std::vector< std::unique_ptr< AtomicReverseRate > > m_atomicReverseRates
 
screening::ScreeningType m_screeningType = screening::ScreeningType::BARE
 Screening type for the reaction network. Default to no screening.
 
std::unique_ptr< screening::ScreeningModelm_screeningModel = screening::selectScreeningModel(m_screeningType)
 
bool m_usePrecomputation = true
 Flag to enable or disable using precomputed reactions for efficiency. Mathematically, this should not change the results. Generally end users should not need to change this.
 
bool m_useReverseReactions = true
 Flag to enable or disable reverse reactions. If false, only forward reactions are considered.
 
BuildDepthType m_depth
 
std::vector< PrecomputedReactionm_precomputedReactions
 Precomputed reactions for efficiency.
 
std::unordered_map< uint64_t, size_t > m_precomputedReactionIndexMap
 Set of hashed precomputed reactions for quick lookup.
 
std::unique_ptr< partition::PartitionFunctionm_partitionFunction
 Partition function for the network.
 

Detailed Description

A reaction network engine that uses a graph-based representation.

The GraphEngine class implements the DynamicEngine interface using a graph-based representation of the reaction network. It uses sparse matrices for efficient storage and computation of the stoichiometry and Jacobian matrices. Automatic differentiation (AD) is used to calculate the Jacobian matrix.

The engine supports:

  • Calculation of the right-hand side (dY/dt) and energy generation rate.
  • Generation and access to the Jacobian matrix.
  • Generation and access to the stoichiometry matrix.
  • Calculation of molar reaction flows.
  • Access to the set of logical reactions in the network.
  • Computation of timescales for each species.
  • Exporting the network to DOT and CSV formats for visualization and analysis.
See also
engine_abstract.h

Member Enumeration Documentation

◆ JacobianMatrixState

Enumerator
UNINITIALIZED 
STALE 
READY_DENSE 
READY_SPARSE 

Constructor & Destructor Documentation

◆ GraphEngine() [1/4]

gridfire::GraphEngine::GraphEngine ( const fourdst::composition::Composition & composition,
BuildDepthType buildDepth = NetworkBuildDepth::Full )
explicit

Constructs a GraphEngine from a composition.

Parameters
compositionThe composition of the material.

This constructor builds the reaction network from the given composition using the build_reaclib_nuclear_network function.

See also
build_reaclib_nuclear_network

◆ GraphEngine() [2/4]

gridfire::GraphEngine::GraphEngine ( const fourdst::composition::Composition & composition,
const partition::PartitionFunction & partitionFunction,
BuildDepthType buildDepth = NetworkBuildDepth::Full )
explicit

◆ GraphEngine() [3/4]

gridfire::GraphEngine::GraphEngine ( const fourdst::composition::Composition & composition,
const partition::PartitionFunction & partitionFunction,
BuildDepthType buildDepth,
NetworkConstructionFlags reactionTypes )
explicit

◆ GraphEngine() [4/4]

gridfire::GraphEngine::GraphEngine ( const reaction::ReactionSet & reactions)
explicit

Constructs a GraphEngine from a set of reactions.

Parameters
reactionsThe set of reactions to use in the network.

This constructor uses the given set of reactions to construct the reaction network.

Member Function Documentation

◆ calculateAllDerivatives() [1/2]

template<IsArithmeticOrAD T>
StepDerivatives< T > gridfire::GraphEngine::calculateAllDerivatives ( const std::vector< T > & Y_in,
const T T9,
const T rho,
const T Ye,
const T mue,
const std::function< std::optional< size_t >(const fourdst::atomic::Species &)> speciesLookup,
const std::function< bool(const reaction::Reaction &)> & reactionLookup ) const

◆ calculateAllDerivatives() [2/2]

template<IsArithmeticOrAD T>
StepDerivatives< T > gridfire::GraphEngine::calculateAllDerivatives ( const std::vector< T > & Y_in,
T T9,
T rho,
T Ye,
T mue,
std::function< std::optional< size_t >(const fourdst::atomic::Species &)> speciesLookup,
const std::function< bool(const reaction::Reaction &)> & reactionLookup ) const
nodiscardprivate

Calculates all derivatives (dY/dt) and the energy generation rate.

Template Parameters
TThe numeric type to use for the calculation.
Parameters
Y_inVector of molar abundances for all species in the network.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
Ye
mue
speciesLookup
reactionLookup
Returns
StepDerivatives<T> containing dY/dt and energy generation rate.

This method calculates the time derivatives of all species and the specific nuclear energy generation rate for the current state.

◆ calculateAllDerivativesUsingPrecomputation()

StepDerivatives< double > gridfire::GraphEngine::calculateAllDerivativesUsingPrecomputation ( const fourdst::composition::Composition & comp,
const std::vector< double > & bare_rates,
const std::vector< double > & bare_reverse_rates,
double T9,
double rho,
const reaction::ReactionSet & activeReactions ) const
nodiscardprivate

◆ calculateEpsDerivatives() [1/2]

EnergyDerivatives gridfire::GraphEngine::calculateEpsDerivatives ( const fourdst::composition::Composition & comp,
double T9,
double rho ) const
nodiscardoverridevirtual

Calculates the derivatives of the energy generation rate with respect to temperature and density.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
Returns
EnergyDerivatives struct containing the derivatives.

This method computes the partial derivatives of the specific nuclear energy generation rate with respect to temperature (∂ε/∂T) and density (∂ε/∂ρ)

See also
EnergyDerivatives

Implements gridfire::DynamicEngine.

◆ calculateEpsDerivatives() [2/2]

EnergyDerivatives gridfire::GraphEngine::calculateEpsDerivatives ( const fourdst::composition::Composition & comp,
double T9,
double rho,
const reaction::ReactionSet & activeReactions ) const
nodiscard

Calculates the derivatives of the energy generation rate with respect to temperature and density for a subset of reactions.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
activeReactionsThe set of reactions to include in the calculation.
Returns
EnergyDerivatives struct containing the derivatives.

This method computes the partial derivatives of the specific nuclear energy generation rate with respect to temperature (∂ε/∂T) and density (∂ε/∂ρ) considering only the specified subset of reactions. This allows for flexible calculations with different reaction sets without modifying the engine's internal state.

See also
EnergyDerivatives

◆ calculateMolarReactionFlow() [1/2]

double gridfire::GraphEngine::calculateMolarReactionFlow ( const reaction::Reaction & reaction,
const fourdst::composition::Composition & comp,
double T9,
double rho ) const
nodiscardoverridevirtual

Calculates the molar reaction flow for a given reaction.

Parameters
reactionThe reaction for which to calculate the flow.
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
Returns
Molar flow rate for the reaction (e.g., mol/g/s).

This method computes the net rate at which the given reaction proceeds under the current state.

Implements gridfire::DynamicEngine.

◆ calculateMolarReactionFlow() [2/2]

template<IsArithmeticOrAD T>
T gridfire::GraphEngine::calculateMolarReactionFlow ( const reaction::Reaction & reaction,
const std::vector< T > & Y,
T T9,
T rho,
T Ye,
T mue,
const std::function< std::optional< size_t >(const fourdst::atomic::Species &)> & speciesIDLookup ) const
private

Calculates the molar reaction flow for a given reaction.

Template Parameters
TThe numeric type to use for the calculation.
Parameters
reactionThe reaction for which to calculate the flow.
YVector of molar abundances for all species in the network.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
Ye
mue
speciesIDLookup
Returns
Molar flow rate for the reaction (e.g., mol/g/s).

This method computes the net rate at which the given reaction proceeds under the current state.

◆ calculateReverseMolarReactionFlow()

template<IsArithmeticOrAD T>
T gridfire::GraphEngine::calculateReverseMolarReactionFlow ( T T9,
T rho,
std::vector< T > screeningFactors,
const std::vector< T > & Y,
size_t reactionIndex,
const reaction::Reaction & reaction ) const
private

◆ calculateReverseRate()

double gridfire::GraphEngine::calculateReverseRate ( const reaction::Reaction & reaction,
double T9,
double rho,
const fourdst::composition::Composition & comp ) const
nodiscard

Calculates the reverse rate for a given reaction.

Parameters
reactionThe reaction for which to calculate the reverse rate.
T9Temperature in units of 10^9 K.
rho
compComposition object containing current abundances.
Returns
Reverse rate for the reaction (e.g., mol/g/s).

This method computes the reverse rate based on the forward rate and thermodynamic properties of the reaction.

◆ calculateReverseRateTwoBody()

double gridfire::GraphEngine::calculateReverseRateTwoBody ( const reaction::Reaction & reaction,
double T9,
double forwardRate,
double expFactor ) const
nodiscard

Calculates the reverse rate for a two-body reaction.

Parameters
reactionThe reaction for which to calculate the reverse rate.
T9Temperature in units of 10^9 K.
forwardRateThe forward rate of the reaction.
expFactorExponential factor for the reaction.
Returns
Reverse rate for the two-body reaction (e.g., mol/g/s).

This method computes the reverse rate using the forward rate and thermodynamic properties of the reaction.

◆ calculateReverseRateTwoBodyDerivative()

double gridfire::GraphEngine::calculateReverseRateTwoBodyDerivative ( const reaction::Reaction & reaction,
double T9,
double rho,
const fourdst::composition::Composition & comp,
double reverseRate ) const
nodiscard

Calculates the derivative of the reverse rate for a two-body reaction with respect to temperature.

Parameters
reactionThe reaction for which to calculate the derivative.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
compComposition object containing current abundances.
reverseRateThe reverse rate of the reaction.
Returns
Derivative of the reverse rate with respect to temperature.

This method computes the derivative of the reverse rate using automatic differentiation.

◆ calculateRHSAndEnergy() [1/2]

std::expected< StepDerivatives< double >, expectations::StaleEngineError > gridfire::GraphEngine::calculateRHSAndEnergy ( const fourdst::composition::Composition & comp,
double T9,
double rho ) const
nodiscardoverridevirtual

Calculates the right-hand side (dY/dt) and energy generation rate.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
Returns
StepDerivatives<double> containing dY/dt and energy generation rate.

This method calculates the time derivatives of all species and the specific nuclear energy generation rate for the current state.

See also
StepDerivatives

Implements gridfire::Engine.

◆ calculateRHSAndEnergy() [2/2]

std::expected< StepDerivatives< double >, expectations::StaleEngineError > gridfire::GraphEngine::calculateRHSAndEnergy ( const fourdst::composition::Composition & comp,
double T9,
double rho,
const reaction::ReactionSet & activeReactions ) const
nodiscard

Calculates the right-hand side (dY/dt) and energy generation rate for a subset of reactions.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
activeReactionsThe set of reactions to include in the calculation.
Returns
StepDerivatives<double> containing dY/dt and energy generation rate.

This method calculates the time derivatives of all species and the specific nuclear energy generation rate considering only the specified subset of reactions. This allows for flexible calculations with different reaction sets without modifying the engine's internal state.

See also
StepDerivatives

◆ collectAtomicReverseRateAtomicBases()

void gridfire::GraphEngine::collectAtomicReverseRateAtomicBases ( )
private

◆ collectComposition()

fourdst::composition::Composition gridfire::GraphEngine::collectComposition ( fourdst::composition::Composition & comp) const
overridevirtual

This will return the input comp with the molar abundances of any species not registered in that but registered in the engine active species set to 0.0.

Note
Effectively this method does not change input composition; rather it ensures that all species which can be tracked by an instance of GraphEngine are registered in the composition object.
If a species is in the input comp but not in the network
Parameters
compInput Composition
Returns
A new composition where all members of the active species set are registered. And any members not in comp have a molar abundance set to 0.
Exceptions
BadCollectionErrorIf the input composition contains species not present in the network species set

Implements gridfire::DynamicEngine.

◆ collectNetworkSpecies()

void gridfire::GraphEngine::collectNetworkSpecies ( )
private

Collects the unique species in the network.

This method collects the unique species in the network from the reactants and products of all reactions.

◆ exportToCSV()

void gridfire::GraphEngine::exportToCSV ( const std::string & filename) const

Exports the network to a CSV file for analysis.

Parameters
filenameThe name of the CSV file to create.

This method generates a CSV file containing information about the reactions in the network, including the reactants, products, Q-value, and reaction rate coefficients.

Exceptions
std::runtime_errorIf the file cannot be opened for writing.

Example usage:

engine.exportToCSV("network.csv");

◆ exportToDot()

void gridfire::GraphEngine::exportToDot ( const std::string & filename) const

Exports the network to a DOT file for visualization.

Parameters
filenameThe name of the DOT file to create.

This method generates a DOT file that can be used to visualize the reaction network as a graph. The DOT file can be converted to a graphical image using Graphviz.

Exceptions
std::runtime_errorIf the file cannot be opened for writing.

Example usage:

engine.exportToDot("network.dot");

◆ generateJacobianMatrix() [1/3]

void gridfire::GraphEngine::generateJacobianMatrix ( const fourdst::composition::Composition & comp,
double T9,
double rho ) const
overridevirtual

Generates the Jacobian matrix for the current state.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.

This method computes and stores the Jacobian matrix (∂(dY/dt)_i/∂Y_j) for the current state using automatic differentiation. The matrix can then be accessed via getJacobianMatrixEntry().

See also
getJacobianMatrixEntry()

Implements gridfire::DynamicEngine.

◆ generateJacobianMatrix() [2/3]

void gridfire::GraphEngine::generateJacobianMatrix ( const fourdst::composition::Composition & comp,
double T9,
double rho,
const SparsityPattern & sparsityPattern ) const
overridevirtual

Generates the Jacobian matrix for the current state with a specified sparsity pattern.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
sparsityPatternThe sparsity pattern to use for the Jacobian matrix.

This method computes and stores the Jacobian matrix (∂(dY/dt)_i/∂Y_j) for the current state using automatic differentiation, taking into account the provided sparsity pattern. The matrix can then be accessed via getJacobianMatrixEntry().

See also
getJacobianMatrixEntry()

Implements gridfire::DynamicEngine.

◆ generateJacobianMatrix() [3/3]

void gridfire::GraphEngine::generateJacobianMatrix ( const fourdst::composition::Composition & comp,
double T9,
double rho,
const std::vector< fourdst::atomic::Species > & activeSpecies ) const
overridevirtual

Generates the Jacobian matrix for the current state with a specified set of active species. generally this will be much faster than the full matrix generation. Here we use forward mode to generate the Jacobian only for the active species.

Parameters
compThe Composition object containing current abundances.
T9The temperature in units of 10^9 K.
rhoThe density in g/cm^3.
activeSpeciesA vector of Species objects representing the active species.
See also
getJacobianMatrixEntry()
generateJacobianMatrix()

Implements gridfire::DynamicEngine.

◆ generateStoichiometryMatrix()

void gridfire::GraphEngine::generateStoichiometryMatrix ( )
overridevirtual

Generates the stoichiometry matrix for the network.

This method computes and stores the stoichiometry matrix, which encodes the net change of each species in each reaction.

Implements gridfire::DynamicEngine.

◆ getDepth()

BuildDepthType gridfire::GraphEngine::getDepth ( ) const
nodiscardoverridevirtual

Gets the depth of the network.

Returns
The build depth of the network.

This method returns the current build depth of the reaction network, which indicates how many levels of reactions are included in the network.

Reimplemented from gridfire::DynamicEngine.

◆ getJacobianMatrixEntry()

double gridfire::GraphEngine::getJacobianMatrixEntry ( const fourdst::atomic::Species & rowSpecies,
const fourdst::atomic::Species & colSpecies ) const
nodiscardoverridevirtual

Gets an entry from the previously generated Jacobian matrix.

Parameters
rowSpeciesSpecies corresponding to the row index (i).
colSpeciesSpecies corresponding to the column index (j).
Returns
Value of the Jacobian matrix at (i, j).

The Jacobian must have been generated by generateJacobianMatrix() before calling this.

See also
generateJacobianMatrix()

Implements gridfire::DynamicEngine.

◆ getNetReactionStoichiometry()

std::unordered_map< fourdst::atomic::Species, int > gridfire::GraphEngine::getNetReactionStoichiometry ( const reaction::Reaction & reaction)
staticnodiscard

Gets the net stoichiometry for a given reaction.

Parameters
reactionThe reaction for which to get the stoichiometry.
Returns
Map of species to their stoichiometric coefficients.

◆ getNetworkReactions()

const reaction::ReactionSet & gridfire::GraphEngine::getNetworkReactions ( ) const
nodiscardoverridevirtual

Gets the set of logical reactions in the network.

Returns
Reference to the LogicalReactionSet containing all reactions.

Implements gridfire::DynamicEngine.

◆ getNetworkSpecies()

const std::vector< fourdst::atomic::Species > & gridfire::GraphEngine::getNetworkSpecies ( ) const
nodiscardoverridevirtual

Gets the list of species in the network.

Returns
Vector of Species objects representing all network species.

Implements gridfire::Engine.

◆ getPartitionFunction()

const partition::PartitionFunction & gridfire::GraphEngine::getPartitionFunction ( ) const
nodiscard

Gets the partition function used for reaction rate calculations.

Returns
Reference to the PartitionFunction object.

This method provides access to the partition function used in the engine, which is essential for calculating thermodynamic properties and reaction rates.

◆ getScreeningModel()

screening::ScreeningType gridfire::GraphEngine::getScreeningModel ( ) const
nodiscardoverridevirtual

Gets the current electron screening model.

Returns
The currently active screening model type.

Example usage:

screening::ScreeningType currentModel = engine.getScreeningModel();
ScreeningType
Enumerates the available plasma screening models.
Definition screening_types.h:15

Implements gridfire::DynamicEngine.

◆ getSpeciesDestructionTimescales() [1/2]

std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > gridfire::GraphEngine::getSpeciesDestructionTimescales ( const fourdst::composition::Composition & comp,
double T9,
double rho ) const
nodiscardoverridevirtual

Computes destruction timescales for all species in the network.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
Returns
Map from Species to their destruction timescales (s).

This method estimates the destruction timescale for each species, which can be useful for understanding reaction flows and equilibrium states.

Implements gridfire::DynamicEngine.

◆ getSpeciesDestructionTimescales() [2/2]

std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > gridfire::GraphEngine::getSpeciesDestructionTimescales ( const fourdst::composition::Composition & comp,
double T9,
double rho,
const reaction::ReactionSet & activeReactions ) const
nodiscard

Computes destruction timescales for all species in the network considering a subset of reactions.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
activeReactionsThe set of reactions to include in the calculation.
Returns
Map from Species to their destruction timescales (s).

This method estimates the destruction timescale for each species, considering only the specified subset of reactions. This allows for flexible calculations with different reaction sets without modifying the engine's internal state.

◆ getSpeciesIndex()

size_t gridfire::GraphEngine::getSpeciesIndex ( const fourdst::atomic::Species & species) const
nodiscardoverridevirtual

Gets the index of a species in the network.

Parameters
speciesThe species for which to get the index.
Returns
Index of the species in the network, or -1 if not found.

This method returns the index of the given species in the network's species vector. If the species is not found, it returns -1.

Implements gridfire::DynamicEngine.

◆ getSpeciesTimescales() [1/2]

std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > gridfire::GraphEngine::getSpeciesTimescales ( const fourdst::composition::Composition & comp,
double T9,
double rho ) const
nodiscardoverridevirtual

Computes timescales for all species in the network.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
Returns
Map from Species to their characteristic timescales (s).

This method estimates the timescale for abundance change of each species, which can be used for timestep control or diagnostics.

Implements gridfire::DynamicEngine.

◆ getSpeciesTimescales() [2/2]

std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > gridfire::GraphEngine::getSpeciesTimescales ( const fourdst::composition::Composition & comp,
double T9,
double rho,
const reaction::ReactionSet & activeReactions ) const
nodiscard

Computes timescales for all species in the network considering a subset of reactions.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
activeReactionsThe set of reactions to include in the calculation.
Returns
Map from Species to their characteristic timescales (s).

This method estimates the timescale for abundance change of each species, considering only the specified subset of reactions. This allows for flexible calculations with different reaction sets without modifying the engine's internal state.

◆ getStoichiometryMatrixEntry()

int gridfire::GraphEngine::getStoichiometryMatrixEntry ( const fourdst::atomic::Species & species,
const reaction::Reaction & reaction ) const
nodiscardoverridevirtual

Gets an entry from the stoichiometry matrix.

Parameters
speciesSpecies to look up stoichiometry for.
reactionReaction to find.
Returns
Stoichiometric coefficient for the species in the reaction.

The stoichiometry matrix must have been generated by generateStoichiometryMatrix().

See also
generateStoichiometryMatrix()

Implements gridfire::DynamicEngine.

◆ involvesSpecies()

bool gridfire::GraphEngine::involvesSpecies ( const fourdst::atomic::Species & species) const
nodiscard

Checks if a given species is involved in the network.

Parameters
speciesThe species to check.
Returns
True if the species is involved in the network, false otherwise.

◆ isPrecomputationEnabled()

bool gridfire::GraphEngine::isPrecomputationEnabled ( ) const
nodiscard

Checks if precomputation of reaction rates is enabled.

Returns
True if precomputation is enabled, false otherwise.

This method allows checking the current state of precomputation for reaction rates in the engine.

◆ isStale()

bool gridfire::GraphEngine::isStale ( const NetIn & netIn)
overridevirtual

Checks if the engine view is stale and needs to be updated.

Parameters
netInThe current network input (unused).
Returns
True if the view is stale, false otherwise.
Deprecated
This method is deprecated and will be removed in future versions. Stale states are returned as part of the results of methods that require the ability to report them.

Implements gridfire::DynamicEngine.

◆ isUsingReverseReactions()

bool gridfire::GraphEngine::isUsingReverseReactions ( ) const
nodiscard

Checks if reverse reactions are enabled.

Returns
True if reverse reactions are enabled, false otherwise.

This method allows checking whether the engine is configured to use reverse reactions in its calculations.

◆ mapNetInToMolarAbundanceVector()

std::vector< double > gridfire::GraphEngine::mapNetInToMolarAbundanceVector ( const NetIn & netIn) const
nodiscardoverridevirtual

Maps the NetIn object to a vector of molar abundances.

Parameters
netInThe NetIn object containing the input conditions.
Returns
Vector of molar abundances corresponding to the species in the network.

This method converts the NetIn object into a vector of molar abundances for each species in the network, which can be used for further calculations.

Implements gridfire::DynamicEngine.

◆ populateReactionIDMap()

void gridfire::GraphEngine::populateReactionIDMap ( )
private

Populates the reaction ID map.

This method populates the reaction ID map, which maps reaction IDs to REACLIBReaction objects.

◆ populateSpeciesToIndexMap()

void gridfire::GraphEngine::populateSpeciesToIndexMap ( )
private

Populates the species-to-index map.

This method populates the species-to-index map, which maps species to their index in the stoichiometry matrix.

◆ precomputeNetwork()

void gridfire::GraphEngine::precomputeNetwork ( )
private

◆ primeEngine()

PrimingReport gridfire::GraphEngine::primeEngine ( const NetIn & netIn)
nodiscardoverridevirtual

Prepares the engine for calculations with initial conditions.

Parameters
netInThe input conditions for the network.
Returns
PrimingReport containing information about the priming process.

This method initializes the engine with the provided input conditions, setting up reactions, species, and precomputing necessary data.

Implements gridfire::DynamicEngine.

◆ rebuild()

void gridfire::GraphEngine::rebuild ( const fourdst::composition::Composition & comp,
BuildDepthType depth )
overridevirtual

Rebuilds the reaction network based on a new composition.

Parameters
compThe new composition to use for rebuilding the network.
depthThe build depth to use for the network.

This method rebuilds the reaction network using the provided composition and build depth. It updates all internal data structures accordingly.

Reimplemented from gridfire::DynamicEngine.

◆ recordADTape()

void gridfire::GraphEngine::recordADTape ( ) const
private

Records the AD tape for the right-hand side of the ODE.

This method records the AD tape for the right-hand side of the ODE, which is used to calculate the Jacobian matrix using automatic differentiation.

Exceptions
std::runtime_errorIf there are no species in the network.

◆ reserveJacobianMatrix()

void gridfire::GraphEngine::reserveJacobianMatrix ( ) const
private

Reserves space for the Jacobian matrix.

This method reserves space for the Jacobian matrix, which is used to store the partial derivatives of the right-hand side of the ODE with respect to the species abundances.

◆ setNetworkReactions()

void gridfire::GraphEngine::setNetworkReactions ( const reaction::ReactionSet & reactions)
overridevirtual

Sets the reactions for the network.

Parameters
reactionsThe set of reactions to use in the network.

This method replaces the current set of reactions in the network with the provided set. It marks the engine as stale, requiring regeneration of matrices and recalculation of rates.

Implements gridfire::DynamicEngine.

◆ setPrecomputation()

void gridfire::GraphEngine::setPrecomputation ( bool precompute)

Sets whether to precompute reaction rates.

Parameters
precomputeTrue to enable precomputation, false to disable.

This method allows enabling or disabling precomputation of reaction rates for performance optimization. When enabled, reaction rates are computed once and stored for later use.

Postcondition
If precomputation is enabled, reaction rates will be precomputed and cached. If disabled, reaction rates will be computed on-the-fly as needed.

◆ setScreeningModel()

void gridfire::GraphEngine::setScreeningModel ( screening::ScreeningType model)
overridevirtual

Sets the electron screening model for reaction rate calculations.

Parameters
modelThe type of screening model to use.

This method allows changing the screening model at runtime. Screening corrections account for the electrostatic shielding of nuclei by electrons, which affects reaction rates in dense stellar plasmas.

Implements gridfire::DynamicEngine.

◆ setUseReverseReactions()

void gridfire::GraphEngine::setUseReverseReactions ( bool useReverse)

Sets whether to use reverse reactions in the engine.

Parameters
useReverseTrue to enable reverse reactions, false to disable.

This method allows enabling or disabling reverse reactions in the engine. If disabled, only forward reactions will be considered in calculations.

Postcondition
If reverse reactions are enabled, the engine will consider both forward and reverse reactions in its calculations. If disabled, only forward reactions will be considered.

◆ syncInternalMaps()

void gridfire::GraphEngine::syncInternalMaps ( )
private

Synchronizes the internal maps.

This method synchronizes the internal maps used by the engine, including the species map, reaction ID map, and species-to-index map. It also generates the stoichiometry matrix and records the AD tape.

◆ update()

fourdst::composition::Composition gridfire::GraphEngine::update ( const NetIn & netIn)
overridevirtual

Updates the state of the network and the composition to be usable for the current network.

For graph engine all this does is ensure that the returned composition has all the species in the network registered. if a species was already in the composition is will keep its abundance, otherwise it will be added with zero abundance.

Parameters
netInThe input netIn to use, this includes the composition, temperature, and density
Returns
The updated composition that includes all species in the network.

Implements gridfire::DynamicEngine.

◆ validateConservation()

bool gridfire::GraphEngine::validateConservation ( ) const
nodiscardprivate

Validates mass and charge conservation across all reactions.

Returns
True if all reactions conserve mass and charge, false otherwise.

This method checks that all reactions in the network conserve mass and charge. If any reaction does not conserve mass or charge, an error message is logged and false is returned.

Member Data Documentation

◆ m_atomicReverseRates

std::vector<std::unique_ptr<AtomicReverseRate> > gridfire::GraphEngine::m_atomicReverseRates
private

◆ m_config

Config& gridfire::GraphEngine::m_config = Config::getInstance()
private

◆ m_constants

constants gridfire::GraphEngine::m_constants
private

◆ m_depth

BuildDepthType gridfire::GraphEngine::m_depth
private

◆ m_epsADFun

CppAD::ADFun<double> gridfire::GraphEngine::m_epsADFun
mutableprivate

CppAD function for the energy generation rate.

◆ m_full_jacobian_sparsity_pattern

CppAD::sparse_rc<std::vector<size_t> > gridfire::GraphEngine::m_full_jacobian_sparsity_pattern
private

Full sparsity pattern for the Jacobian matrix.

◆ m_full_sparsity_set

std::set<std::pair<size_t, size_t> > gridfire::GraphEngine::m_full_sparsity_set
private

For quick lookups of the base sparsity pattern.

◆ m_indexToSpeciesMap

std::unordered_map<size_t, fourdst::atomic::Species> gridfire::GraphEngine::m_indexToSpeciesMap
private

Map from index to species in the stoichiometry matrix.

◆ m_jac_work

CppAD::sparse_jac_work gridfire::GraphEngine::m_jac_work
mutableprivate

Work object for sparse Jacobian calculations.

◆ m_jacobianMatrix

boost::numeric::ublas::compressed_matrix<double> gridfire::GraphEngine::m_jacobianMatrix
mutableprivate

Jacobian matrix (species x species).

◆ m_jacobianMatrixState

JacobianMatrixState gridfire::GraphEngine::m_jacobianMatrixState = JacobianMatrixState::UNINITIALIZED
mutableprivate

◆ m_jacobianMatrixStateNameMap

std::unordered_map<JacobianMatrixState, std::string> gridfire::GraphEngine::m_jacobianMatrixStateNameMap
private
Initial value:
= {
{JacobianMatrixState::READY_SPARSE, "Ready (sparse)"},
}
@ READY_DENSE
Definition engine_graph.h:782
@ STALE
Definition engine_graph.h:781
@ READY_SPARSE
Definition engine_graph.h:783
@ UNINITIALIZED
Definition engine_graph.h:780

◆ m_logger

quill::Logger* gridfire::GraphEngine::m_logger = LogManager::getInstance().getLogger("log")
private

◆ m_networkSpecies

std::vector<fourdst::atomic::Species> gridfire::GraphEngine::m_networkSpecies
private

Vector of unique species in the network.

◆ m_networkSpeciesMap

std::unordered_map<std::string_view, fourdst::atomic::Species> gridfire::GraphEngine::m_networkSpeciesMap
private

Map from species name to Species object.

◆ m_partitionFunction

std::unique_ptr<partition::PartitionFunction> gridfire::GraphEngine::m_partitionFunction
private

Partition function for the network.

◆ m_precomputedReactionIndexMap

std::unordered_map<uint64_t, size_t> gridfire::GraphEngine::m_precomputedReactionIndexMap
private

Set of hashed precomputed reactions for quick lookup.

◆ m_precomputedReactions

std::vector<PrecomputedReaction> gridfire::GraphEngine::m_precomputedReactions
private

Precomputed reactions for efficiency.

◆ m_reactionIDMap

std::unordered_map<std::string_view, reaction::Reaction*> gridfire::GraphEngine::m_reactionIDMap
private

Map from reaction ID to REACLIBReaction. //PERF: This makes copies of REACLIBReaction and could be a performance bottleneck.

◆ m_reactions

reaction::ReactionSet gridfire::GraphEngine::m_reactions
private

Set of REACLIB reactions in the network.

◆ m_rhsADFun

CppAD::ADFun<double> gridfire::GraphEngine::m_rhsADFun
mutableprivate

CppAD function for the right-hand side of the ODE.

◆ m_screeningModel

std::unique_ptr<screening::ScreeningModel> gridfire::GraphEngine::m_screeningModel = screening::selectScreeningModel(m_screeningType)
private

◆ m_screeningType

screening::ScreeningType gridfire::GraphEngine::m_screeningType = screening::ScreeningType::BARE
private

Screening type for the reaction network. Default to no screening.

◆ m_speciesToIndexMap

std::unordered_map<fourdst::atomic::Species, size_t> gridfire::GraphEngine::m_speciesToIndexMap
private

Map from species to their index in the stoichiometry matrix.

◆ m_stoichiometryMatrix

boost::numeric::ublas::compressed_matrix<int> gridfire::GraphEngine::m_stoichiometryMatrix
private

Stoichiometry matrix (species x reactions).

◆ m_usePrecomputation

bool gridfire::GraphEngine::m_usePrecomputation = true
private

Flag to enable or disable using precomputed reactions for efficiency. Mathematically, this should not change the results. Generally end users should not need to change this.

◆ m_useReverseReactions

bool gridfire::GraphEngine::m_useReverseReactions = true
private

Flag to enable or disable reverse reactions. If false, only forward reactions are considered.

◆ m_weakRateInterpolator

rates::weak::WeakRateInterpolator gridfire::GraphEngine::m_weakRateInterpolator
private

Interpolator for weak reaction rates.


The documentation for this class was generated from the following files: