GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
engine_abstract.h
Go to the documentation of this file.
1#pragma once
2
4#include "gridfire/network.h"
7
10
12
13#include <vector>
14#include <unordered_map>
15#include <utility>
16#include <expected>
17
30
31namespace gridfire {
32
39 template<typename T>
40 concept IsArithmeticOrAD = std::is_same_v<T, double> || std::is_same_v<T, CppAD::AD<double>>;
41
59 template <IsArithmeticOrAD T>
61 std::vector<T> dydt;
63 };
64
65 using SparsityPattern = std::vector<std::pair<size_t, size_t>>;
66
84 class Engine {
85 public:
89 virtual ~Engine() = default;
90
95 [[nodiscard]] virtual const std::vector<fourdst::atomic::Species>& getNetworkSpecies() const = 0;
96
109 [[nodiscard]] virtual std::expected<StepDerivatives<double>, expectations::StaleEngineError> calculateRHSAndEnergy(
110 const std::vector<double>& Y,
111 double T9,
112 double rho
113 ) const = 0;
114 };
115
130 class DynamicEngine : public Engine {
131 public:
143 const std::vector<double>& Y_dynamic,
144 double T9,
145 double rho
146 ) const = 0;
147
149 const std::vector<double>& Y_dynamic,
150 double T9,
151 double rho,
152 const SparsityPattern& sparsityPattern
153 ) const {
154 throw std::logic_error("Sparsity pattern not supported by this engine.");
155 }
156
166 [[nodiscard]] virtual double getJacobianMatrixEntry(
167 int i,
168 int j
169 ) const = 0;
170
171
178 virtual void generateStoichiometryMatrix() = 0;
179
189 [[nodiscard]] virtual int getStoichiometryMatrixEntry(
190 int speciesIndex,
191 int reactionIndex
192 ) const = 0;
193
206 [[nodiscard]] virtual double calculateMolarReactionFlow(
208 const std::vector<double>& Y,
209 double T9,
210 double rho
211 ) const = 0;
212
218 [[nodiscard]] virtual const reaction::LogicalReactionSet& getNetworkReactions() const = 0;
219
220 virtual void setNetworkReactions(const reaction::LogicalReactionSet& reactions) = 0;
221
233 [[nodiscard]] virtual std::expected<std::unordered_map<fourdst::atomic::Species, double>, expectations::StaleEngineError> getSpeciesTimescales(
234 const std::vector<double>& Y,
235 double T9,
236 double rho
237 ) const = 0;
238
239 [[nodiscard]] virtual std::expected<std::unordered_map<fourdst::atomic::Species, double>, expectations::StaleEngineError> getSpeciesDestructionTimescales(
240 const std::vector<double>& Y,
241 double T9,
242 double rho
243 ) const = 0;
244
265 virtual fourdst::composition::Composition update(const NetIn &netIn) = 0;
266
267 virtual bool isStale(const NetIn& netIn) = 0;
268
286
297 [[nodiscard]] virtual screening::ScreeningType getScreeningModel() const = 0;
298
299 [[nodiscard]] virtual int getSpeciesIndex(const fourdst::atomic::Species &species) const = 0;
300
301 [[nodiscard]] virtual std::vector<double> mapNetInToMolarAbundanceVector(const NetIn &netIn) const = 0;
302
303 [[nodiscard]] virtual PrimingReport primeEngine(const NetIn &netIn) = 0;
304
305 [[nodiscard]] virtual BuildDepthType getDepth() const {
306 throw std::logic_error("Network depth not supported by this engine.");
307 }
308
309 virtual void rebuild(const fourdst::composition::Composition& comp, BuildDepthType depth) {
310 throw std::logic_error("Setting network depth not supported by this engine.");
311 }
312
313 };
314}
Abstract class for engines supporting Jacobian and stoichiometry operations.
virtual BuildDepthType getDepth() const
virtual double getJacobianMatrixEntry(int i, int j) const =0
Get an entry from the previously generated Jacobian matrix.
virtual PrimingReport primeEngine(const NetIn &netIn)=0
virtual void generateJacobianMatrix(const std::vector< double > &Y_dynamic, double T9, double rho) const =0
Generate the Jacobian matrix for the current state.
virtual void setScreeningModel(screening::ScreeningType model)=0
Set the electron screening model.
virtual void rebuild(const fourdst::composition::Composition &comp, BuildDepthType depth)
virtual std::vector< double > mapNetInToMolarAbundanceVector(const NetIn &netIn) const =0
virtual double calculateMolarReactionFlow(const reaction::Reaction &reaction, const std::vector< double > &Y, double T9, double rho) const =0
Calculate the molar reaction flow for a given reaction.
virtual std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > getSpeciesTimescales(const std::vector< double > &Y, double T9, double rho) const =0
Compute timescales for all species in the network.
virtual std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > getSpeciesDestructionTimescales(const std::vector< double > &Y, double T9, double rho) const =0
virtual screening::ScreeningType getScreeningModel() const =0
Get the current electron screening model.
virtual void generateJacobianMatrix(const std::vector< double > &Y_dynamic, double T9, double rho, const SparsityPattern &sparsityPattern) const
virtual bool isStale(const NetIn &netIn)=0
virtual fourdst::composition::Composition update(const NetIn &netIn)=0
Update the internal state of the engine.
virtual const reaction::LogicalReactionSet & getNetworkReactions() const =0
Get the set of logical reactions in the network.
virtual int getSpeciesIndex(const fourdst::atomic::Species &species) const =0
virtual void generateStoichiometryMatrix()=0
Generate the stoichiometry matrix for the network.
virtual int getStoichiometryMatrixEntry(int speciesIndex, int reactionIndex) const =0
Get an entry from the stoichiometry matrix.
virtual void setNetworkReactions(const reaction::LogicalReactionSet &reactions)=0
Abstract base class for a reaction network engine.
virtual const std::vector< fourdst::atomic::Species > & getNetworkSpecies() const =0
Get the list of species in the network.
virtual ~Engine()=default
Virtual destructor.
virtual std::expected< StepDerivatives< double >, expectations::StaleEngineError > calculateRHSAndEnergy(const std::vector< double > &Y, double T9, double rho) const =0
Calculate the right-hand side (dY/dt) and energy generation.
Represents a single nuclear reaction from a specific data source.
Definition reaction.h:72
Concept for types allowed in engine calculations.
TemplatedReactionSet< LogicalReaction > LogicalReactionSet
A set of logical reactions.
Definition reaction.h:563
ScreeningType
Enumerates the available plasma screening models.
std::variant< NetworkBuildDepth, int > BuildDepthType
Variant specifying either a predefined NetworkBuildDepth or a custom integer depth.
Definition building.h:37
std::vector< std::pair< size_t, size_t > > SparsityPattern
Defines classes for representing and managing nuclear reactions.
Captures the result of a network priming operation.
Definition reporting.h:67
Structure holding derivatives and energy generation for a network step.
T nuclearEnergyGenerationRate
Specific energy generation rate (e.g., erg/g/s).
std::vector< T > dydt
Derivatives of abundances (dY/dt for each species).