GridFire 0.6.0
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
 Get the index of a species in the network.
 
std::vector< double > mapNetInToMolarAbundanceVector (const gridfire::NetIn &netIn) const override
 Map a NetIn object to a vector of molar abundances.
 
gridfire::PrimingReport primeEngine (const gridfire::NetIn &netIn) override
 Prime the engine with initial conditions.
 
gridfire::BuildDepthType getDepth () const override
 Get the depth of the network.
 
void rebuild (const fourdst::composition::Composition &comp, gridfire::BuildDepthType depth) override
 Rebuild the network with a specified depth.
 
- Public Member Functions inherited from gridfire::Engine
virtual ~Engine ()=default
 Virtual destructor.
 

Private Attributes

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

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.

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

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

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

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

◆ getDepth()

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

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 from gridfire::DynamicEngine.

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

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

◆ getNetworkSpecies()

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

PyDynamicEngine Implementation ///.

Implements gridfire::Engine.

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

◆ 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

◆ getSpeciesIndex()

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

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.

Implements gridfire::DynamicEngine.

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

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

◆ isStale()

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

◆ mapNetInToMolarAbundanceVector()

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

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.

Implements gridfire::DynamicEngine.

◆ primeEngine()

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

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.

Implements gridfire::DynamicEngine.

◆ rebuild()

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

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 from gridfire::DynamicEngine.

◆ setNetworkReactions()

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

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

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

Member Data Documentation

◆ m_species_cache

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

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