GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
gridfire::DynamicEngine Class Referenceabstract

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

#include <engine_abstract.h>

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

Public Member Functions

virtual void generateJacobianMatrix (const std::vector< double > &Y, double T9, double rho)=0
 Generate the Jacobian matrix for the current state.
 
virtual double getJacobianMatrixEntry (int i, int j) const =0
 Get an entry from the previously generated Jacobian matrix.
 
virtual void generateStoichiometryMatrix ()=0
 Generate the stoichiometry matrix for the network.
 
virtual int getStoichiometryMatrixEntry (int speciesIndex, int reactionIndex) const =0
 Get an entry from the stoichiometry matrix.
 
virtual double calculateMolarReactionFlow (const reaction::Reaction &reaction, const std::vector< double > &Y, double T9, double rho) const =0
 Calculate the molar reaction flow for a given reaction.
 
virtual const reaction::LogicalReactionSetgetNetworkReactions () const =0
 Get the set of logical reactions in the network.
 
virtual std::unordered_map< fourdst::atomic::Species, double > getSpeciesTimescales (const std::vector< double > &Y, double T9, double rho) const =0
 Compute timescales for all species in the network.
 
virtual void update (const NetIn &netIn)=0
 Update the internal state of the engine.
 
virtual void setScreeningModel (screening::ScreeningType model)=0
 Set the electron screening model.
 
virtual screening::ScreeningType getScreeningModel () const =0
 Get the current electron screening model.
 
- Public Member Functions inherited from gridfire::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 StepDerivatives< double > calculateRHSAndEnergy (const std::vector< double > &Y, 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.

Definition at line 121 of file engine_abstract.h.

Member Function Documentation

◆ calculateMolarReactionFlow()

virtual double gridfire::DynamicEngine::calculateMolarReactionFlow ( const reaction::Reaction & reaction,
const std::vector< double > & Y,
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.
YVector of 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::AdaptiveEngineView, gridfire::FileDefinedEngineView, and gridfire::GraphEngine.

◆ generateJacobianMatrix()

virtual void gridfire::DynamicEngine::generateJacobianMatrix ( const std::vector< double > & Y,
double T9,
double rho )
pure virtual

Generate the Jacobian matrix for the current state.

Parameters
YVector of 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::AdaptiveEngineView, gridfire::FileDefinedEngineView, and gridfire::GraphEngine.

◆ generateStoichiometryMatrix()

virtual void gridfire::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::AdaptiveEngineView, gridfire::FileDefinedEngineView, and gridfire::GraphEngine.

◆ getJacobianMatrixEntry()

virtual double gridfire::DynamicEngine::getJacobianMatrixEntry ( int i,
int j ) const
nodiscardpure virtual

Get an entry from the previously generated Jacobian matrix.

Parameters
iRow index (species index).
jColumn index (species index).
Returns
Value of the Jacobian matrix at (i, j).

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

Implemented in gridfire::AdaptiveEngineView, gridfire::FileDefinedEngineView, and gridfire::GraphEngine.

◆ getNetworkReactions()

virtual const reaction::LogicalReactionSet & gridfire::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::AdaptiveEngineView, gridfire::FileDefinedEngineView, and gridfire::GraphEngine.

◆ getScreeningModel()

virtual screening::ScreeningType gridfire::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.

Implemented in gridfire::AdaptiveEngineView, gridfire::FileDefinedEngineView, and gridfire::GraphEngine.

◆ getSpeciesTimescales()

virtual std::unordered_map< fourdst::atomic::Species, double > gridfire::DynamicEngine::getSpeciesTimescales ( const std::vector< double > & Y,
double T9,
double rho ) const
nodiscardpure virtual

Compute timescales for all species in the network.

Parameters
YVector of 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::AdaptiveEngineView, gridfire::FileDefinedEngineView, and gridfire::GraphEngine.

◆ getStoichiometryMatrixEntry()

virtual int gridfire::DynamicEngine::getStoichiometryMatrixEntry ( int speciesIndex,
int reactionIndex ) const
nodiscardpure virtual

Get an entry from the stoichiometry matrix.

Parameters
speciesIndexIndex of the species.
reactionIndexIndex of the reaction.
Returns
Stoichiometric coefficient for the species in the reaction.

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

Implemented in gridfire::AdaptiveEngineView, gridfire::FileDefinedEngineView, and gridfire::GraphEngine.

◆ setScreeningModel()

virtual void gridfire::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).
Postcondition
The engine will use the specified screening model for subsequent rate calculations.

Implemented in gridfire::AdaptiveEngineView, gridfire::FileDefinedEngineView, and gridfire::GraphEngine.

◆ update()

virtual void gridfire::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::AdaptiveEngineView, gridfire::FileDefinedEngineView, and gridfire::GraphEngine.


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