GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
gridfire::solver::QSENetworkSolver Class Referencefinal

A network solver that uses a Quasi-Steady-State Equilibrium (QSE) approach. More...

#include <solver.h>

Inheritance diagram for gridfire::solver::QSENetworkSolver:
gridfire::solver::NetworkSolverStrategy< DynamicEngine >

Classes

struct  EigenFunctor
 Functor for calculating the residual and Jacobian for the QSE species using Eigen. More...
 
struct  JacobianFunctor
 Functor for calculating the Jacobian matrix of the ODEs for the dynamic species. More...
 
struct  RHSFunctor
 Functor for calculating the right-hand side of the ODEs for the dynamic species. More...
 

Public Member Functions

NetOut evaluate (const NetIn &netIn) override
 Evaluates the network for a given timestep using the QSE approach.
 
- Public Member Functions inherited from gridfire::solver::NetworkSolverStrategy< DynamicEngine >
 NetworkSolverStrategy (DynamicEngine &engine)
 Constructor for the NetworkSolverStrategy.
 
virtual ~NetworkSolverStrategy ()=default
 Virtual destructor.
 
 NetworkSolverStrategy (DynamicEngine &engine)
 Constructor for the NetworkSolverStrategy.
 
virtual ~NetworkSolverStrategy ()=default
 Virtual destructor.
 

Private Member Functions

dynamicQSESpeciesIndices packSpeciesTypeIndexVectors (const std::vector< double > &Y, const double T9, const double rho) const
 Packs the species indices into vectors based on their type (dynamic or QSE).
 
Eigen::VectorXd calculateSteadyStateAbundances (const std::vector< double > &Y, const double T9, const double rho, const dynamicQSESpeciesIndices &indices) const
 Calculates the steady-state abundances of the QSE species.
 
NetOut initializeNetworkWithShortIgnition (const NetIn &netIn) const
 Initializes the network with a short ignition phase.
 
bool shouldUpdateView (const NetIn &conditions) const
 Determines whether the adaptive engine view should be updated.
 

Private Attributes

quill::Logger * m_logger = fourdst::logging::LogManager::getInstance().getLogger("log")
 Logger instance.
 
fourdst::config::Config & m_config = fourdst::config::Config::getInstance()
 Configuration instance.
 
bool m_isViewInitialized = false
 Flag indicating whether the adaptive engine view has been initialized.
 
NetIn m_lastSeenConditions
 The last seen input conditions.
 

Additional Inherited Members

- Protected Attributes inherited from gridfire::solver::NetworkSolverStrategy< DynamicEngine >
DynamicEnginem_engine
 The engine used by this solver strategy.
 
DynamicEnginem_engine
 The engine used by this solver strategy.
 

Detailed Description

A network solver that uses a Quasi-Steady-State Equilibrium (QSE) approach.

This solver partitions the network into "fast" species in QSE and "slow" (dynamic) species. The abundances of the fast species are determined by solving a system of algebraic equations, while the abundances of the slow species are integrated using an ODE solver. This hybrid approach is highly effective for stiff networks with disparate timescales.

The QSE solver uses an AdaptiveEngineView to dynamically cull unimportant species and reactions, which significantly improves performance for large networks.

See also
AdaptiveEngineView
DynamicEngine::getSpeciesTimescales()

Definition at line 98 of file solver.h.

Member Function Documentation

◆ calculateSteadyStateAbundances()

Eigen::VectorXd gridfire::solver::QSENetworkSolver::calculateSteadyStateAbundances ( const std::vector< double > & Y,
const double T9,
const double rho,
const dynamicQSESpeciesIndices & indices ) const
private

Calculates the steady-state abundances of the QSE species.

Parameters
YVector of current abundances for all species.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
indicesA dynamicQSESpeciesIndices struct containing the indices of the dynamic and QSE species.
Returns
An Eigen::VectorXd containing the steady-state abundances of the QSE species.

