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 "unsupported/Eigen/NonLinearOptimization" // Required for LevenbergMarquardt
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 DynamicNetworkSolverStrategy::DynamicNetworkSolverStrategy;
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 enum {
314 InputsAtCompileTime = Eigen::Dynamic,
315 ValuesAtCompileTime = Eigen::Dynamic
316 };
317
319 const std::vector<double>& m_YFull;
320 const std::vector<size_t>& m_dynamicSpeciesIndices;
321 const std::vector<size_t>& m_QSESpeciesIndices;
322 const double m_T9;
323 const double m_rho;
324
335 DynamicEngine& engine,
336 const std::vector<double>& YFull,
337 const std::vector<size_t>& dynamicSpeciesIndices,
338 const std::vector<size_t>& QSESpeciesIndices,
339 const double T9,
340 const double rho
341 ) :
342 m_engine(engine),
343 m_YFull(YFull),
344 m_dynamicSpeciesIndices(dynamicSpeciesIndices),
345 m_QSESpeciesIndices(QSESpeciesIndices),
346 m_T9(T9),
347 m_rho(rho) {}
348
349 int values() const { return m_QSESpeciesIndices.size(); }
350 int inputs() const { return m_QSESpeciesIndices.size(); }
351
358 int operator()(const InputType& v_QSE_log, OutputType& f_QSE) const;
359
366 int df(const InputType& v_QSE_log, JacobianType& J_QSE) const;
367 };
368 private:
369 quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
370 fourdst::config::Config& m_config = fourdst::config::Config::getInstance();
371
372 bool m_isViewInitialized = false;
374 };
375
387 public:
392 using DynamicNetworkSolverStrategy::DynamicNetworkSolverStrategy;
393
399 NetOut evaluate(const NetIn& netIn) override;
400 private:
409 struct RHSFunctor {
411 const double m_T9;
412 const double m_rho;
413 const size_t m_numSpecies;
414 quill::Logger* m_logger = LogManager::getInstance().getLogger("log");
415
423 DynamicEngine& engine,
424 const double T9,
425 const double rho
426 ) :
427 m_engine(engine),
428 m_T9(T9),
429 m_rho(rho),
430 m_numSpecies(engine.getNetworkSpecies().size()) {}
431
438 void operator()(
439 const boost::numeric::ublas::vector<double>& Y,
440 boost::numeric::ublas::vector<double>& dYdt,
441 double t
442 ) const;
443 };
444
454 const double m_T9;
455 const double m_rho;
456 const size_t m_numSpecies;
457
465 DynamicEngine& engine,
466 const double T9,
467 const double rho
468 ) :
469 m_engine(engine),
470 m_T9(T9),
471 m_rho(rho),
472 m_numSpecies(engine.getNetworkSpecies().size()) {}
473
481 void operator()(
482 const boost::numeric::ublas::vector<double>& Y,
483 boost::numeric::ublas::matrix<double>& J,
484 double t,
485 boost::numeric::ublas::vector<double>& dfdt
486 ) const;
487
488 };
489
490 private:
491 quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
492 fourdst::config::Config& m_config = fourdst::config::Config::getInstance();
493 };
494
495 template<typename T>
497 std::vector<double> y = m_YFull; // Full vector of species abundances
498 Eigen::VectorXd Y_QSE = v_QSE_log.array().exp();
499
500 for (size_t i = 0; i < m_QSESpeciesIndices.size(); ++i) {
501 y[m_QSESpeciesIndices[i]] = Y_QSE(i);
502 }
503
504
505 const auto [dydt, specificEnergyRate] = m_engine.calculateRHSAndEnergy(y, m_T9, m_rho);
506 f_QSE.resize(m_QSESpeciesIndices.size());
507 for (size_t i = 0; i < m_QSESpeciesIndices.size(); ++i) {
508 f_QSE(i) = dydt[m_QSESpeciesIndices[i]];
509 }
510 return 0; // Success
511 }
512
513 template <typename T>
515 std::vector<double> y = m_YFull;
516 Eigen::VectorXd Y_QSE = v_QSE_log.array().exp();
517
518 for (size_t i = 0; i < m_QSESpeciesIndices.size(); ++i) {
519 y[m_QSESpeciesIndices[i]] = Y_QSE(i);
520 }
521
522 m_engine.generateJacobianMatrix(y, m_T9, m_rho);
523
524 J_QSE.resize(m_QSESpeciesIndices.size(), m_QSESpeciesIndices.size());
525 for (size_t i = 0; i < m_QSESpeciesIndices.size(); ++i) {
526 for (size_t j = 0; j < m_QSESpeciesIndices.size(); ++j) {
527 J_QSE(i, j) = m_engine.getJacobianMatrixEntry(m_QSESpeciesIndices[i], m_QSESpeciesIndices[j]);
528 }
529 }
530
531 for (long j = 0; j < J_QSE.cols(); ++j) {
532 J_QSE(j, j) *= Y_QSE(j); // chain rule for log space transformation
533 }
534 return 0;
535 }
536
537
538}
Abstract class for engines supporting Jacobian and stoichiometry operations.
A network solver that directly integrates the reaction network ODEs.
Definition solver.h:386
quill::Logger * m_logger
Logger instance.
Definition solver.h:491
fourdst::config::Config & m_config
Configuration instance.
Definition solver.h:492
NetOut evaluate(const NetIn &netIn) override
Evaluates the network for a given timestep using direct integration.
Definition solver.cpp:375
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:205
bool shouldUpdateView(const NetIn &conditions) const
Determines whether the adaptive engine view should be updated.
Definition solver.cpp:308
NetIn m_lastSeenConditions
The last seen input conditions.
Definition solver.h:373
quill::Logger * m_logger
Logger instance.
Definition solver.h:369
NetOut evaluate(const NetIn &netIn) override
Evaluates the network for a given timestep using the QSE approach.
Definition solver.cpp:26
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:137
fourdst::config::Config & m_config
Configuration instance.
Definition solver.h:370
bool m_isViewInitialized
Flag indicating whether the adaptive engine view has been initialized.
Definition solver.h:372
NetOut initializeNetworkWithShortIgnition(const NetIn &netIn) const
Initializes the network with a short ignition phase.
Definition solver.cpp:264
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:456
DynamicEngine & m_engine
The engine used to evaluate the network.
Definition solver.h:453
const double m_T9
Temperature in units of 10^9 K.
Definition solver.h:454
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:463
JacobianFunctor(DynamicEngine &engine, const double T9, const double rho)
Constructor for the JacobianFunctor.
Definition solver.h:464
DynamicEngine & m_engine
The engine used to evaluate the network.
Definition solver.h:410
const double m_T9
Temperature in units of 10^9 K.
Definition solver.h:411
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:442
const double m_rho
Density in g/cm^3.
Definition solver.h:412
quill::Logger * m_logger
Logger instance.
Definition solver.h:414
const size_t m_numSpecies
The number of species in the network.
Definition solver.h:413
RHSFunctor(DynamicEngine &engine, const double T9, const double rho)
Constructor for the RHSFunctor.
Definition solver.h:422
const std::vector< size_t > & m_dynamicSpeciesIndices
Indices of the dynamic species.
Definition solver.h:320
int operator()(const InputType &v_QSE_log, OutputType &f_QSE) const
Calculates the residual vector for the QSE species.
Definition solver.h:496
const std::vector< double > & m_YFull
The full, initial abundance vector.
Definition solver.h:319
int df(const InputType &v_QSE_log, JacobianType &J_QSE) const
Calculates the Jacobian matrix for the QSE species.
Definition solver.h:514
const double m_T9
Temperature in units of 10^9 K.
Definition solver.h:322
DynamicEngine & m_engine
The engine used to evaluate the network.
Definition solver.h:318
EigenFunctor(DynamicEngine &engine, const std::vector< double > &YFull, 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:334
const double m_rho
Density in g/cm^3.
Definition solver.h:323
const std::vector< size_t > & m_QSESpeciesIndices
Indices of the QSE species.
Definition solver.h:321
Eigen::Matrix< T, Eigen::Dynamic, 1 > OutputType
Definition solver.h:311
Eigen::Matrix< T, Eigen::Dynamic, 1 > InputType
Definition solver.h:310
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > JacobianType
Definition solver.h:312
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:350
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