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
5#include <vector>
6#include <unordered_map>
7
20
21namespace gridfire {
22
29 template<typename T>
30 concept IsArithmeticOrAD = std::is_same_v<T, double> || std::is_same_v<T, CppAD::AD<double>>;
31
49 template <IsArithmeticOrAD T>
51 std::vector<T> dydt;
53 };
54
72 class Engine {
73 public:
77 virtual ~Engine() = default;
78
83 virtual const std::vector<fourdst::atomic::Species>& getNetworkSpecies() const = 0;
84
98 const std::vector<double>& Y,
99 double T9,
100 double rho
101 ) const = 0;
102 };
103
118 class DynamicEngine : public Engine {
119 public:
131 const std::vector<double>& Y,
132 double T9, double rho
133 ) = 0;
134
145 int i,
146 int j
147 ) const = 0;
148
155 virtual void generateStoichiometryMatrix() = 0;
156
167 int speciesIndex,
168 int reactionIndex
169 ) const = 0;
170
185 const std::vector<double>& Y,
186 double T9,
187 double rho
188 ) const = 0;
189
196
208 virtual std::unordered_map<fourdst::atomic::Species, double> getSpeciesTimescales(
209 const std::vector<double>& Y,
210 double T9,
211 double rho
212 ) const = 0;
213 };
214}
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 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 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.
A collection of LogicalReaction objects.
Definition reaction.h:554
Represents a single nuclear reaction from a specific data source.
Definition reaction.h:71
Concept for types allowed in engine calculations.
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).