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
8#include <vector>
9#include <unordered_map>
10
23
24namespace gridfire {
25
32 template<typename T>
33 concept IsArithmeticOrAD = std::is_same_v<T, double> || std::is_same_v<T, CppAD::AD<double>>;
34
52 template <IsArithmeticOrAD T>
54 std::vector<T> dydt;
56 };
57
75 class Engine {
76 public:
80 virtual ~Engine() = default;
81
86 [[nodiscard]] virtual const std::vector<fourdst::atomic::Species>& getNetworkSpecies() const = 0;
87
101 const std::vector<double>& Y,
102 double T9,
103 double rho
104 ) const = 0;
105 };
106
121 class DynamicEngine : public Engine {
122 public:
134 const std::vector<double>& Y,
135 double T9, double rho
136 ) = 0;
137
147 [[nodiscard]] virtual double getJacobianMatrixEntry(
148 int i,
149 int j
150 ) const = 0;
151
158 virtual void generateStoichiometryMatrix() = 0;
159
169 [[nodiscard]] virtual int getStoichiometryMatrixEntry(
170 int speciesIndex,
171 int reactionIndex
172 ) const = 0;
173
186 [[nodiscard]] virtual double calculateMolarReactionFlow(
188 const std::vector<double>& Y,
189 double T9,
190 double rho
191 ) const = 0;
192
198 [[nodiscard]] virtual const reaction::LogicalReactionSet& getNetworkReactions() const = 0;
199
211 [[nodiscard]] virtual std::unordered_map<fourdst::atomic::Species, double> getSpeciesTimescales(
212 const std::vector<double>& Y,
213 double T9,
214 double rho
215 ) const = 0;
216
237 virtual void update(const NetIn& netIn) = 0;
238
256
267 [[nodiscard]] virtual screening::ScreeningType getScreeningModel() const = 0;
268 };
269}
Abstract class for engines supporting Jacobian and stoichiometry operations.
virtual double getJacobianMatrixEntry(int i, int j) const =0
Get an entry from the previously generated Jacobian matrix.
virtual void generateJacobianMatrix(const std::vector< double > &Y, double T9, double rho)=0
Generate the Jacobian matrix for the current state.
virtual void setScreeningModel(screening::ScreeningType model)=0
Set the electron screening model.
virtual std::unordered_map< fourdst::atomic::Species, double > getSpeciesTimescales(const std::vector< double > &Y, double T9, double rho) const =0
Compute timescales for all species in the network.
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 screening::ScreeningType getScreeningModel() const =0
Get the current electron screening model.
virtual void 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 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.
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 StepDerivatives< double > 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:557
ScreeningType
Enumerates the available plasma screening models.
Defines classes for representing and managing nuclear reactions.
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).