GridFire 0.6.0
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 PyDynamicEngine gridfire::AdaptiveEngineView gridfire::DefinedEngineView gridfire::GraphEngine gridfire::MultiscalePartitioningEngineView gridfire::FileDefinedEngineView gridfire::NetworkPrimingEngineView

Public Member Functions

virtual void generateJacobianMatrix (const std::vector< double > &Y_dynamic, double T9, double rho) const =0
 Generate the Jacobian matrix for the current state.
 
virtual void generateJacobianMatrix (const std::vector< double > &Y_dynamic, double T9, double rho, const SparsityPattern &sparsityPattern) const
 
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 void setNetworkReactions (const reaction::LogicalReactionSet &reactions)=0
 
virtual std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineErrorgetSpeciesTimescales (const std::vector< double > &Y, double T9, double rho) const =0
 Compute timescales for all species in the network.
 
virtual std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineErrorgetSpeciesDestructionTimescales (const std::vector< double > &Y, double T9, double rho) const =0
 
virtual fourdst::composition::Composition update (const NetIn &netIn)=0
 Update the internal state of the engine.
 
virtual bool isStale (const NetIn &netIn)=0
 
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 int 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::Composition &comp, BuildDepthType depth)
 Rebuild the network with a specified depth.
 
- 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 std::expected< StepDerivatives< double >, expectations::StaleEngineErrorcalculateRHSAndEnergy (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.

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::DefinedEngineView, gridfire::GraphEngine, gridfire::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ generateJacobianMatrix() [1/2]

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

Generate the Jacobian matrix for the current state.

Parameters
Y_dynamicVector 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::DefinedEngineView, gridfire::GraphEngine, gridfire::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ generateJacobianMatrix() [2/2]

virtual void gridfire::DynamicEngine::generateJacobianMatrix ( const std::vector< double > & Y_dynamic,
double T9,
double rho,
const SparsityPattern & sparsityPattern ) const
inlinevirtual

Reimplemented in gridfire::GraphEngine, and PyDynamicEngine.

◆ 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::DefinedEngineView, gridfire::GraphEngine, gridfire::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ getDepth()

virtual BuildDepthType gridfire::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::GraphEngine, and PyDynamicEngine.

◆ 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::DefinedEngineView, gridfire::GraphEngine, gridfire::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ 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::DefinedEngineView, gridfire::GraphEngine, gridfire::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ 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.
Definition screening_types.h:15

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

◆ getSpeciesDestructionTimescales()

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

◆ getSpeciesIndex()

virtual int gridfire::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::AdaptiveEngineView, gridfire::DefinedEngineView, gridfire::GraphEngine, gridfire::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ getSpeciesTimescales()

virtual std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > 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::DefinedEngineView, gridfire::GraphEngine, gridfire::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ 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::DefinedEngineView, gridfire::GraphEngine, gridfire::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ isStale()

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

◆ mapNetInToMolarAbundanceVector()

virtual std::vector< double > gridfire::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::AdaptiveEngineView, gridfire::DefinedEngineView, gridfire::GraphEngine, gridfire::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ primeEngine()

virtual PrimingReport gridfire::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::AdaptiveEngineView, gridfire::DefinedEngineView, gridfire::GraphEngine, gridfire::MultiscalePartitioningEngineView, and PyDynamicEngine.

◆ rebuild()

virtual void gridfire::DynamicEngine::rebuild ( const fourdst::composition::Composition & 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::GraphEngine, and PyDynamicEngine.

◆ setNetworkReactions()

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

◆ 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).
Definition screening_types.h:35
Postcondition
The engine will use the specified screening model for subsequent rate calculations.

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

◆ update()

virtual fourdst::composition::Composition 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);
Definition network.h:53
Postcondition
The internal state of the engine is updated to reflect the new conditions.

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


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