GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
PyDynamicEngine Class Referencefinal

#include <py_engine.h>

Inheritance diagram for PyDynamicEngine:
gridfire::DynamicEngine gridfire::Engine

Public Member Functions

const std::vector< fourdst::atomic::Species > & getNetworkSpecies () const override
 PyDynamicEngine Implementation ///.
 
std::expected< gridfire::StepDerivatives< double >, gridfire::expectations::StaleEngineErrorcalculateRHSAndEnergy (const std::vector< double > &Y, double T9, double rho) const override
 Calculate the right-hand side (dY/dt) and energy generation.
 
void generateJacobianMatrix (const std::vector< double > &Y_dynamic, double T9, double rho) const override
 Generate the Jacobian matrix for the current state.
 
void generateJacobianMatrix (const std::vector< double > &Y_dynamic, double T9, double rho, const gridfire::SparsityPattern &sparsityPattern) const override
 
double getJacobianMatrixEntry (int i, int j) const override
 Get an entry from the previously generated Jacobian matrix.
 
void generateStoichiometryMatrix () override
 Generate the stoichiometry matrix for the network.
 
int getStoichiometryMatrixEntry (int speciesIndex, int reactionIndex) const override
 Get an entry from the stoichiometry matrix.
 
double calculateMolarReactionFlow (const gridfire::reaction::Reaction &reaction, const std::vector< double > &Y, double T9, double rho) const override
 Calculate the molar reaction flow for a given reaction.
 
const gridfire::reaction::LogicalReactionSetgetNetworkReactions () const override
 Get the set of logical reactions in the network.
 
void setNetworkReactions (const gridfire::reaction::LogicalReactionSet &reactions) override
 
std::expected< std::unordered_map< fourdst::atomic::Species, double >, gridfire::expectations::StaleEngineErrorgetSpeciesTimescales (const std::vector< double > &Y, double T9, double rho) const override
 Compute timescales for all species in the network.
 
std::expected< std::unordered_map< fourdst::atomic::Species, double >, gridfire::expectations::StaleEngineErrorgetSpeciesDestructionTimescales (const std::vector< double > &Y, double T9, double rho) const override
 
fourdst::composition::Composition update (const gridfire::NetIn &netIn) override
 Update the internal state of the engine.
 
bool isStale (const gridfire::NetIn &netIn) override
 
void setScreeningModel (gridfire::screening::ScreeningType model) override
 Set the electron screening model.
 
gridfire::screening::ScreeningType getScreeningModel () const override
 Get the current electron screening model.
 
int getSpeciesIndex (const fourdst::atomic::Species &species) const override
 
std::vector< double > mapNetInToMolarAbundanceVector (const gridfire::NetIn &netIn) const override
 
gridfire::PrimingReport primeEngine (const gridfire::NetIn &netIn) override
 
gridfire::BuildDepthType getDepth () const override
 
void rebuild (const fourdst::composition::Composition &comp, gridfire::BuildDepthType depth) override
 
- Public Member Functions inherited from gridfire::Engine
virtual ~Engine ()=default
 Virtual destructor.
 

Private Attributes

std::vector< fourdst::atomic::Species > m_species_cache
 

Detailed Description

Definition at line 20 of file py_engine.h.

Member Function Documentation

◆ calculateMolarReactionFlow()

double PyDynamicEngine::calculateMolarReactionFlow ( const gridfire::reaction::Reaction & reaction,
const std::vector< double > & Y,
double T9,
double rho ) const
overridevirtual

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.

Implements gridfire::DynamicEngine.

Definition at line 123 of file py_engine.cpp.

◆ calculateRHSAndEnergy()

std::expected< gridfire::StepDerivatives< double >, gridfire::expectations::StaleEngineError > PyDynamicEngine::calculateRHSAndEnergy ( const std::vector< double > & Y,
double T9,
double rho ) const
overridevirtual

Calculate the right-hand side (dY/dt) and energy generation.

Parameters
YVector of current abundances for all species.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
Returns
StepDerivatives<double> containing dY/dt and energy generation rate.

This function must be implemented by derived classes to compute the time derivatives of all species and the specific nuclear energy generation rate for the current state.

Implements gridfire::Engine.

Definition at line 70 of file py_engine.cpp.

◆ generateJacobianMatrix() [1/2]

void PyDynamicEngine::generateJacobianMatrix ( const std::vector< double > & Y_dynamic,
double T9,
double rho ) const
overridevirtual

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

Implements gridfire::DynamicEngine.

Definition at line 79 of file py_engine.cpp.

◆ generateJacobianMatrix() [2/2]

void PyDynamicEngine::generateJacobianMatrix ( const std::vector< double > & Y_dynamic,
double T9,
double rho,
const gridfire::SparsityPattern & sparsityPattern ) const
overridevirtual

Reimplemented from gridfire::DynamicEngine.

