|
GridFire v0.7.6rc4.0
General Purpose Nuclear Network
|
A reaction network engine that uses a graph-based representation. More...
#include <engine_graph.h>
Classes | |
| class | AtomicReverseRate |
| struct | constants |
| struct | PrecomputationKernelResults |
| 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. | |
| void | addReaction (const reaction::Reaction &reaction) |
| void | addReaction (const std::string &reaction_id) |
| std::unique_ptr< scratch::StateBlob > | constructStateBlob (const scratch::StateBlob *blob=nullptr) const override |
| std::expected< StepDerivatives< double >, engine::EngineStatus > | calculateRHSAndEnergy (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &comp, double T9, double rho, bool trust) const override |
| Calculates the right-hand side (dY/dt) and energy generation rate. | |
| std::expected< StepDerivatives< double >, EngineStatus > | calculateRHSAndEnergy (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &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 (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const override |
| Calculates the derivatives of the energy generation rate with respect to temperature and density. | |
| EnergyDerivatives | calculateEpsDerivatives (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &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 | generate_jacobian_sparsity_pattern () |
| NetworkJacobian | generateJacobianMatrix (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const override |
| Generates the Jacobian matrix for the current state. | |
| NetworkJacobian | generateJacobianMatrix (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &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. | |
| NetworkJacobian | generateJacobianMatrix (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &comp, double T9, double rho, const SparsityPattern &sparsityPattern) const override |
| Generates the Jacobian matrix for the current state with a specified sparsity pattern. | |
| double | calculateMolarReactionFlow (scratch::StateBlob &, const reaction::Reaction &reaction, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const override |
| Calculates the molar reaction flow for a given reaction. | |
| const std::vector< fourdst::atomic::Species > & | getNetworkSpecies (scratch::StateBlob &ctx) const override |
| Gets the list of species in the network. | |
| const reaction::ReactionSet & | getNetworkReactions (scratch::StateBlob &) const override |
| Gets the set of logical reactions in the network. | |
| std::expected< std::unordered_map< fourdst::atomic::Species, double >, EngineStatus > | getSpeciesTimescales (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const override |
| Computes timescales for all species in the network. | |
| std::expected< std::unordered_map< fourdst::atomic::Species, double >, EngineStatus > | getSpeciesTimescales (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &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 >, EngineStatus > | getSpeciesDestructionTimescales (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &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 >, EngineStatus > | getSpeciesDestructionTimescales (scratch::StateBlob &ctx, const fourdst::composition::CompositionAbstract &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 | project (scratch::StateBlob &ctx, const NetIn &netIn) const override |
| Updates the state of the network and the composition to be usable for the current network. | |
| bool | involvesSpecies (scratch::StateBlob &ctx, const fourdst::atomic::Species &species) const |
| Checks if a given species is involved in the network. | |
| void | exportToDot (scratch::StateBlob &ctx, const std::string &filename) const |
| Exports the network to a DOT file for visualization. | |
| void | exportToCSV (scratch::StateBlob &ctx, const std::string &filename) const |
| Exports the network to a CSV file for analysis. | |
| screening::ScreeningType | getScreeningModel (scratch::StateBlob &ctx) const override |
| Gets the current electron screening model. | |
| bool | isPrecomputationEnabled (scratch::StateBlob &ctx) const |
| Checks if precomputation of reaction rates is enabled. | |
| const partition::PartitionFunction & | getPartitionFunction (scratch::StateBlob &ctx) const |
| Gets the partition function used for reaction rate calculations. | |
| double | calculateReverseRate (const reaction::Reaction &reaction, double T9, double rho, const fourdst::composition::CompositionAbstract &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 (scratch::StateBlob &ctx) const |
| Checks if reverse reactions are enabled. | |
| size_t | getSpeciesIndex (scratch::StateBlob &ctx, const fourdst::atomic::Species &species) const override |
| Gets the index of a species in the network. | |
| PrimingReport | primeEngine (scratch::StateBlob &ctx, const NetIn &netIn) const override |
| Prepares the engine for calculations with initial conditions. | |
| fourdst::composition::Composition | collectComposition (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) 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. | |
| SpeciesStatus | getSpeciesStatus (scratch::StateBlob &, const fourdst::atomic::Species &species) const override |
| Gets the status of a species in the network. | |
| bool | get_store_intermediate_reaction_contributions () const |
| void | set_store_intermediate_reaction_contributions (const bool value) |
| std::optional< StepDerivatives< double > > | getMostRecentRHSCalculation (scratch::StateBlob &) const override |
| const CppAD::ADFun< double > & | getAuthoritativeADFun () const |
| 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::DynamicEngine | |
| virtual reaction::ReactionSet | getInactiveNetworkReactions (scratch::StateBlob &ctx) const |
| Get the set of inactive reactions in the network. | |
| virtual double | getInactiveReactionMolarReactionFlow (scratch::StateBlob &ctx, const reaction::Reaction &reaction, const fourdst::composition::CompositionAbstract &comp, const double T9, const double rho) const |
Public Member Functions inherited from gridfire::engine::Engine | |
| virtual | ~Engine ()=default |
| Virtual destructor. | |
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 | recordADTape () |
| Records the AD tape for the right-hand side of the ODE. | |
| void | collectAtomicReverseRateAtomicBases () |
| void | precomputeNetwork () |
| double | compute_reaction_flow (scratch::StateBlob &ctx, const std::vector< double > &local_abundances, const std::vector< double > &screening_factors, const std::vector< double > &bare_rates, const std::vector< double > &bare_reverse_rates, double rho, size_t reactionCounter, const reaction::Reaction &reaction, size_t reactionIndex, const PrecomputedReaction &precomputedReaction) const |
| std::pair< double, double > | compute_neutrino_fluxes (scratch::StateBlob &ctx, double netFlow, const reaction::Reaction &reaction) const |
| PrecomputationKernelResults | accumulate_flows_serial (scratch::StateBlob &ctx, const std::vector< double > &local_abundances, const std::vector< double > &screening_factors, const std::vector< double > &bare_rates, const std::vector< double > &bare_reverse_rates, double rho, const reaction::ReactionSet &activeReactions) const |
| StepDerivatives< double > | calculateAllDerivativesUsingPrecomputation (scratch::StateBlob &ctx, const fourdst::composition::CompositionAbstract &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> | |
| 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> | |
| 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< config::GridFireConfig > | m_config |
| 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. | |
| std::unique_ptr< partition::PartitionFunction > | m_partitionFunction |
| Partition function for the network. | |
| 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::ScreeningModel > | m_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. | |
| std::vector< PrecomputedReaction > | m_precomputed_reactions |
| std::unordered_map< uint64_t, size_t > | m_precomputed_reaction_index_map |
| bool | m_useReverseReactions = false |
| Flag to enable or disable reverse reactions. If false, only forward reactions are considered. | |
| bool | m_store_intermediate_reaction_contributions = false |
| Flag to enable or disable storing intermediate reaction contributions for debugging. | |
| BuildDepthType | m_depth |
| CppAD::ADFun< double > | m_authoritativeADFun |
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:
|
strongprivate |
|
explicit |
Constructs a GraphEngine from a composition.
| composition | The composition of the material. |
This constructor builds the reaction network from the given composition using the build_reaclib_nuclear_network function.
|
explicit |
|
explicit |
|
explicit |
Constructs a GraphEngine from a set of reactions.
| reactions | The set of reactions to use in the network. |
This constructor uses the given set of reactions to construct the reaction network.
|
private |
| void gridfire::engine::GraphEngine::addReaction | ( | const reaction::Reaction & | reaction | ) |
| void gridfire::engine::GraphEngine::addReaction | ( | const std::string & | reaction_id | ) |
| StepDerivatives< T > gridfire::engine::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 |
|
nodiscardprivate |
Calculates all derivatives (dY/dt) and the energy generation rate.
| T | The numeric type to use for the calculation. |
| Y_in | Vector of molar abundances for all species in the network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
| Ye | |
| mue | |
| speciesLookup | |
| reactionLookup |
This method calculates the time derivatives of all species and the specific nuclear energy generation rate for the current state.
|
nodiscardprivate |
|
nodiscardoverridevirtual |
Calculates the derivatives of the energy generation rate with respect to temperature and density.
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
This method computes the partial derivatives of the specific nuclear energy generation rate with respect to temperature (∂ε/∂T) and density (∂ε/∂ρ)
Implements gridfire::engine::DynamicEngine.
|
nodiscard |
Calculates the derivatives of the energy generation rate with respect to temperature and density for a subset of reactions.
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
| activeReactions | The set of reactions to include in the calculation. |
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.
|
private |
Calculates the molar reaction flow for a given reaction.
| T | The numeric type to use for the calculation. |
| reaction | The reaction for which to calculate the flow. |
| Y | Vector of molar abundances for all species in the network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
| Ye | |
| mue | |
| speciesIDLookup |
This method computes the net rate at which the given reaction proceeds under the current state.
|
nodiscardoverridevirtual |
Calculates the molar reaction flow for a given reaction.
| reaction | The reaction for which to calculate the flow. |
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
This method computes the net rate at which the given reaction proceeds under the current state.
Implements gridfire::engine::DynamicEngine.
|
private |
|
nodiscard |
Calculates the reverse rate for a given reaction.
| reaction | The reaction for which to calculate the reverse rate. |
| T9 | Temperature in units of 10^9 K. |
| rho | |
| comp | Composition object containing current abundances. |
This method computes the reverse rate based on the forward rate and thermodynamic properties of the reaction.
|
nodiscard |
Calculates the reverse rate for a two-body reaction.
| reaction | The reaction for which to calculate the reverse rate. |
| T9 | Temperature in units of 10^9 K. |
| forwardRate | The forward rate of the reaction. |
| expFactor | Exponential factor for the reaction. |
This method computes the reverse rate using the forward rate and thermodynamic properties of the reaction.
|
nodiscard |
Calculates the derivative of the reverse rate for a two-body reaction with respect to temperature.
| reaction | The reaction for which to calculate the derivative. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
| comp | Composition object containing current abundances. |
| reverseRate | The reverse rate of the reaction. |
This method computes the derivative of the reverse rate using automatic differentiation.
|
nodiscardoverridevirtual |
Calculates the right-hand side (dY/dt) and energy generation rate.
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
| trust |
This method calculates the time derivatives of all species and the specific nuclear energy generation rate for the current state.
Implements gridfire::engine::Engine.
|
nodiscard |
Calculates the right-hand side (dY/dt) and energy generation rate for a subset of reactions.
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
| activeReactions | The set of reactions to include in the calculation. |
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.
|
private |
|
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.
| comp | Input Composition |
| T9 | |
| rho | |
| T9 | |
| rho |
| BadCollectionError | If the input composition contains species not present in the network species set |
Implements gridfire::engine::DynamicEngine.
|
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.
|
private |
|
private |
|
overridevirtual |
Implements gridfire::engine::DynamicEngine.
| void gridfire::engine::GraphEngine::exportToCSV | ( | scratch::StateBlob & | ctx, |
| const std::string & | filename ) const |
Exports the network to a CSV file for analysis.
| filename | The 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.
| std::runtime_error | If the file cannot be opened for writing. |
Example usage:
| void gridfire::engine::GraphEngine::exportToDot | ( | scratch::StateBlob & | ctx, |
| const std::string & | filename ) const |
Exports the network to a DOT file for visualization.
| filename | The 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.
| std::runtime_error | If the file cannot be opened for writing. |
Example usage:
| void gridfire::engine::GraphEngine::generate_jacobian_sparsity_pattern | ( | ) |
|
nodiscardoverridevirtual |
Generates the Jacobian matrix for the current state.
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density 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().
Implements gridfire::engine::DynamicEngine.
|
nodiscardoverridevirtual |
Generates the Jacobian matrix for the current state with a specified sparsity pattern.
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
| sparsityPattern | The 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().
Implements gridfire::engine::DynamicEngine.
|
nodiscardoverridevirtual |
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.
| comp | The Composition object containing current abundances. |
| T9 | The temperature in units of 10^9 K. |
| rho | The density in g/cm^3. |
| activeSpecies | A vector of Species objects representing the active species. |
Implements gridfire::engine::DynamicEngine.
|
inlinenodiscard |
|
inlinenodiscard |
|
nodiscardoverridevirtual |
Implements gridfire::engine::DynamicEngine.
|
nodiscardoverridevirtual |
Gets the set of logical reactions in the network.
Implements gridfire::engine::DynamicEngine.
|
nodiscardoverridevirtual |
Gets the list of species in the network.
Implements gridfire::engine::Engine.
|
nodiscard |
Gets the partition function used for reaction rate calculations.
This method provides access to the partition function used in the engine, which is essential for calculating thermodynamic properties and reaction rates.
|
nodiscardoverridevirtual |
Gets the current electron screening model.
Example usage:
Implements gridfire::engine::DynamicEngine.
|
nodiscardoverridevirtual |
Computes destruction timescales for all species in the network.
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
This method estimates the destruction timescale for each species, which can be useful for understanding reaction flows and equilibrium states.
Implements gridfire::engine::DynamicEngine.
|
nodiscard |
Computes destruction timescales for all species in the network considering a subset of reactions.
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
| activeReactions | The set of reactions to include in the calculation. |
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.
|
nodiscardoverridevirtual |
Gets the index of a species in the network.
| species | The species for which to get the index. |
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::engine::DynamicEngine.
|
nodiscardoverridevirtual |
Gets the status of a species in the network.
| species | The species for which to get the status. |
This method checks if the given species is part of the network and returns its status (e.g., Active, Inactive, NotFound).
Implements gridfire::engine::DynamicEngine.
|
nodiscardoverridevirtual |
Computes timescales for all species in the network.
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
This method estimates the timescale for abundance change of each species, which can be used for timestep control or diagnostics.
Implements gridfire::engine::DynamicEngine.
|
nodiscard |
Computes timescales for all species in the network considering a subset of reactions.
| comp | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
| activeReactions | The set of reactions to include in the calculation. |
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.
|
nodiscard |
Checks if a given species is involved in the network.
| species | The species to check. |
|
nodiscard |
Checks if precomputation of reaction rates is enabled.
This method allows checking the current state of precomputation for reaction rates in the engine.
|
nodiscard |
Checks if reverse reactions are enabled.
This method allows checking whether the engine is configured to use reverse reactions in its calculations.
|
private |
Populates the reaction ID map.
This method populates the reaction ID map, which maps reaction IDs to REACLIBReaction objects.
|
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.
|
private |
|
nodiscardoverridevirtual |
Prepares the engine for calculations with initial conditions.
| netIn | The input conditions for the network. |
This method initializes the engine with the provided input conditions, setting up reactions, species, and precomputing necessary data.
Implements gridfire::engine::DynamicEngine.
|
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.
| netIn | The input netIn to use, this includes the composition, temperature, and density |
Implements gridfire::engine::DynamicEngine.
|
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.
| std::runtime_error | If there are no species in the network. |
|
inline |
|
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.
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
Full sparsity pattern for the Jacobian matrix.
|
private |
For quick lookups of the base sparsity pattern.
|
private |
Map from index to species in the stoichiometry matrix.
|
private |
|
private |
|
private |
Vector of unique species in the network.
|
private |
Map from species name to Species object.
|
private |
Partition function for the network.
|
private |
|
private |
|
private |
Map from reaction ID to REACLIBReaction. //PERF: This makes copies of REACLIBReaction and could be a performance bottleneck.
|
private |
Set of REACLIB reactions in the network.
|
private |
|
private |
Screening type for the reaction network. Default to no screening.
|
private |
Map from species to their index in the stoichiometry matrix.
|
private |
Flag to enable or disable storing intermediate reaction contributions for debugging.
|
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.
|
private |
Flag to enable or disable reverse reactions. If false, only forward reactions are considered.
|
private |
Interpolator for weak reaction rates.