GridFire v0.7.0_rc2
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 (const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const =0
 Generate the Jacobian matrix for the current state.
 
virtual NetworkJacobian generateJacobianMatrix (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 (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 void generateStoichiometryMatrix ()=0
 Generate the stoichiometry matrix for the network.
 
virtual int getStoichiometryMatrixEntry (const fourdst::atomic::Species &species, const reaction::Reaction &reaction) const =0
 Get an entry from the stoichiometry matrix.
 
virtual double calculateMolarReactionFlow (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 (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 () const =0
 Get the set of logical reactions in the network.
 
virtual void setNetworkReactions (const reaction::ReactionSet &reactions)=0
 Set the reactions for the network.
 
virtual std::expected< std::unordered_map< fourdst::atomic::Species, double >, EngineStatusgetSpeciesTimescales (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 (const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const =0
 Compute destruction timescales for all species in the network.
 
virtual fourdst::composition::Composition update (const NetIn &netIn)=0
 Update the internal state of the engine.
 
virtual bool isStale (const NetIn &netIn)=0
 Check if the engine's internal state is stale.
 
virtual void setScreeningModel (screening::ScreeningType model)=0
 Set the electron screening model.
 
virtual screening::ScreeningType getScreeningModel () const =0
 Get the current electron screening model.
 
virtual size_t getSpeciesIndex (const fourdst::atomic::Species &species) const =0
 Get the index of a species in the network.
 
virtual std::vector< double > mapNetInToMolarAbundanceVector (const NetIn &netIn) const =0
 Map a NetIn object to a vector of molar abundances.
 
virtual PrimingReport primeEngine (const NetIn &netIn)=0
 Prime the engine with initial conditions.
 
virtual BuildDepthType getDepth () const
 Get the depth of the network.
 
virtual void rebuild (const fourdst::composition::CompositionAbstract &comp, BuildDepthType depth)
 Rebuild the network with a specified depth.
 
virtual fourdst::composition::Composition collectComposition (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 (const fourdst::atomic::Species &species) const =0
 Get the status of a species in the network.
 
- Public Member Functions inherited from gridfire::engine::Engine
virtual ~Engine ()=default
 Virtual destructor.
 
virtual const std::vector< fourdst::atomic::Species > & getNetworkSpecies () const =0
 Get the list of species in the network.
 
virtual std::expected< StepDerivatives< double >, EngineStatuscalculateRHSAndEnergy (const fourdst::composition::CompositionAbstract &comp, double T9, double rho) 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 ( 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 ( 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
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 ( const fourdst::composition::CompositionAbstract & comp,
double T9,
double rho ) const
pure 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
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.

◆ generateJacobianMatrix() [1/3]

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

Generate 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 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 ( 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
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 ( 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
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.

◆ generateStoichiometryMatrix()

virtual void gridfire::engine::DynamicEngine::generateStoichiometryMatrix ( )
pure virtual

Generate the stoichiometry matrix for the network.

This method must compute and store the stoichiometry matrix, which encodes the net change of each species in each reaction.

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

◆ getDepth()

virtual BuildDepthType gridfire::engine::DynamicEngine::getDepth ( ) const
inlinenodiscardvirtual

Get the depth of the network.

Returns
The depth of the network, which may indicate the level of detail or complexity in the reaction network.

This method is intended to provide information about the network's structure, such as how many layers of reactions or species are present. It can be useful for diagnostics and understanding the network's complexity.

Reimplemented in gridfire::engine::GraphEngine, and PyDynamicEngine.

◆ getNetworkReactions()

virtual const reaction::ReactionSet & gridfire::engine::DynamicEngine::getNetworkReactions ( ) 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 ( ) const
nodiscardpure virtual

Get the current electron screening model.

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 ( const fourdst::composition::CompositionAbstract & comp,
double T9,
double rho ) const
nodiscardpure virtual

Compute 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.

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

◆ getSpeciesIndex()

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

Get the index of a species in the network.

Parameters
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 ( 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 ( 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.

◆ getStoichiometryMatrixEntry()

virtual int gridfire::engine::DynamicEngine::getStoichiometryMatrixEntry ( const fourdst::atomic::Species & species,
const reaction::Reaction & reaction ) const
nodiscardpure virtual

Get 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().

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

◆ isStale()

virtual bool gridfire::engine::DynamicEngine::isStale ( const NetIn & netIn)
nodiscardpure virtual

Check if the engine's internal state is stale.

Parameters
netInA struct containing the current network input, such as temperature, density, and composition.
Returns
True if the engine's state is stale and needs to be updated; false otherwise.

This method allows derived classes to determine if their internal state is out-of-date with respect to the provided network conditions. If the engine is stale, it may require a call to update() before performing calculations.

Usage Example:
NetIn input = { ... };
if (myEngine.isStale(input)) {
// Update the engine before proceeding
}
Definition types.h:27

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

◆ mapNetInToMolarAbundanceVector()

virtual std::vector< double > gridfire::engine::DynamicEngine::mapNetInToMolarAbundanceVector ( const NetIn & netIn) const
nodiscardpure virtual

Map a NetIn object to a vector of molar abundances.

Parameters
netInThe input conditions for the network.
Returns
A vector of molar abundances corresponding to the species in the network.

This method converts the input conditions into a vector of molar abundances, which can be used for further calculations or diagnostics.

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

◆ primeEngine()

virtual PrimingReport gridfire::engine::DynamicEngine::primeEngine ( const NetIn & netIn)
nodiscardpure virtual

Prime the engine with initial conditions.

Parameters
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.

◆ rebuild()

virtual void gridfire::engine::DynamicEngine::rebuild ( const fourdst::composition::CompositionAbstract & comp,
BuildDepthType depth )
inlinevirtual

Rebuild the network with a specified depth.

Parameters
compThe composition to rebuild the network with.
depthThe desired depth of the network.

This method is intended to allow dynamic adjustment of the network's depth, which may involve adding or removing species and reactions based on the specified depth. However, not all engines support this operation.

Reimplemented in gridfire::engine::GraphEngine, and PyDynamicEngine.

◆ setNetworkReactions()

virtual void gridfire::engine::DynamicEngine::setNetworkReactions ( const reaction::ReactionSet & reactions)
pure virtual

Set 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.

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

◆ setScreeningModel()

virtual void gridfire::engine::DynamicEngine::setScreeningModel ( screening::ScreeningType model)
pure virtual

Set the electron screening model.

Parameters
modelThe type of screening model to use for reaction rate calculations.

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.

Usage Example:
myEngine.setScreeningModel(screening::ScreeningType::WEAK);
@ WEAK
Weak screening model (Salpeter, 1954).
Definition screening_types.h:35
Postcondition
The engine will use the specified screening model for subsequent rate calculations.

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

◆ update()

virtual fourdst::composition::Composition gridfire::engine::DynamicEngine::update ( const NetIn & netIn)
pure virtual

Update the internal state of the engine.

Parameters
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);
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: