GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
solver.h
Go to the documentation of this file.
1#pragma once
2
6#include "gridfire/network.h"
7
8#include "fourdst/logging/logging.h"
9#include "fourdst/config/config.h"
10
11#include "quill/Logger.h"
12
13#include "Eigen/Dense"
14
15#include <vector>
16
18
28 std::vector<size_t> dynamicSpeciesIndices;
29 std::vector<size_t> QSESpeciesIndices;
30 };
31
42 template <typename EngineT>
44 public:
49 explicit NetworkSolverStrategy(EngineT& engine) : m_engine(engine) {};
50
54 virtual ~NetworkSolverStrategy() = default;
55
61 virtual NetOut evaluate(const NetIn& netIn) = 0;
62 protected:
63 EngineT& m_engine;
64 };
65
70
75
80
99 public:
104 using AdaptiveNetworkSolverStrategy::AdaptiveNetworkSolverStrategy;
105
124 NetOut evaluate(const NetIn& netIn) override;
125 private: // methods
140 const std::vector<double>& Y,
141 const double T9,
142 const double rho
143 ) const;
144
158 Eigen::VectorXd calculateSteadyStateAbundances(
159 const std::vector<double>& Y,
160 const double T9,
161 const double rho,
162 const dynamicQSESpeciesIndices& indices
163 ) const;
164
177 const NetIn& netIn
178 ) const;
179
191 bool shouldUpdateView(const NetIn& conditions) const;
192 private: // Nested functors for ODE integration
201 struct RHSFunctor {
203 const std::vector<size_t>& m_dynamicSpeciesIndices;
204 const std::vector<size_t>& m_QSESpeciesIndices;
205 const Eigen::VectorXd& m_Y_QSE;
206 const double m_T9;
207 const double m_rho;
208
210
221 DynamicEngine& engine,
222 const std::vector<size_t>& dynamicSpeciesIndices,
223 const std::vector<size_t>& QSESpeciesIndices,
224 const Eigen::VectorXd& Y_QSE,
225 const double T9,
226 const double rho
227 ) :
228 m_engine(engine),
229 m_dynamicSpeciesIndices(dynamicSpeciesIndices),
230 m_QSESpeciesIndices(QSESpeciesIndices),
231 m_Y_QSE(Y_QSE),
232 m_T9(T9),
233 m_rho(rho) {}
234
241 void operator()(
242 const boost::numeric::ublas::vector<double>& YDynamic,
243 boost::numeric::ublas::vector<double>& dYdtDynamic,
244 double t
245 ) const;
246
247 };
248
261 const std::vector<size_t>& m_dynamicSpeciesIndices;
262 const std::vector<size_t>& m_QSESpeciesIndices;
263 const double m_T9;
264 const double m_rho;
265
275 DynamicEngine& engine,
276 const std::vector<size_t>& dynamicSpeciesIndices,
277 const std::vector<size_t>& QSESpeciesIndices,
278 const double T9,
279 const double rho
280 ) :
281 m_engine(engine),
282 m_dynamicSpeciesIndices(dynamicSpeciesIndices),
283 m_QSESpeciesIndices(QSESpeciesIndices),
284 m_T9(T9),
285 m_rho(rho) {}
286
295 const boost::numeric::ublas::vector<double>& YDynamic,
296 boost::numeric::ublas::matrix<double>& JDynamic,
297 double t,
298 boost::numeric::ublas::vector<double>& dfdt
299 ) const;
300 };
301
308 template<typename T>
310 using InputType = Eigen::Matrix<T, Eigen::Dynamic, 1>;
311 using OutputType = Eigen::Matrix<T, Eigen::Dynamic, 1>;
312 using JacobianType = Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic>;
313
315 const std::vector<double>& m_YDynamic;
316 const std::vector<size_t>& m_dynamicSpeciesIndices;
317 const std::vector<size_t>& m_QSESpeciesIndices;
318 const double m_T9;
319 const double m_rho;
320
331 DynamicEngine& engine,
332 const std::vector<double>& YDynamic,
333 const std::vector<size_t>& dynamicSpeciesIndices,
334 const std::vector<size_t>& QSESpeciesIndices,
335 const double T9,
336 const double rho
337 ) :
338 m_engine(engine),
339 m_YDynamic(YDynamic),
340 m_dynamicSpeciesIndices(dynamicSpeciesIndices),
341 m_QSESpeciesIndices(QSESpeciesIndices),
342 m_T9(T9),
343 m_rho(rho) {}
344
351 int operator()(const InputType& v_QSE, OutputType& f_QSE) const;
352
359 int df(const InputType& v_QSE, JacobianType& J_QSE) const;
360 };
361 private:
362 quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
363 fourdst::config::Config& m_config = fourdst::config::Config::getInstance();
364
365 bool m_isViewInitialized = false;
367 };
368
380 public:
385 using DynamicNetworkSolverStrategy::DynamicNetworkSolverStrategy;
386
392 NetOut evaluate(const NetIn& netIn) override;
393 private:
402 struct RHSFunctor {
404 const double m_T9;
405 const double m_rho;
406 const size_t m_numSpecies;
407
415 DynamicEngine& engine,
416 const double T9,
417 const double rho
418 ) :
419 m_engine(engine),
420 m_T9(T9),
421 m_rho(rho),
422 m_numSpecies(engine.getNetworkSpecies().size()) {}
423
430 void operator()(
431 const boost::numeric::ublas::vector<double>& Y,
432 boost::numeric::ublas::vector<double>& dYdt,
433 double t
434 ) const;
435 };
436
446 const double m_T9;
447 const double m_rho;
448 const size_t m_numSpecies;
449
457 DynamicEngine& engine,
458 const double T9,
459 const double rho
460 ) :
461 m_engine(engine),
462 m_T9(T9),
463 m_rho(rho),
464 m_numSpecies(engine.getNetworkSpecies().size()) {}
465
473 void operator()(
474 const boost::numeric::ublas::vector<double>& Y,
475 boost::numeric::ublas::matrix<double>& J,
476 double t,
477 boost::numeric::ublas::vector<double>& dfdt
478 ) const;
479
480 };
481
482 private:
483 quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
484 fourdst::config::Config& m_config = fourdst::config::Config::getInstance();
485 };
486
487 template<typename T>
489 std::vector<double> YFull(m_engine.getNetworkSpecies().size(), 0.0);
490 Eigen::VectorXd Y_QSE(v_QSE.array().exp());
491 for (size_t i = 0; i < m_dynamicSpeciesIndices.size(); ++i) {
492 YFull[m_dynamicSpeciesIndices[i]] = m_YDynamic[i];
493 }
494 for (size_t i = 0; i < m_QSESpeciesIndices.size(); ++i) {
495 YFull[m_QSESpeciesIndices[i]] = Y_QSE(i);
496 }
497 const auto [full_dYdt, specificEnergyGenerationRate] = m_engine.calculateRHSAndEnergy(YFull, m_T9, m_rho);
498 f_QSE.resize(m_QSESpeciesIndices.size());
499 for (size_t i = 0; i < m_QSESpeciesIndices.size(); ++i) {
500 f_QSE(i) = full_dYdt[m_QSESpeciesIndices[i]];
501 }
502
503 return 0; // Success
504 }
505
506 template <typename T>
508 std::vector<double> YFull(m_engine.getNetworkSpecies().size(), 0.0);
509 Eigen::VectorXd Y_QSE(v_QSE.array().exp());
510 for (size_t i = 0; i < m_dynamicSpeciesIndices.size(); ++i) {
511 YFull[m_dynamicSpeciesIndices[i]] = m_YDynamic[i];
512 }
513 for (size_t i = 0; i < m_QSESpeciesIndices.size(); ++i) {
514 YFull[m_QSESpeciesIndices[i]] = Y_QSE(i);
515 }
516
517 m_engine.generateJacobianMatrix(YFull, m_T9, m_rho);
518
519 Eigen::MatrixXd J_orig(m_QSESpeciesIndices.size(), m_QSESpeciesIndices.size());
520 for (size_t i = 0; i < m_QSESpeciesIndices.size(); ++i) {
521 for (size_t j = 0; j < m_QSESpeciesIndices.size(); ++j) {
522 J_orig(i, j) = m_engine.getJacobianMatrixEntry(m_QSESpeciesIndices[i], m_QSESpeciesIndices[j]);
523 }
524 }
525
526 J_QSE = J_orig;
527 for (long j = 0; j < J_QSE.cols(); ++j) {
528 J_QSE.col(j) *= Y_QSE(j); // Chain rule for log space
529 }
530 return 0; // Success
531 }
532
533
534}
Abstract class for engines supporting Jacobian and stoichiometry operations.
A network solver that directly integrates the reaction network ODEs.
Definition solver.h:379
quill::Logger * m_logger
Logger instance.
Definition solver.h:483
fourdst::config::Config & m_config
Configuration instance.
Definition solver.h:484
NetOut evaluate(const NetIn &netIn) override
Evaluates the network for a given timestep using direct integration.
Definition solver.cpp:308
Abstract base class for network solver strategies.
Definition solver.h:43
NetworkSolverStrategy(EngineT &engine)
Constructor for the NetworkSolverStrategy.
Definition solver.h:49
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.
Definition solver.h:98
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.
Definition solver.cpp:187
bool shouldUpdateView(const NetIn &conditions) const
Determines whether the adaptive engine view should be updated.
Definition solver.cpp:241
NetIn m_lastSeenConditions
The last seen input conditions.
Definition solver.h:366
quill::Logger * m_logger
Logger instance.
Definition solver.h:362
NetOut evaluate(const NetIn &netIn) override
Evaluates the network for a given timestep using the QSE approach.
Definition solver.cpp:23
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).
Definition solver.cpp:119
fourdst::config::Config & m_config
Configuration instance.
Definition solver.h:363
bool m_isViewInitialized
Flag indicating whether the adaptive engine view has been initialized.
Definition solver.h:365
NetOut initializeNetworkWithShortIgnition(const NetIn &netIn) const
Initializes the network with a short ignition phase.
Definition solver.cpp:202
Abstract interfaces for reaction network engines in GridFire.
NetworkSolverStrategy< Engine > StaticNetworkSolverStrategy
Type alias for a network solver strategy that uses a static Engine.
Definition solver.h:79
NetworkSolverStrategy< DynamicEngine > DynamicNetworkSolverStrategy
Type alias for a network solver strategy that uses a DynamicEngine.
Definition solver.h:69
NetworkSolverStrategy< AdaptiveEngineView > AdaptiveNetworkSolverStrategy
Type alias for a network solver strategy that uses an AdaptiveEngineView.
Definition solver.h:74
const size_t m_numSpecies
The number of species in the network.
Definition solver.h:448
DynamicEngine & m_engine
The engine used to evaluate the network.
Definition solver.h:445
const double m_T9
Temperature in units of 10^9 K.
Definition solver.h:446
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.
Definition solver.cpp:386
JacobianFunctor(DynamicEngine &engine, const double T9, const double rho)
Constructor for the JacobianFunctor.
Definition solver.h:456
DynamicEngine & m_engine
The engine used to evaluate the network.
Definition solver.h:403
const double m_T9
Temperature in units of 10^9 K.
Definition solver.h:404
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.
Definition solver.cpp:374
const double m_rho
Density in g/cm^3.
Definition solver.h:405
const size_t m_numSpecies
The number of species in the network.
Definition solver.h:406
RHSFunctor(DynamicEngine &engine, const double T9, const double rho)
Constructor for the RHSFunctor.
Definition solver.h:414
const std::vector< size_t > & m_dynamicSpeciesIndices
Indices of the dynamic species.
Definition solver.h:316
const double m_T9
Temperature in units of 10^9 K.
Definition solver.h:318
int operator()(const InputType &v_QSE, OutputType &f_QSE) const
Calculates the residual vector for the QSE species.
Definition solver.h:488
DynamicEngine & m_engine
The engine used to evaluate the network.
Definition solver.h:314
const double m_rho
Density in g/cm^3.
Definition solver.h:319
const std::vector< size_t > & m_QSESpeciesIndices
Indices of the QSE species.
Definition solver.h:317
Eigen::Matrix< T, Eigen::Dynamic, 1 > OutputType
Definition solver.h:311
int df(const InputType &v_QSE, JacobianType &J_QSE) const
Calculates the Jacobian matrix for the QSE species.
Definition solver.h:507
Eigen::Matrix< T, Eigen::Dynamic, 1 > InputType
Definition solver.h:310
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > JacobianType
Definition solver.h:312
const std::vector< double > & m_YDynamic
Abundances of the dynamic species.
Definition solver.h:315
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.
Definition solver.h:330
const double m_rho
Density in g/cm^3.
Definition solver.h:264
const std::vector< size_t > & m_QSESpeciesIndices
Indices of the QSE species.
Definition solver.h:262
const std::vector< size_t > & m_dynamicSpeciesIndices
Indices of the dynamic species.
Definition solver.h:261
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.
Definition solver.h:263
DynamicEngine & m_engine
The engine used to evaluate the network.
Definition solver.h:260
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.
Definition solver.h:274
const Eigen::VectorXd & m_Y_QSE
Steady-state abundances of the QSE species.
Definition solver.h:205
DynamicEngine & m_engine
The engine used to evaluate the network.
Definition solver.h:202
const double m_T9
Temperature in units of 10^9 K.
Definition solver.h:206
const std::vector< size_t > & m_dynamicSpeciesIndices
Indices of the dynamic species.
Definition solver.h:203
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.
Definition solver.h:220
const std::vector< size_t > & m_QSESpeciesIndices
Indices of the QSE species.
Definition solver.h:204
const double m_rho
Density in g/cm^3.
Definition solver.h:207
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.
Definition solver.cpp:283
Structure to hold indices of dynamic and QSE species.
Definition solver.h:27
std::vector< size_t > QSESpeciesIndices
Indices of fast species that are in QSE.
Definition solver.h:29
std::vector< size_t > dynamicSpeciesIndices
Indices of slow species that are not in QSE.
Definition solver.h:28