8#include "fourdst/logging/logging.h"
9#include "fourdst/config/config.h"
11#include "quill/Logger.h"
42 template <
typename EngineT>
104 using AdaptiveNetworkSolverStrategy::AdaptiveNetworkSolverStrategy;
140 const std::vector<double>& Y,
159 const std::vector<double>& Y,
222 const std::vector<size_t>& dynamicSpeciesIndices,
223 const std::vector<size_t>& QSESpeciesIndices,
224 const Eigen::VectorXd& Y_QSE,
242 const boost::numeric::ublas::vector<double>& YDynamic,
243 boost::numeric::ublas::vector<double>& dYdtDynamic,
276 const std::vector<size_t>& dynamicSpeciesIndices,
277 const std::vector<size_t>& QSESpeciesIndices,
295 const boost::numeric::ublas::vector<double>& YDynamic,
296 boost::numeric::ublas::matrix<double>& JDynamic,
298 boost::numeric::ublas::vector<double>& dfdt
332 const std::vector<double>& YDynamic,
333 const std::vector<size_t>& dynamicSpeciesIndices,
334 const std::vector<size_t>& QSESpeciesIndices,
362 quill::Logger*
m_logger = fourdst::logging::LogManager::getInstance().getLogger(
"log");
363 fourdst::config::Config&
m_config = fourdst::config::Config::getInstance();
385 using DynamicNetworkSolverStrategy::DynamicNetworkSolverStrategy;
431 const boost::numeric::ublas::vector<double>& Y,
432 boost::numeric::ublas::vector<double>& dYdt,
474 const boost::numeric::ublas::vector<double>& Y,
475 boost::numeric::ublas::matrix<double>& J,
477 boost::numeric::ublas::vector<double>& dfdt
483 quill::Logger*
m_logger = fourdst::logging::LogManager::getInstance().getLogger(
"log");
484 fourdst::config::Config&
m_config = fourdst::config::Config::getInstance();
489 std::vector<double> YFull(
m_engine.getNetworkSpecies().size(), 0.0);
490 Eigen::VectorXd Y_QSE(v_QSE.array().exp());
497 const auto [full_dYdt, specificEnergyGenerationRate] =
m_engine.calculateRHSAndEnergy(YFull,
m_T9,
m_rho);
506 template <
typename T>
508 std::vector<double> YFull(
m_engine.getNetworkSpecies().size(), 0.0);
509 Eigen::VectorXd Y_QSE(v_QSE.array().exp());
527 for (
long j = 0; j < J_QSE.cols(); ++j) {
528 J_QSE.col(j) *= Y_QSE(j);
Abstract class for engines supporting Jacobian and stoichiometry operations.
A network solver that directly integrates the reaction network ODEs.
quill::Logger * m_logger
Logger instance.
fourdst::config::Config & m_config
Configuration instance.
NetOut evaluate(const NetIn &netIn) override
Evaluates the network for a given timestep using direct integration.
Abstract base class for network solver strategies.
NetworkSolverStrategy(EngineT &engine)
Constructor for the NetworkSolverStrategy.
virtual ~NetworkSolverStrategy()=default
Virtual destructor.
virtual NetOut evaluate(const NetIn &netIn)=0
Evaluates the network for a given timestep.
A network solver that uses a Quasi-Steady-State Equilibrium (QSE) approach.
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.
bool shouldUpdateView(const NetIn &conditions) const
Determines whether the adaptive engine view should be updated.
NetIn m_lastSeenConditions
The last seen input conditions.
quill::Logger * m_logger
Logger instance.
NetOut evaluate(const NetIn &netIn) override
Evaluates the network for a given timestep using the QSE approach.
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).
fourdst::config::Config & m_config
Configuration instance.
bool m_isViewInitialized
Flag indicating whether the adaptive engine view has been initialized.
NetOut initializeNetworkWithShortIgnition(const NetIn &netIn) const
Initializes the network with a short ignition phase.
Abstract interfaces for reaction network engines in GridFire.
NetworkSolverStrategy< Engine > StaticNetworkSolverStrategy
Type alias for a network solver strategy that uses a static Engine.
NetworkSolverStrategy< DynamicEngine > DynamicNetworkSolverStrategy
Type alias for a network solver strategy that uses a DynamicEngine.
NetworkSolverStrategy< AdaptiveEngineView > AdaptiveNetworkSolverStrategy
Type alias for a network solver strategy that uses an AdaptiveEngineView.
const size_t m_numSpecies
The number of species in the network.
DynamicEngine & m_engine
The engine used to evaluate the network.
const double m_T9
Temperature in units of 10^9 K.
const double m_rho
Density in g/cm^3.
void operator()(const boost::numeric::ublas::vector< double > &Y, boost::numeric::ublas::matrix< double > &J, double t, boost::numeric::ublas::vector< double > &dfdt) const
Calculates the Jacobian matrix.
JacobianFunctor(DynamicEngine &engine, const double T9, const double rho)
Constructor for the JacobianFunctor.
DynamicEngine & m_engine
The engine used to evaluate the network.
const double m_T9
Temperature in units of 10^9 K.
void operator()(const boost::numeric::ublas::vector< double > &Y, boost::numeric::ublas::vector< double > &dYdt, double t) const
Calculates the time derivatives of the species abundances.
const double m_rho
Density in g/cm^3.
const size_t m_numSpecies
The number of species in the network.
RHSFunctor(DynamicEngine &engine, const double T9, const double rho)
Constructor for the RHSFunctor.
const std::vector< size_t > & m_dynamicSpeciesIndices
Indices of the dynamic species.
const double m_T9
Temperature in units of 10^9 K.
int operator()(const InputType &v_QSE, OutputType &f_QSE) const
Calculates the residual vector for the QSE species.
DynamicEngine & m_engine
The engine used to evaluate the network.
const double m_rho
Density in g/cm^3.
const std::vector< size_t > & m_QSESpeciesIndices
Indices of the QSE species.
Eigen::Matrix< T, Eigen::Dynamic, 1 > OutputType
int df(const InputType &v_QSE, JacobianType &J_QSE) const
Calculates the Jacobian matrix for the QSE species.
Eigen::Matrix< T, Eigen::Dynamic, 1 > InputType
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > JacobianType
const std::vector< double > & m_YDynamic
Abundances of the dynamic species.
EigenFunctor(DynamicEngine &engine, const std::vector< double > &YDynamic, const std::vector< size_t > &dynamicSpeciesIndices, const std::vector< size_t > &QSESpeciesIndices, const double T9, const double rho)
Constructor for the EigenFunctor.
const double m_rho
Density in g/cm^3.
const std::vector< size_t > & m_QSESpeciesIndices
Indices of the QSE species.
const std::vector< size_t > & m_dynamicSpeciesIndices
Indices of the dynamic species.
void operator()(const boost::numeric::ublas::vector< double > &YDynamic, boost::numeric::ublas::matrix< double > &JDynamic, double t, boost::numeric::ublas::vector< double > &dfdt) const
Calculates the Jacobian matrix of the ODEs for the dynamic species.
const double m_T9
Temperature in units of 10^9 K.
DynamicEngine & m_engine
The engine used to evaluate the network.
JacobianFunctor(DynamicEngine &engine, const std::vector< size_t > &dynamicSpeciesIndices, const std::vector< size_t > &QSESpeciesIndices, const double T9, const double rho)
Constructor for the JacobianFunctor.
const Eigen::VectorXd & m_Y_QSE
Steady-state abundances of the QSE species.
DynamicEngine & m_engine
The engine used to evaluate the network.
const double m_T9
Temperature in units of 10^9 K.
const std::vector< size_t > & m_dynamicSpeciesIndices
Indices of the dynamic species.
RHSFunctor(DynamicEngine &engine, const std::vector< size_t > &dynamicSpeciesIndices, const std::vector< size_t > &QSESpeciesIndices, const Eigen::VectorXd &Y_QSE, const double T9, const double rho)
Constructor for the RHSFunctor.
const std::vector< size_t > & m_QSESpeciesIndices
Indices of the QSE species.
const double m_rho
Density in g/cm^3.
void operator()(const boost::numeric::ublas::vector< double > &YDynamic, boost::numeric::ublas::vector< double > &dYdtDynamic, double t) const
Calculates the time derivatives of the dynamic species.
Structure to hold indices of dynamic and QSE species.
std::vector< size_t > QSESpeciesIndices
Indices of fast species that are in QSE.
std::vector< size_t > dynamicSpeciesIndices
Indices of slow species that are not in QSE.