GridFire v0.7.6rc4.0
General Purpose Nuclear Network
Loading...
Searching...
No Matches
gridfire::engine::DynamicEngine Class Referenceabstract

Abstract class for engines supporting Jacobian and stoichiometry operations. More...

#include <engine_abstract.h>

Inheritance diagram for gridfire::engine::DynamicEngine:
[legend]
Collaboration diagram for gridfire::engine::DynamicEngine:
[legend]

Public Member Functions

virtual NetworkJacobian generateJacobianMatrix (scratch::StateBlob &ctx, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const =0
 Generate the Jacobian matrix for the current state.
 
virtual NetworkJacobian generateJacobianMatrix (scratch::StateBlob &ctx, const fourdst::composition::CompositionAbstract &comp, double T9, double rho, const std::vector< fourdst::atomic::Species > &activeSpecies) const =0
 Generate the Jacobian matrix for the current state using a subset of active species.
 
virtual NetworkJacobian generateJacobianMatrix (scratch::StateBlob &ctx, const fourdst::composition::CompositionAbstract &comp, double T9, double rho, const SparsityPattern &sparsityPattern) const =0
 Generate the Jacobian matrix for the current state with a specified sparsity pattern.
 
virtual double calculateMolarReactionFlow (scratch::StateBlob &ctx, const reaction::Reaction &reaction, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const =0
 Calculate the molar reaction flow for a given reaction.
 
virtual EnergyDerivatives calculateEpsDerivatives (scratch::StateBlob &ctx, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const =0
 Calculate the derivatives of the energy generation rate with respect to T and rho.
 
virtual const reaction::ReactionSetgetNetworkReactions (scratch::StateBlob &ctx) const =0
 Get the set of logical reactions in the network.
 
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
 
virtual std::expected< std::unordered_map< fourdst::atomic::Species, double >, EngineStatusgetSpeciesTimescales (scratch::StateBlob &ctx, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const =0
 Compute timescales for all species in the network.
 
virtual std::expected< std::unordered_map< fourdst::atomic::Species, double >, EngineStatusgetSpeciesDestructionTimescales (scratch::StateBlob &ctx, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const =0
 Compute destruction timescales for all species in the network.
 
virtual fourdst::composition::Composition project (scratch::StateBlob &ctx, const NetIn &netIn) const =0
 Update the thread local scratch pad state of a network.
 
virtual screening::ScreeningType getScreeningModel (scratch::StateBlob &ctx) const =0
 Get the current electron screening model.
 
virtual size_t getSpeciesIndex (scratch::StateBlob &ctx, const fourdst::atomic::Species &species) const =0
 Get the index of a species in the network.
 
virtual PrimingReport primeEngine (scratch::StateBlob &ctx, const NetIn &netIn) const =0
 Prime the engine with initial conditions.
 
virtual fourdst::composition::Composition collectComposition (scratch::StateBlob &ctx, const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const =0
 Recursively collect composition from current engine and any sub engines if they exist.
 
virtual SpeciesStatus getSpeciesStatus (scratch::StateBlob &ctx, const fourdst::atomic::Species &species) const =0
 Get the status of a species in the network.
 
virtual std::optional< StepDerivatives< double > > getMostRecentRHSCalculation (scratch::StateBlob &ctx) const =0
 
virtual std::unique_ptr< scratch::StateBlobconstructStateBlob (const scratch::StateBlob *blob) const =0
 
- Public Member Functions inherited from gridfire::engine::Engine
virtual ~Engine ()=default
 Virtual destructor.
 
virtual const std::vector< fourdst::atomic::Species > & getNetworkSpecies (scratch::StateBlob &ctx) const =0
 Get the list of species in the network.
 
virtual std::expected< StepDerivatives< double >, EngineStatuscalculateRHSAndEnergy (scratch::StateBlob &, const fourdst::composition::CompositionAbstract &comp, double T9, double rho, bool trust) const =0
 Calculate the right-hand side (dY/dt) and energy generation.
 

Detailed Description

Abstract class for engines supporting Jacobian and stoichiometry operations.

Extends Engine with additional methods for:

  • Generating and accessing the Jacobian matrix (for implicit solvers).
  • Generating and accessing the stoichiometry matrix.
  • Calculating molar reaction flows for individual reactions.
  • Accessing the set of logical reactions in the network.
  • Computing timescales for each species.

Intended usage: Derive from this class to implement engines that support advanced solver features such as implicit integration, sensitivity analysis, QSE (Quasi-Steady-State Equilibrium) handling, and more. Generally this will be the main engine type

Member Function Documentation

◆ calculateEpsDerivatives()

virtual EnergyDerivatives gridfire::engine::DynamicEngine::calculateEpsDerivatives ( scratch::StateBlob & ctx,
const fourdst::composition::CompositionAbstract & comp,
double T9,
double rho ) const
nodiscardpure virtual

Calculate the derivatives of the energy generation rate with respect to T and rho.

Parameters
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
Returns
EnergyDerivatives containing dEps/dT and dEps/dRho.

This method computes the partial derivatives of the specific nuclear energy generation rate with respect to temperature and density for the current state.

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ calculateMolarReactionFlow()

virtual double gridfire::engine::DynamicEngine::calculateMolarReactionFlow ( scratch::StateBlob & ctx,
const reaction::Reaction & reaction,
const fourdst::composition::CompositionAbstract & comp,
double T9,
double rho ) const
nodiscardpure virtual

Calculate the molar reaction flow for a given reaction.

Parameters
ctxThe scratchpad context for the current state.
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.

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ collectComposition()

virtual fourdst::composition::Composition gridfire::engine::DynamicEngine::collectComposition ( scratch::StateBlob & ctx,
const fourdst::composition::CompositionAbstract & comp,
double T9,
double rho ) const
nodiscardpure virtual

Recursively collect composition from current engine and any sub engines if they exist.

If species i is defined in comp and in any sub engine or self composition then the molar abundance of species i in the returned composition will be that defined in comp. If there are species defined in sub engine compositions which are not defined in comp then their molar abundances will be based on the reported values from each sub engine.

Note
It is up to each engine to decide how to handle filling in the return composition.
These methods return an unfinalized composition which must then be finalized by the caller
Parameters
ctxThe scratchpad context for the current state.
compInput composition to "normalize".
T9
rho
Returns
An updated composition which is a superset of comp. This may contain species which were culled, for example, by either QSE partitioning or reaction flow rate culling

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ constructStateBlob()

virtual std::unique_ptr< scratch::StateBlob > gridfire::engine::DynamicEngine::constructStateBlob ( const scratch::StateBlob * blob) const
nodiscardpure virtual

◆ generateJacobianMatrix() [1/3]

virtual NetworkJacobian gridfire::engine::DynamicEngine::generateJacobianMatrix ( scratch::StateBlob & ctx,
const fourdst::composition::CompositionAbstract & comp,
double T9,
double rho ) const
nodiscardpure virtual

Generate the Jacobian matrix for the current state.

Parameters
ctxThe scratchpad context for the current state.
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.

This method must compute and store the Jacobian matrix (∂(dY/dt)_i/∂Y_j) for the current state. The matrix can then be accessed via getJacobianMatrixEntry().

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ generateJacobianMatrix() [2/3]

virtual NetworkJacobian gridfire::engine::DynamicEngine::generateJacobianMatrix ( scratch::StateBlob & ctx,
const fourdst::composition::CompositionAbstract & comp,
double T9,
double rho,
const SparsityPattern & sparsityPattern ) const
nodiscardpure virtual

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

Parameters
ctxGet the scratchpad context for the current state.
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 must compute and store 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()

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ generateJacobianMatrix() [3/3]

virtual NetworkJacobian gridfire::engine::DynamicEngine::generateJacobianMatrix ( scratch::StateBlob & ctx,
const fourdst::composition::CompositionAbstract & comp,
double T9,
double rho,
const std::vector< fourdst::atomic::Species > & activeSpecies ) const
nodiscardpure virtual

Generate the Jacobian matrix for the current state using a subset of active species.

Parameters
ctxThe scratchpad context for the current state.
compComposition object containing current abundances.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
activeSpeciesThe set of species to include in the Jacobian calculation.

This method must compute and store the Jacobian matrix (∂(dY/dt)_i/∂Y_j) for the current state, considering only the specified subset of active species. The matrix can then be accessed via getJacobianMatrixEntry().

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ getInactiveNetworkReactions()

virtual reaction::ReactionSet gridfire::engine::DynamicEngine::getInactiveNetworkReactions ( scratch::StateBlob & ctx) const
inlinenodiscardvirtual

Get the set of inactive reactions in the network.

Returns
ReactionSet containing all inactive reactions.

By default, this method returns an empty set. Derived classes can override this method to provide the actual set of inactive reactions based on their internal logic (e.g., reaction flow culling, QSE partitioning).

Reimplemented in gridfire::engine::AdaptiveEngineView.

◆ getInactiveReactionMolarReactionFlow()

virtual double gridfire::engine::DynamicEngine::getInactiveReactionMolarReactionFlow ( scratch::StateBlob & ctx,
const reaction::Reaction & reaction,
const fourdst::composition::CompositionAbstract & comp,
const double T9,
const double rho ) const
inlinenodiscardvirtual

◆ getMostRecentRHSCalculation()

virtual std::optional< StepDerivatives< double > > gridfire::engine::DynamicEngine::getMostRecentRHSCalculation ( scratch::StateBlob & ctx) const
nodiscardpure virtual

◆ getNetworkReactions()

virtual const reaction::ReactionSet & gridfire::engine::DynamicEngine::getNetworkReactions ( scratch::StateBlob & ctx) const
nodiscardpure virtual

Get the set of logical reactions in the network.

Returns
Reference to the LogicalReactionSet containing all reactions.

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ getScreeningModel()

virtual screening::ScreeningType gridfire::engine::DynamicEngine::getScreeningModel ( scratch::StateBlob & ctx) const
nodiscardpure virtual

Get the current electron screening model.

Parameters
ctxThe scratchpad context for the current state.
Returns
The currently active screening model type.
Usage Example:
screening::ScreeningType currentModel = myEngine.getScreeningModel();
ScreeningType
Enumerates the available plasma screening models.
Definition screening_types.h:15

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ getSpeciesDestructionTimescales()

virtual std::expected< std::unordered_map< fourdst::atomic::Species, double >, EngineStatus > gridfire::engine::DynamicEngine::getSpeciesDestructionTimescales ( scratch::StateBlob & ctx,
const fourdst::composition::CompositionAbstract & comp,
double T9,
double rho ) const
nodiscardpure virtual

Compute destruction timescales for all species in the network.

Parameters
ctxThe scratchpad context for the current state.
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.

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ getSpeciesIndex()

virtual size_t gridfire::engine::DynamicEngine::getSpeciesIndex ( scratch::StateBlob & ctx,
const fourdst::atomic::Species & species ) const
nodiscardpure virtual

Get the index of a species in the network.

Parameters
ctxThe scratchpad context for the current state.
speciesThe species to look up.

This method allows querying the index of a specific species in the engine's internal representation. It is useful for accessing species data efficiently.

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ getSpeciesStatus()

virtual SpeciesStatus gridfire::engine::DynamicEngine::getSpeciesStatus ( scratch::StateBlob & ctx,
const fourdst::atomic::Species & species ) const
nodiscardpure virtual

Get the status of a species in the network.

Parameters
speciesThe species to check.
Returns
SpeciesStatus indicating whether the species is active, inactive, or culled.

This method allows querying the current status of a specific species within the engine's network.

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ getSpeciesTimescales()

virtual std::expected< std::unordered_map< fourdst::atomic::Species, double >, EngineStatus > gridfire::engine::DynamicEngine::getSpeciesTimescales ( scratch::StateBlob & ctx,
const fourdst::composition::CompositionAbstract & comp,
double T9,
double rho ) const
nodiscardpure virtual

Compute 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, diagnostics, and reaction network culling.

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ primeEngine()

virtual PrimingReport gridfire::engine::DynamicEngine::primeEngine ( scratch::StateBlob & ctx,
const NetIn & netIn ) const
nodiscardpure virtual

Prime the engine with initial conditions.

Parameters
ctxThe scratchpad context for the current state.
netInThe input conditions for the network.
Returns
PrimingReport containing information about the priming process.

This method is used to prepare the engine for calculations by setting up initial conditions, reactions, and species. It may involve compiling reaction rates, initializing internal data structures, and performing any necessary pre-computation.

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ project()

virtual fourdst::composition::Composition gridfire::engine::DynamicEngine::project ( scratch::StateBlob & ctx,
const NetIn & netIn ) const
nodiscardpure virtual

Update the thread local scratch pad state of a network.

Parameters
ctxThe scratchpad context for the current state.
netInA struct containing the current network input, such as temperature, density, and composition.

This method is intended to be implemented by derived classes to update their internal state based on the provided network conditions. For example, an adaptive engine might use this to re-evaluate which reactions and species are active. For other engines that do not support manually updating, this method might do nothing.

Usage Example:
NetIn input = { ... };
myEngine.update(input);
Definition types.h:27
Postcondition
The internal state of the engine is updated to reflect the new conditions.

Implemented in gridfire::engine::AdaptiveEngineView, gridfire::engine::DefinedEngineView, gridfire::engine::GraphEngine, gridfire::engine::MultiscalePartitioningEngineView, and PyDynamicEngine.


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