|
GridFire 0.0.1a
General Purpose Nuclear Network
|
An engine view that partitions the reaction network into multiple groups based on timescales. More...
#include <engine_multiscale.h>
Classes | |
| struct | CacheStats |
| Struct for tracking cache statistics. More... | |
| struct | EigenFunctor |
| Functor for solving QSE abundances using Eigen's nonlinear optimization. More... | |
| struct | QSEGroup |
| Struct representing a QSE group. More... | |
Public Member Functions | |
| MultiscalePartitioningEngineView (GraphEngine &baseEngine) | |
| Constructs a MultiscalePartitioningEngineView. | |
| const std::vector< fourdst::atomic::Species > & | getNetworkSpecies () const override |
| Gets the list of species in the network. | |
| std::expected< StepDerivatives< double >, expectations::StaleEngineError > | calculateRHSAndEnergy (const std::vector< double > &Y_full, double T9, double rho) const override |
| Calculates the right-hand side (dY/dt) and energy generation. | |
| void | generateJacobianMatrix (const std::vector< double > &Y_full, double T9, double rho) const override |
| Generates the Jacobian matrix for the current state. | |
| double | getJacobianMatrixEntry (int i_full, int j_full) const override |
| Gets an entry from the previously generated Jacobian matrix. | |
| void | generateStoichiometryMatrix () override |
| Generates the stoichiometry matrix for the network. | |
| int | getStoichiometryMatrixEntry (int speciesIndex, int reactionIndex) const override |
| Gets an entry from the stoichiometry matrix. | |
| double | calculateMolarReactionFlow (const reaction::Reaction &reaction, const std::vector< double > &Y_full, double T9, double rho) const override |
| Calculates the molar reaction flow for a given reaction. | |
| const reaction::LogicalReactionSet & | getNetworkReactions () const override |
| Gets the set of logical reactions in the network. | |
| void | setNetworkReactions (const reaction::LogicalReactionSet &reactions) override |
| Sets the set of logical reactions in the network. | |
| std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > | getSpeciesTimescales (const std::vector< double > &Y, double T9, double rho) const override |
| Computes timescales for all species in the network. | |
| std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > | getSpeciesDestructionTimescales (const std::vector< double > &Y, double T9, double rho) const override |
| Computes destruction timescales for all species in the network. | |
| fourdst::composition::Composition | update (const NetIn &netIn) override |
| Updates the internal state of the engine, performing partitioning and QSE equilibration. | |
| bool | isStale (const NetIn &netIn) override |
| Checks if the engine's internal state is stale relative to the provided conditions. | |
| void | setScreeningModel (screening::ScreeningType model) override |
| Sets the electron screening model. | |
| screening::ScreeningType | getScreeningModel () const override |
| Gets the current electron screening model. | |
| const DynamicEngine & | getBaseEngine () const override |
| Gets the base engine. | |
| std::vector< std::vector< size_t > > | analyzeTimescalePoolConnectivity (const std::vector< std::vector< size_t > > ×cale_pools, const std::vector< double > &Y, double T9, double rho) const |
| Analyzes the connectivity of timescale pools. | |
| void | partitionNetwork (const std::vector< double > &Y, double T9, double rho) |
| Partitions the network into dynamic and algebraic (QSE) groups based on timescales. | |
| void | partitionNetwork (const NetIn &netIn) |
Partitions the network based on timescales from a NetIn struct. | |
| void | exportToDot (const std::string &filename, const std::vector< double > &Y, const double T9, const double rho) const |
| Exports the network to a DOT file for visualization. | |
| int | getSpeciesIndex (const fourdst::atomic::Species &species) const override |
| Gets the index of a species in the full network. | |
| std::vector< double > | mapNetInToMolarAbundanceVector (const NetIn &netIn) const override |
Maps a NetIn struct to a molar abundance vector for the full network. | |
| PrimingReport | primeEngine (const NetIn &netIn) override |
| Primes the engine with a specific species. | |
| std::vector< fourdst::atomic::Species > | getFastSpecies () const |
| Gets the fast species in the network. | |
| const std::vector< fourdst::atomic::Species > & | getDynamicSpecies () const |
| Gets the dynamic species in the network. | |
| fourdst::composition::Composition | equilibrateNetwork (const std::vector< double > &Y, double T9, double rho) |
| Equilibrates the network by partitioning and solving for QSE abundances. | |
| fourdst::composition::Composition | equilibrateNetwork (const NetIn &netIn) |
Equilibrates the network using QSE from a NetIn struct. | |
Public Member Functions inherited from gridfire::DynamicEngine | |
| virtual void | generateJacobianMatrix (const std::vector< double > &Y_dynamic, double T9, double rho, const SparsityPattern &sparsityPattern) const |
| virtual BuildDepthType | getDepth () const |
| virtual void | rebuild (const fourdst::composition::Composition &comp, BuildDepthType depth) |
Public Member Functions inherited from gridfire::Engine | |
| virtual | ~Engine ()=default |
| Virtual destructor. | |
Public Member Functions inherited from gridfire::EngineView< DynamicEngine > | |
| virtual | ~EngineView ()=default |
| Virtual destructor. | |
Private Types | |
| typedef std::tuple< std::vector< fourdst::atomic::Species >, std::vector< size_t >, std::vector< fourdst::atomic::Species >, std::vector< size_t > > | QSEPartition |
| Type alias for a QSE partition. | |
Private Member Functions | |
| std::vector< std::vector< size_t > > | partitionByTimescale (const std::vector< double > &Y_full, double T9, double rho) const |
| Partitions the network by timescale. | |
| std::unordered_map< size_t, std::vector< size_t > > | buildConnectivityGraph (const std::unordered_set< size_t > &fast_reaction_indices) const |
| Builds a connectivity graph from a set of fast reaction indices. | |
| std::vector< QSEGroup > | validateGroupsWithFluxAnalysis (const std::vector< QSEGroup > &candidate_groups, const std::vector< double > &Y, double T9, double rho) const |
| Validates candidate QSE groups using flux analysis. | |
| std::vector< double > | solveQSEAbundances (const std::vector< double > &Y_full, double T9, double rho) |
| Solves for the QSE abundances of the algebraic species in a given state. | |
| size_t | identifyMeanSlowestPool (const std::vector< std::vector< size_t > > &pools, const std::vector< double > &Y, double T9, double rho) const |
| Identifies the pool with the slowest mean timescale. | |
| std::unordered_map< size_t, std::vector< size_t > > | buildConnectivityGraph (const std::vector< size_t > &species_pool) const |
| Builds a connectivity graph from a species pool. | |
| std::vector< QSEGroup > | constructCandidateGroups (const std::vector< std::vector< size_t > > &candidate_pools, const std::vector< double > &Y, double T9, double rho) const |
| Constructs candidate QSE groups from connected timescale pools. | |
Private Attributes | |
| quill::Logger * | m_logger = LogManager::getInstance().getLogger("log") |
| Logger instance for logging messages. | |
| GraphEngine & | m_baseEngine |
| The base engine to which this view delegates calculations. | |
| std::vector< QSEGroup > | m_qse_groups |
| The list of identified equilibrium groups. | |
| std::vector< fourdst::atomic::Species > | m_dynamic_species |
| The simplified set of species presented to the solver (the "slow" species). | |
| std::vector< size_t > | m_dynamic_species_indices |
| Indices mapping the dynamic species back to the base engine's full species list. | |
| std::vector< fourdst::atomic::Species > | m_algebraic_species |
| Species that are treated as algebraic (in QSE) in the QSE groups. | |
| std::vector< size_t > | m_algebraic_species_indices |
| Indices of algebraic species in the full network. | |
| std::vector< size_t > | m_activeSpeciesIndices |
| Indices of all species considered active in the current partition (dynamic + algebraic). | |
| std::vector< size_t > | m_activeReactionIndices |
| Indices of all reactions involving only active species. | |
| std::unordered_map< QSECacheKey, std::vector< double > > | m_qse_abundance_cache |
| Cache for QSE abundances based on T9, rho, and Y. | |
| CacheStats | m_cacheStats |
| Statistics for the QSE abundance cache. | |
An engine view that partitions the reaction network into multiple groups based on timescales.
@purpose This class is designed to accelerate the integration of stiff nuclear reaction networks. It identifies species that react on very short timescales ("fast" species) and treats them as being in Quasi-Steady-State Equilibrium (QSE). Their abundances are solved for algebraically, removing their stiff differential equations from the system. The remaining "slow" or "dynamic" species are integrated normally. This significantly improves the stability and performance of the solver.
@how The core logic resides in the partitionNetwork() and equilibrateNetwork() methods. The partitioning process involves:
getSpeciesDestructionTimescales from the base engine, all species are sorted by their characteristic timescales.solveQSEAbundances uses a Levenberg-Marquardt nonlinear solver (Eigen::LevenbergMarquardt) to find the equilibrium abundances of the "algebraic" species, holding the "seed" species constant.All calculations are cached using QSECacheKey to avoid re-partitioning and re-solving for similar thermodynamic conditions.
<DynamicEngine>
|
private |
Type alias for a QSE partition.
A QSE partition is a tuple containing the fast species, their indices, the slow species, and their indices.
|
explicit |
Constructs a MultiscalePartitioningEngineView.
| baseEngine | The underlying GraphEngine to which this view delegates calculations. It must be a GraphEngine and not a more general DynamicEngine because this view relies on its specific implementation details. |
| std::vector< std::vector< size_t > > gridfire::MultiscalePartitioningEngineView::analyzeTimescalePoolConnectivity | ( | const std::vector< std::vector< size_t > > & | timescale_pools, |
| const std::vector< double > & | Y, | ||
| double | T9, | ||
| double | rho ) const |
Analyzes the connectivity of timescale pools.
| timescale_pools | A vector of vectors of species indices, where each inner vector represents a timescale pool. |
| Y | Vector of current molar abundances for the full network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
@purpose To merge timescale pools that are strongly connected by reactions, forming cohesive groups for QSE analysis.
@how For each pool, it builds a reaction connectivity graph using buildConnectivityGraph. It then finds the connected components within that graph using a Breadth-First Search (BFS). The resulting components from all pools are collected and returned.
|
private |
Builds a connectivity graph from a set of fast reaction indices.
| fast_reaction_indices | A set of indices for reactions considered "fast". |
@purpose To represent the reaction pathways among a subset of reactions.
@how It iterates through the specified fast reactions. For each reaction, it creates a two-way edge in the graph between every reactant and every product, signifying that mass can flow between them.
|
private |
Builds a connectivity graph from a species pool.
| species_pool | A vector of species indices representing a species pool. |
@purpose To find reaction connections within a specific group of species.
@how It iterates through all reactions in the base engine. If a reaction involves at least two distinct species from the input species_pool (one as a reactant and one as a product), it adds edges between all reactants and products from that reaction that are also in the pool.
|
nodiscardoverridevirtual |
Calculates the molar reaction flow for a given reaction.
| reaction | The reaction for which to calculate the flow. |
| Y_full | Vector of current molar abundances for the full network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
@purpose To compute the net rate of a single reaction.
@how It first checks the QSE cache. On a hit, it retrieves the cached equilibrium abundances for the algebraic species. It creates a mutable copy of Y_full, overwrites the algebraic species abundances with the cached equilibrium values, and then calls the base engine's calculateMolarReactionFlow with this modified abundance vector.
| StaleEngineError | If the QSE cache misses. |
Implements gridfire::DynamicEngine.
|
nodiscardoverridevirtual |
Calculates the right-hand side (dY/dt) and energy generation.
| Y_full | Vector of current molar abundances for all species in the base engine. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
std::expected containing StepDerivatives<double> on success, or a StaleEngineError if the engine's QSE cache does not contain a solution for the given state.@purpose To compute the time derivatives for the ODE solver. This implementation modifies the derivatives from the base engine to enforce the QSE condition.
@how It first performs a lookup in the QSE abundance cache (m_qse_abundance_cache). If a cache hit occurs, it calls the base engine's calculateRHSAndEnergy. It then manually sets the time derivatives (dydt) of all identified algebraic species to zero, effectively removing their differential equations from the system being solved.
update() or equilibrateNetwork() for the current thermodynamic conditions, so that a valid entry exists in the QSE cache. dydt=0 for all algebraic species.| StaleEngineError | If the QSE cache does not contain an entry for the given (T9, rho, Y_full). This indicates update() was not called recently enough. |
Implements gridfire::Engine.
|
private |
Constructs candidate QSE groups from connected timescale pools.
| candidate_pools | A vector of vectors of species indices, where each inner vector represents a connected pool of species with similar fast timescales. |
| Y | Vector of current molar abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
QSEGroup structs, ready for flux validation.@how For each input pool, it identifies "bridge" reactions that connect the pool to species outside the pool. The reactants of these bridge reactions that are not in the pool are identified as "seed" species. The original pool members are the "algebraic" species. It then bundles the seed and algebraic species into a QSEGroup struct.
candidate_pools should be connected components from analyzeTimescalePoolConnectivity. QSEGroup objects is returned. | fourdst::composition::Composition gridfire::MultiscalePartitioningEngineView::equilibrateNetwork | ( | const NetIn & | netIn | ) |
Equilibrates the network using QSE from a NetIn struct.
| netIn | A struct containing the current network input. |
@purpose A convenience overload for equilibrateNetwork.
@how It unpacks the netIn struct into Y, T9, and rho and then calls the primary equilibrateNetwork method.
| fourdst::composition::Composition gridfire::MultiscalePartitioningEngineView::equilibrateNetwork | ( | const std::vector< double > & | Y, |
| double | T9, | ||
| double | rho ) |
Equilibrates the network by partitioning and solving for QSE abundances.
| Y | Vector of current molar abundances for the full network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
@purpose A convenience method to run the full QSE analysis and get an equilibrated composition object as a result.
@how It first calls partitionNetwork() with the given state to define the QSE groups. Then, it calls solveQSEAbundances() to compute the new equilibrium abundances for the algebraic species. Finally, it packs the resulting full abundance vector into a new fourdst::composition::Composition object and returns it.
| void gridfire::MultiscalePartitioningEngineView::exportToDot | ( | const std::string & | filename, |
| const std::vector< double > & | Y, | ||
| const double | T9, | ||
| const double | rho ) const |
Exports the network to a DOT file for visualization.
| filename | The name of the DOT file to create. |
| Y | Vector of current molar abundances for the full network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
@purpose To visualize the partitioned network graph.
@how This method delegates the DOT file export to the base engine. It does not currently add any partitioning information to the output graph.
|
overridevirtual |
Generates the Jacobian matrix for the current state.
| Y_full | Vector of current molar abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
@purpose To compute the Jacobian matrix required by implicit ODE solvers.
@how It first performs a QSE cache lookup. On a hit, it delegates the full Jacobian calculation to the base engine. While this view could theoretically return a modified, sparser Jacobian reflecting the QSE constraints, the current implementation returns the full Jacobian from the base engine. The solver is expected to handle the algebraic constraints (e.g., via dydt=0 from calculateRHSAndEnergy).
| exceptions::StaleEngineError | If the QSE cache misses, as it cannot proceed without a valid partition. |
Implements gridfire::DynamicEngine.
|
overridevirtual |
Generates the stoichiometry matrix for the network.
@purpose To prepare the stoichiometry matrix for later queries.
@how This method delegates directly to the base engine's generateStoichiometryMatrix(). The stoichiometry is based on the full, unpartitioned network.
Implements gridfire::DynamicEngine.
|
overridevirtual |
Gets the base engine.
Implements gridfire::EngineView< DynamicEngine >.
|
nodiscard |
Gets the dynamic species in the network.
@purpose To allow external queries of the partitioning results.
@how It returns a const reference to the m_dynamic_species member vector.
partitionNetwork() must have been called.
|
nodiscard |
Gets the fast species in the network.
@purpose To allow external queries of the partitioning results.
@how It returns a copy of the m_algebraic_species member vector.
partitionNetwork() must have been called.
|
nodiscardoverridevirtual |
Gets an entry from the previously generated Jacobian matrix.
| i_full | Row index (species index) in the full network. |
| j_full | Column index (species index) in the full network. |
@purpose To provide Jacobian entries to an implicit solver.
@how This method directly delegates to the base engine's getJacobianMatrixEntry. It does not currently modify the Jacobian to reflect the QSE algebraic constraints, as these are handled by setting dY/dt = 0 in calculateRHSAndEnergy.
generateJacobianMatrix() must have been called for the current state. Implements gridfire::DynamicEngine.
|
nodiscardoverridevirtual |
Gets the set of logical reactions in the network.
LogicalReactionSet from the base engine, containing all reactions in the full network. Implements gridfire::DynamicEngine.
|
nodiscardoverridevirtual |
Gets the list of species in the network.
Species objects representing all species in the underlying base engine. This view does not alter the species list itself, only how their abundances are evolved. Implements gridfire::Engine.
|
nodiscardoverridevirtual |
Gets the current electron screening model.
@how This method delegates directly to the base engine's getScreeningModel().
Implements gridfire::DynamicEngine.
|
nodiscardoverridevirtual |
Computes destruction timescales for all species in the network.
| Y | Vector of current molar abundances for the full network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
std::expected containing a map from Species to their characteristic destruction timescales (s) on success, or a StaleEngineError on failure.@purpose To get the timescale for species destruction, which is used as the primary metric for network partitioning.
@how It delegates the calculation to the base engine. For any species identified as algebraic (in QSE), it manually sets their timescale to 0.0.
| StaleEngineError | If the QSE cache misses. |
Implements gridfire::DynamicEngine.
|
nodiscardoverridevirtual |
Gets the index of a species in the full network.
| species | The species to get the index of. |
@how This method delegates directly to the base engine's getSpeciesIndex().
Implements gridfire::DynamicEngine.
|
nodiscardoverridevirtual |
Computes timescales for all species in the network.
| Y | Vector of current molar abundances for the full network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
std::expected containing a map from Species to their characteristic timescales (s) on success, or a StaleEngineError on failure.@purpose To get the characteristic timescale Y / (dY/dt) for each species.
@how It delegates the calculation to the base engine. For any species identified as algebraic (in QSE), it manually sets their timescale to 0.0 to signify that they equilibrate instantaneously on the timescale of the solver.
| StaleEngineError | If the QSE cache misses. |
Implements gridfire::DynamicEngine.
|
nodiscardoverridevirtual |
Gets an entry from the stoichiometry matrix.
| speciesIndex | Index of the species in the full network. |
| reactionIndex | Index of the reaction in the full network. |
@purpose To query the stoichiometric relationship between a species and a reaction.
@how This method delegates directly to the base engine's getStoichiometryMatrixEntry().
generateStoichiometryMatrix() must have been called. Implements gridfire::DynamicEngine.
|
private |
Identifies the pool with the slowest mean timescale.
| pools | A vector of vectors of species indices, where each inner vector represents a timescale pool. |
| Y | Vector of current molar abundances for the full network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
@purpose To identify the core set of dynamic species that will not be part of any QSE group.
@how It calculates the geometric mean of the destruction timescales for all species in each pool and returns the index of the pool with the maximum mean timescale.
|
overridevirtual |
Checks if the engine's internal state is stale relative to the provided conditions.
| netIn | A struct containing the current network input. |
true if the engine is stale, false otherwise.@purpose To determine if update() needs to be called.
@how It creates a QSECacheKey from the netIn data and checks for its existence in the m_qse_abundance_cache. A cache miss indicates the engine is stale because it does not have a valid QSE partition for the current conditions. It also queries the base engine's isStale() method.
Implements gridfire::DynamicEngine.
|
nodiscardoverridevirtual |
Maps a NetIn struct to a molar abundance vector for the full network.
| netIn | A struct containing the current network input. |
@how This method delegates directly to the base engine's mapNetInToMolarAbundanceVector().
Implements gridfire::DynamicEngine.
|
private |
Partitions the network by timescale.
| Y_full | Vector of current molar abundances for all species. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
@purpose To group species into "pools" based on their destruction timescales.
@how It retrieves all species destruction timescales from the base engine, sorts them, and then iterates through the sorted list, creating a new pool whenever it detects a gap between consecutive timescales that is larger than a predefined threshold (e.g., a factor of 100).
| void gridfire::MultiscalePartitioningEngineView::partitionNetwork | ( | const NetIn & | netIn | ) |
Partitions the network based on timescales from a NetIn struct.
| netIn | A struct containing the current network input. |
@purpose A convenience overload for partitionNetwork.
@how It unpacks the netIn struct into Y, T9, and rho and then calls the primary partitionNetwork method.
| void gridfire::MultiscalePartitioningEngineView::partitionNetwork | ( | const std::vector< double > & | Y, |
| double | T9, | ||
| double | rho ) |
Partitions the network into dynamic and algebraic (QSE) groups based on timescales.
| Y | Vector of current molar abundances for the full network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
@purpose To perform the core partitioning logic that identifies which species are "fast" (and can be treated algebraically) and which are "slow" (and must be integrated dynamically).
@how
partitionByTimescale: Gets species destruction timescales from the base engine, sorts them, and looks for large gaps to create timescale "pools".identifyMeanSlowestPool: The pool with the slowest average timescale is designated as the core set of dynamic species.analyzeTimescalePoolConnectivity: The other (faster) pools are analyzed for reaction connectivity to form cohesive groups.constructCandidateGroups: These connected groups are processed to identify "seed" species (dynamic species that feed the group) and "algebraic" species (the rest).validateGroupsWithFluxAnalysis: The groups are validated by ensuring their internal reaction flux is much larger than the flux connecting them to the outside network.m_qse_groups, m_dynamic_species, and m_algebraic_species (and their index maps) are populated with the results of the partitioning.
|
nodiscardoverridevirtual |
Primes the engine with a specific species.
| netIn | A struct containing the current network input. |
PrimingReport struct containing information about the priming process.@purpose To prepare the network for ignition or specific pathway studies.
@how This method delegates directly to the base engine's primeEngine(). The multiscale view does not currently interact with the priming process.
Implements gridfire::DynamicEngine.
|
overridevirtual |
Sets the set of logical reactions in the network.
| reactions | The set of logical reactions to use. |
@purpose To modify the reaction network.
@how This operation is not supported by the MultiscalePartitioningEngineView as it would invalidate the partitioning logic. It logs a critical error and throws an exception. Network modifications should be done on the base engine before it is wrapped by this view.
| exceptions::UnableToSetNetworkReactionsError | Always. |
Implements gridfire::DynamicEngine.
|
overridevirtual |
Sets the electron screening model.
| model | The type of screening model to use for reaction rate calculations. |
@how This method delegates directly to the base engine's setScreeningModel().
Implements gridfire::DynamicEngine.
|
private |
Solves for the QSE abundances of the algebraic species in a given state.
| Y_full | Vector of current molar abundances for all species in the base engine. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
@purpose To find the equilibrium abundances of the algebraic species that satisfy the QSE conditions.
@how It uses the Levenberg-Marquardt algorithm via Eigen's LevenbergMarquardt class. The problem is defined by the EigenFunctor which computes the residuals and Jacobian for the QSE equations.
|
overridevirtual |
Updates the internal state of the engine, performing partitioning and QSE equilibration.
| netIn | A struct containing the current network input: temperature, density, and composition. |
@purpose This is the main entry point for preparing the multiscale engine for use. It triggers the network partitioning and solves for the initial QSE abundances, caching the result.
@how
equilibrateNetwork().equilibrateNetwork() in turn calls partitionNetwork() to define the dynamic and algebraic species sets.solveQSEAbundances() to compute the equilibrium abundances.m_qse_abundance_cache.fourdst::composition::Composition object reflecting the equilibrated state is created and returned.netIn struct should contain a valid physical state. m_dynamic_species, m_algebraic_species, etc. are populated). The m_qse_abundance_cache is populated with the QSE solution for the given state. The returned composition reflects the new equilibrium. Implements gridfire::DynamicEngine.
|
private |
Validates candidate QSE groups using flux analysis.
| candidate_groups | A vector of candidate QSE groups. |
| Y | Vector of current molar abundances for the full network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
@purpose To ensure that a candidate QSE group is truly in equilibrium by checking that the reaction fluxes within the group are much larger than the fluxes leaving the group.
@how For each candidate group, it calculates the sum of all internal reaction fluxes and the sum of all external (bridge) reaction fluxes. If the ratio of internal to external flux exceeds a configurable threshold, the group is considered valid and is added to the returned vector.
|
private |
Indices of all reactions involving only active species.
|
private |
Indices of all species considered active in the current partition (dynamic + algebraic).
|
private |
Species that are treated as algebraic (in QSE) in the QSE groups.
|
private |
Indices of algebraic species in the full network.
|
private |
The base engine to which this view delegates calculations.
|
mutableprivate |
Statistics for the QSE abundance cache.
|
private |
The simplified set of species presented to the solver (the "slow" species).
|
private |
Indices mapping the dynamic species back to the base engine's full species list.
|
private |
Logger instance for logging messages.
|
mutableprivate |
Cache for QSE abundances based on T9, rho, and Y.
@purpose This is the core of the caching mechanism. It stores the results of QSE solves to avoid re-computation. The key is a QSECacheKey which hashes the thermodynamic state, and the value is the vector of solved molar abundances for the algebraic species.
|
private |
The list of identified equilibrium groups.