This method solves a system of algebraic equations to determine the steady-state abundances of the QSE species.

Exceptions
std::runtime_errorIf the steady-state abundances cannot be calculated.

Definition at line 205 of file solver.cpp.

◆ evaluate()

NetOut gridfire::solver::QSENetworkSolver::evaluate ( const NetIn & netIn)
overridevirtual

Evaluates the network for a given timestep using the QSE approach.

Parameters
netInThe input conditions for the network.
Returns
The output conditions after the timestep.

This method performs the following steps:

  1. Updates the adaptive engine view (if necessary).
  2. Partitions the species into dynamic and QSE species based on their timescales.
  3. Calculates the steady-state abundances of the QSE species.
  4. Integrates the ODEs for the dynamic species using a Runge-Kutta solver.
  5. Marshals the output variables into a NetOut struct.
Exceptions
std::runtime_errorIf the steady-state abundances cannot be calculated.
See also
AdaptiveEngineView::update()
packSpeciesTypeIndexVectors()
calculateSteadyStateAbundances()

Implements gridfire::solver::NetworkSolverStrategy< DynamicEngine >.

Definition at line 26 of file solver.cpp.

◆ initializeNetworkWithShortIgnition()

NetOut gridfire::solver::QSENetworkSolver::initializeNetworkWithShortIgnition ( const NetIn & netIn) const
private

Initializes the network with a short ignition phase.

Parameters
netInThe input conditions for the network.
Returns
The output conditions after the ignition phase.

This method performs a short integration of the network at a high temperature and density to ignite the network and bring it closer to equilibrium. This can improve the convergence of the QSE solver.

See also
DirectNetworkSolver::evaluate()

Definition at line 264 of file solver.cpp.

◆ packSpeciesTypeIndexVectors()

dynamicQSESpeciesIndices gridfire::solver::QSENetworkSolver::packSpeciesTypeIndexVectors ( const std::vector< double > & Y,
const double T9,
const double rho ) const
private

Packs the species indices into vectors based on their type (dynamic or QSE).

Parameters
YVector of current abundances for all species.
T9Temperature in units of 10^9 K.
rhoDensity in g/cm^3.
Returns
A dynamicQSESpeciesIndices struct containing the indices of the dynamic and QSE species.

This method determines whether each species should be treated dynamically or as being in QSE based on its timescale and abundance. Species with short timescales or low abundances are assumed to be in QSE.

See also
DynamicEngine::getSpeciesTimescales()

Definition at line 137 of file solver.cpp.

◆ shouldUpdateView()

bool gridfire::solver::QSENetworkSolver::shouldUpdateView ( const NetIn & conditions) const
private

Determines whether the adaptive engine view should be updated.

Parameters
conditionsThe current input conditions.
Returns
True if the view should be updated, false otherwise.

This method implements a policy for determining when the adaptive engine view should be updated. The view is updated if the temperature or density has changed significantly, or if a primary fuel source has been depleted.

See also
AdaptiveEngineView::update()

Definition at line 308 of file solver.cpp.

Member Data Documentation

◆ m_config

fourdst::config::Config& gridfire::solver::QSENetworkSolver::m_config = fourdst::config::Config::getInstance()
private

Configuration instance.

Definition at line 370 of file solver.h.

◆ m_isViewInitialized

bool gridfire::solver::QSENetworkSolver::m_isViewInitialized = false
private

Flag indicating whether the adaptive engine view has been initialized.

Definition at line 372 of file solver.h.

◆ m_lastSeenConditions

NetIn gridfire::solver::QSENetworkSolver::m_lastSeenConditions
private

The last seen input conditions.

Definition at line 373 of file solver.h.

◆ m_logger

quill::Logger* gridfire::solver::QSENetworkSolver::m_logger = fourdst::logging::LogManager::getInstance().getLogger("log")
private

Logger instance.

Definition at line 369 of file solver.h.


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