Definition at line 88 of file py_engine.cpp.

◆ generateStoichiometryMatrix()

void PyDynamicEngine::generateStoichiometryMatrix ( )
overridevirtual

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.

Implements gridfire::DynamicEngine.

Definition at line 106 of file py_engine.cpp.

◆ getDepth()

gridfire::BuildDepthType PyDynamicEngine::getDepth ( ) const
inlineoverridevirtual

Reimplemented from gridfire::DynamicEngine.

Definition at line 41 of file py_engine.h.

◆ getJacobianMatrixEntry()

double PyDynamicEngine::getJacobianMatrixEntry ( int i,
int j ) const
overridevirtual

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.

Implements gridfire::DynamicEngine.

Definition at line 97 of file py_engine.cpp.

◆ getNetworkReactions()

const gridfire::reaction::LogicalReactionSet & PyDynamicEngine::getNetworkReactions ( ) const
overridevirtual

Get the set of logical reactions in the network.

Returns
Reference to the LogicalReactionSet containing all reactions.

Implements gridfire::DynamicEngine.

Definition at line 132 of file py_engine.cpp.

◆ getNetworkSpecies()

const std::vector< fourdst::atomic::Species > & PyDynamicEngine::getNetworkSpecies ( ) const
overridevirtual

PyDynamicEngine Implementation ///.

Implements gridfire::Engine.

Definition at line 50 of file py_engine.cpp.

◆ getScreeningModel()

gridfire::screening::ScreeningType PyDynamicEngine::getScreeningModel ( ) const
overridevirtual

Get the current electron screening model.

Returns
The currently active screening model type.
Usage Example:
screening::ScreeningType currentModel = myEngine.getScreeningModel();

Implements gridfire::DynamicEngine.

Definition at line 194 of file py_engine.cpp.

◆ getSpeciesDestructionTimescales()

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

Implements gridfire::DynamicEngine.

Definition at line 158 of file py_engine.cpp.

◆ getSpeciesIndex()

int PyDynamicEngine::getSpeciesIndex ( const fourdst::atomic::Species & species) const
overridevirtual

Implements gridfire::DynamicEngine.

Definition at line 202 of file py_engine.cpp.

◆ getSpeciesTimescales()

std::expected< std::unordered_map< fourdst::atomic::Species, double >, gridfire::expectations::StaleEngineError > PyDynamicEngine::getSpeciesTimescales ( const std::vector< double > & Y,
double T9,
double rho ) const
overridevirtual

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.

Implements gridfire::DynamicEngine.

Definition at line 149 of file py_engine.cpp.

◆ getStoichiometryMatrixEntry()

int PyDynamicEngine::getStoichiometryMatrixEntry ( int speciesIndex,
int reactionIndex ) const
overridevirtual

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

Implements gridfire::DynamicEngine.

Definition at line 114 of file py_engine.cpp.

◆ isStale()

bool PyDynamicEngine::isStale ( const gridfire::NetIn & netIn)
overridevirtual

Implements gridfire::DynamicEngine.

Definition at line 176 of file py_engine.cpp.

◆ mapNetInToMolarAbundanceVector()

std::vector< double > PyDynamicEngine::mapNetInToMolarAbundanceVector ( const gridfire::NetIn & netIn) const
overridevirtual

Implements gridfire::DynamicEngine.

Definition at line 211 of file py_engine.cpp.

◆ primeEngine()

gridfire::PrimingReport PyDynamicEngine::primeEngine ( const gridfire::NetIn & netIn)
overridevirtual

Implements gridfire::DynamicEngine.

Definition at line 220 of file py_engine.cpp.

◆ rebuild()

void PyDynamicEngine::rebuild ( const fourdst::composition::Composition & comp,
gridfire::BuildDepthType depth )
inlineoverridevirtual

Reimplemented from gridfire::DynamicEngine.

Definition at line 44 of file py_engine.h.

◆ setNetworkReactions()

void PyDynamicEngine::setNetworkReactions ( const gridfire::reaction::LogicalReactionSet & reactions)
overridevirtual

Implements gridfire::DynamicEngine.

Definition at line 140 of file py_engine.cpp.

◆ setScreeningModel()

void PyDynamicEngine::setScreeningModel ( gridfire::screening::ScreeningType model)
overridevirtual

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

Implements gridfire::DynamicEngine.

Definition at line 185 of file py_engine.cpp.

◆ update()

fourdst::composition::Composition PyDynamicEngine::update ( const gridfire::NetIn & netIn)
overridevirtual

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.

Implements gridfire::DynamicEngine.

Definition at line 167 of file py_engine.cpp.

Member Data Documentation

◆ m_species_cache

std::vector<fourdst::atomic::Species> PyDynamicEngine::m_species_cache
mutableprivate

Definition at line 48 of file py_engine.h.


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