GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
engine_adaptive.h
Go to the documentation of this file.
1#pragma once
6#include "gridfire/network.h"
7
8#include "fourdst/composition/atomicSpecies.h"
9#include "fourdst/config/config.h"
10#include "fourdst/logging/logging.h"
11
12#include "quill/Logger.h"
13
14namespace gridfire {
47 class AdaptiveEngineView final : public DynamicEngine, public EngineView<DynamicEngine> {
48 public:
57 explicit AdaptiveEngineView(DynamicEngine& baseEngine);
58
76 void update(const NetIn& netIn) override;
77
82 [[nodiscard]] const std::vector<fourdst::atomic::Species>& getNetworkSpecies() const override;
83
101 const std::vector<double> &Y_culled,
102 const double T9,
103 const double rho
104 ) const override;
105
120 const std::vector<double> &Y_culled,
121 const double T9,
122 const double rho
123 ) override;
124
139 [[nodiscard]] double getJacobianMatrixEntry(
140 const int i_culled,
141 const int j_culled
142 ) const override;
143
153 void generateStoichiometryMatrix() override;
154
169 [[nodiscard]] int getStoichiometryMatrixEntry(
170 const int speciesIndex_culled,
171 const int reactionIndex_culled
172 ) const override;
173
189 [[nodiscard]] double calculateMolarReactionFlow(
191 const std::vector<double> &Y_culled,
192 double T9,
193 double rho
194 ) const override;
195
201 [[nodiscard]] const reaction::LogicalReactionSet& getNetworkReactions() const override;
202
216 [[nodiscard]] std::unordered_map<fourdst::atomic::Species, double> getSpeciesTimescales(
217 const std::vector<double> &Y_culled,
218 double T9,
219 double rho
220 ) const override;
221
226 [[nodiscard]] const DynamicEngine& getBaseEngine() const override { return m_baseEngine; }
227
243 void setScreeningModel(screening::ScreeningType model) override;
244
258 [[nodiscard]] screening::ScreeningType getScreeningModel() const override;
259 private:
260 using Config = fourdst::config::Config;
261 using LogManager = fourdst::logging::LogManager;
263 Config& m_config = Config::getInstance();
265 quill::Logger* m_logger = LogManager::getInstance().getLogger("log");
266
269
271 std::vector<fourdst::atomic::Species> m_activeSpecies;
274
276 std::vector<size_t> m_speciesIndexMap;
278 std::vector<size_t> m_reactionIndexMap;
279
281 bool m_isStale = true;
282
283 private:
291 private:
302 [[nodiscard]] std::vector<size_t> constructSpeciesIndexMap() const;
303
314 [[nodiscard]] std::vector<size_t> constructReactionIndexMap() const;
315
323 [[nodiscard]] std::vector<double> mapCulledToFull(const std::vector<double>& culled) const;
324
332 [[nodiscard]] std::vector<double> mapFullToCulled(const std::vector<double>& full) const;
333
342 [[nodiscard]] size_t mapCulledToFullSpeciesIndex(size_t culledSpeciesIndex) const;
343
352 [[nodiscard]] size_t mapCulledToFullReactionIndex(size_t culledReactionIndex) const;
353
359 void validateState() const;
360
382 std::vector<ReactionFlow> calculateAllReactionFlows(
383 const NetIn& netIn,
384 std::vector<double>& out_Y_Full
385 ) const;
403 [[nodiscard]] std::unordered_set<fourdst::atomic::Species> findReachableSpecies(
404 const NetIn& netIn
405 ) const;
426 [[nodiscard]] std::vector<const reaction::LogicalReaction*> cullReactionsByFlow(
427 const std::vector<ReactionFlow>& allFlows,
428 const std::unordered_set<fourdst::atomic::Species>& reachableSpecies,
429 const std::vector<double>& Y_full,
430 double maxFlow
431 ) const;
448 const std::vector<const reaction::LogicalReaction*>& finalReactions
449 );
450 };
451}
double calculateMolarReactionFlow(const reaction::Reaction &reaction, const std::vector< double > &Y_culled, double T9, double rho) const override
Calculates the molar reaction flow for a given reaction in the active network.
screening::ScreeningType getScreeningModel() const override
Gets the screening model from the base engine.
std::unordered_set< fourdst::atomic::Species > findReachableSpecies(const NetIn &netIn) const
Finds all species that are reachable from the initial fuel through the reaction network.
const reaction::LogicalReactionSet & getNetworkReactions() const override
Gets the set of active logical reactions in the network.
Config & m_config
A reference to the singleton Config instance, used for retrieving configuration parameters.
reaction::LogicalReactionSet m_activeReactions
The set of reactions that are currently active in the network.
std::vector< size_t > m_reactionIndexMap
A map from the indices of the active reactions to the indices of the corresponding reactions in the f...
void generateStoichiometryMatrix() override
Generates the stoichiometry matrix for the active reactions and species.
size_t mapCulledToFullSpeciesIndex(size_t culledSpeciesIndex) const
Maps a culled species index to a full species index.
std::vector< double > mapFullToCulled(const std::vector< double > &full) const
Maps a vector of full abundances to a vector of culled abundances.
std::vector< const reaction::LogicalReaction * > cullReactionsByFlow(const std::vector< ReactionFlow > &allFlows, const std::unordered_set< fourdst::atomic::Species > &reachableSpecies, const std::vector< double > &Y_full, double maxFlow) const
Culls reactions from the network based on their flow rates.
double getJacobianMatrixEntry(const int i_culled, const int j_culled) const override
Gets an entry from the Jacobian matrix for the active species.
DynamicEngine & m_baseEngine
The underlying engine to which this view delegates calculations.
fourdst::logging::LogManager LogManager
std::vector< size_t > m_speciesIndexMap
A map from the indices of the active species to the indices of the corresponding species in the full ...
bool m_isStale
A flag indicating whether the view is stale and needs to be updated.
int getStoichiometryMatrixEntry(const int speciesIndex_culled, const int reactionIndex_culled) const override
Gets an entry from the stoichiometry matrix for the active species and reactions.
std::vector< double > mapCulledToFull(const std::vector< double > &culled) const
Maps a vector of culled abundances to a vector of full abundances.
StepDerivatives< double > calculateRHSAndEnergy(const std::vector< double > &Y_culled, const double T9, const double rho) const override
Calculates the right-hand side (dY/dt) and energy generation for the active species.
void update(const NetIn &netIn) override
Updates the active species and reactions based on the current conditions.
std::vector< size_t > constructReactionIndexMap() const
Constructs the reaction index map.
std::vector< size_t > constructSpeciesIndexMap() const
Constructs the species index map.
size_t mapCulledToFullReactionIndex(size_t culledReactionIndex) const
Maps a culled reaction index to a full reaction index.
std::unordered_map< fourdst::atomic::Species, double > getSpeciesTimescales(const std::vector< double > &Y_culled, double T9, double rho) const override
Computes timescales for all active species in the network.
void finalizeActiveSet(const std::vector< const reaction::LogicalReaction * > &finalReactions)
Finalizes the set of active species and reactions.
void setScreeningModel(screening::ScreeningType model) override
Sets the screening model for the base engine.
std::vector< ReactionFlow > calculateAllReactionFlows(const NetIn &netIn, std::vector< double > &out_Y_Full) const
Calculates the molar reaction flow rate for all reactions in the full network.
quill::Logger * m_logger
A pointer to the logger instance, used for logging messages.
const std::vector< fourdst::atomic::Species > & getNetworkSpecies() const override
Gets the list of active species in the network.
void generateJacobianMatrix(const std::vector< double > &Y_culled, const double T9, const double rho) override
Generates the Jacobian matrix for the active species.
AdaptiveEngineView(DynamicEngine &baseEngine)
Constructs an AdaptiveEngineView.
void validateState() const
Validates that the AdaptiveEngineView is not stale.
const DynamicEngine & getBaseEngine() const override
Gets the base engine.
std::vector< fourdst::atomic::Species > m_activeSpecies
The set of species that are currently active in the network.
fourdst::config::Config Config
Abstract class for engines supporting Jacobian and stoichiometry operations.
Abstract base class for a "view" of a reaction network engine.
Represents a "logical" reaction that aggregates rates from multiple sources.
Definition reaction.h:308
Represents a single nuclear reaction from a specific data source.
Definition reaction.h:72
Abstract interfaces for reaction network engines in GridFire.
Abstract interfaces for engine "views" in GridFire.
TemplatedReactionSet< LogicalReaction > LogicalReactionSet
A set of logical reactions.
Definition reaction.h:557
ScreeningType
Enumerates the available plasma screening models.
A struct to hold a reaction and its flow rate.
const reaction::LogicalReaction * reactionPtr
Structure holding derivatives and energy generation for a network step.