GridFire v0.7.0_rc2
General Purpose Nuclear Network
Loading...
Searching...
No Matches
gridfire::reaction::ReaclibReaction Class Reference

#include <reaction.h>

Inheritance diagram for gridfire::reaction::ReaclibReaction:
[legend]
Collaboration diagram for gridfire::reaction::ReaclibReaction:
[legend]

Public Member Functions

 ~ReaclibReaction () override=default
 
 ReaclibReaction (std::string_view id, std::string_view peName, int chapter, const std::vector< fourdst::atomic::Species > &reactants, const std::vector< fourdst::atomic::Species > &products, double qValue, std::string_view label, const RateCoefficientSet &sets, bool reverse=false)
 Constructs a Reaction object.
 
double calculate_rate (double T9, double rho, double Ye, double mue, const std::vector< double > &Y, const std::unordered_map< size_t, fourdst::atomic::Species > &index_to_species_map) const override
 Calculates the reaction rate for a given temperature.
 
CppAD::AD< double > calculate_rate (CppAD::AD< double > T9, CppAD::AD< double > rho, CppAD::AD< double > Ye, CppAD::AD< double > mue, const std::vector< CppAD::AD< double > > &Y, const std::unordered_map< size_t, fourdst::atomic::Species > &index_to_species_map) const override
 Calculates the reaction rate for a given temperature using CppAD types.
 
double calculate_log_rate_partial_deriv_wrt_T9 (double T9, double rho, double Ye, double mue, const fourdst::composition::Composition &comp) const override
 Logarithmic partial derivative of the rate with respect to temperature.
 
virtual std::string_view peName () const
 Gets the reaction name in (projectile, ejectile) notation.
 
int chapter () const
 Gets the REACLIB chapter number.
 
std::string_view sourceLabel () const
 Gets the source label for the rate data.
 
ReactionType type () const override
 Category of this reaction (e.g., REACLIB, WEAK, LOGICAL_REACLIB).
 
const RateCoefficientSetrateCoefficients () const
 Gets the set of rate coefficients.
 
std::optional< std::vector< RateCoefficientSet > > getRateCoefficients () const override
 
bool contains (const fourdst::atomic::Species &species) const override
 Checks if the reaction involves a given species as a reactant or product.
 
bool contains_reactant (const fourdst::atomic::Species &species) const override
 Checks if the reaction involves a given species as a reactant.
 
bool contains_product (const fourdst::atomic::Species &species) const override
 Checks if the reaction involves a given species as a product.
 
std::unordered_set< fourdst::atomic::Species > all_species () const override
 Gets a set of all unique species involved in the reaction.
 
std::unordered_set< fourdst::atomic::Species > reactant_species () const override
 Gets a set of all unique reactant species.
 
std::unordered_set< fourdst::atomic::Species > product_species () const override
 Gets a set of all unique product species.
 
size_t num_species () const override
 Gets the number of unique species involved in the reaction.
 
int stoichiometry (const fourdst::atomic::Species &species) const override
 Calculates the stoichiometric coefficient for a given species.
 
std::unordered_map< fourdst::atomic::Species, int > stoichiometry () const override
 Gets a map of all species to their stoichiometric coefficients.
 
std::string_view id () const override
 Gets the unique identifier of the reaction.
 
double qValue () const override
 Gets the Q-value of the reaction.
 
const std::vector< fourdst::atomic::Species > & reactants () const override
 Gets the vector of reactant species.
 
const std::vector< fourdst::atomic::Species > & products () const override
 Gets the vector of product species.
 
bool is_reverse () const override
 Checks if this is a reverse reaction rate.
 
double excess_energy () const
 Calculates the excess energy from the mass difference of reactants and products.
 
bool operator== (const ReaclibReaction &other) const
 Compares this reaction with another for equality based on their IDs.
 
bool operator!= (const ReaclibReaction &other) const
 Compares this reaction with another for inequality.
 
uint64_t hash (uint64_t seed) const override
 Computes a hash for the reaction based on its ID.
 
std::unique_ptr< Reactionclone () const override
 Polymorphic deep copy.
 
size_t countReactantOccurrences (const fourdst::atomic::Species &species) const override
 
size_t countProductOccurrences (const fourdst::atomic::Species &species) const override
 
- Public Member Functions inherited from gridfire::reaction::Reaction
virtual ~Reaction ()=default
 Virtual destructor for correct polymorphic cleanup.
 
virtual double calculate_energy_generation_rate (const double T9, const double rho, const double Ye, double mue, const std::vector< double > &Y, const std::unordered_map< size_t, fourdst::atomic::Species > &index_to_species_map) const
 Convenience: energy generation rate from this reaction (double version).
 
virtual CppAD::AD< double > calculate_energy_generation_rate (const CppAD::AD< double > &T9, const CppAD::AD< double > &rho, const CppAD::AD< double > &Ye, const CppAD::AD< double > &mue, const std::vector< CppAD::AD< double > > &Y, const std::unordered_map< size_t, fourdst::atomic::Species > &index_to_species_map) const
 Convenience: AD-enabled energy generation rate (AD version).
 

Protected Attributes

quill::Logger * m_logger = fourdst::logging::LogManager::getInstance().getLogger("log")
 
std::string m_id
 Unique identifier for the reaction (e.g., "h1+h1=>h2+e+nu").
 
std::string m_peName
 Name of the reaction in (projectile, ejectile) notation (e.g. "p(p,g)d").
 
int m_chapter
 Chapter number from the REACLIB database, defining the reaction structure.
 
double m_qValue = 0.0
 Q-value of the reaction in MeV.
 
std::unordered_map< fourdst::atomic::Species, size_t > m_reactants
 Reactants of the reaction.
 
std::unordered_map< fourdst::atomic::Species, size_t > m_products
 Products of the reaction.
 
std::optional< std::vector< fourdst::atomic::Species > > m_reactantsVec
 
std::optional< std::vector< fourdst::atomic::Species > > m_productsVec
 
std::string m_sourceLabel
 Source label for the rate data (e.g., "wc12w", "st08").
 
RateCoefficientSet m_rateCoefficients
 The seven rate coefficients.
 
bool m_reverse = false
 Flag indicating if this is a reverse reaction rate.
 

Private Member Functions

template<typename T>
calculate_rate (const T T9) const
 Template implementation for calculating the reaction rate.
 

Friends

std::ostream & operator<< (std::ostream &os, const ReaclibReaction &r)
 

Constructor & Destructor Documentation

◆ ~ReaclibReaction()

gridfire::reaction::ReaclibReaction::~ReaclibReaction ( )
overridedefault

◆ ReaclibReaction()

gridfire::reaction::ReaclibReaction::ReaclibReaction ( std::string_view id,
std::string_view peName,
int chapter,
const std::vector< fourdst::atomic::Species > & reactants,
const std::vector< fourdst::atomic::Species > & products,
double qValue,
std::string_view label,
const RateCoefficientSet & sets,
bool reverse = false )

Constructs a Reaction object.

Parameters
idA unique identifier for the reaction.
peNameThe name in (projectile, ejectile) notation (e.g., "p(p,g)d").
chapterThe REACLIB chapter number, defining reaction structure.
reactantsA vector of reactant species.
productsA vector of product species.
qValueThe Q-value of the reaction in MeV.
labelThe sources label for the rate data (e.g., "wc12", "st08").
setsThe set of rate coefficients.
reverseTrue if this is a reverse reaction rate.

Member Function Documentation

◆ all_species()

std::unordered_set< Species > gridfire::reaction::ReaclibReaction::all_species ( ) const
nodiscardoverridevirtual

Gets a set of all unique species involved in the reaction.

Returns
An unordered_set of all reactant and product species.

Implements gridfire::reaction::Reaction.

◆ calculate_log_rate_partial_deriv_wrt_T9()

double gridfire::reaction::ReaclibReaction::calculate_log_rate_partial_deriv_wrt_T9 ( double T9,
double rho,
double Ye,
double mue,
const fourdst::composition::Composition & comp ) const
nodiscardoverridevirtual

Logarithmic partial derivative of the rate with respect to temperature.

Implementations return d(ln rate)/d(ln T9) or an equivalent measure (as documented by the concrete class), evaluated at the provided state.

Parameters
T9Temperature in GK (10^9 K).
rhoMass density (g cm^-3).
YeElectron fraction.
mueElectron chemical potential.
compComposition object providing composition in a convenient form.
Returns
The logarithmic temperature derivative of the rate.

Implements gridfire::reaction::Reaction.

◆ calculate_rate() [1/3]

template<typename T>
T gridfire::reaction::ReaclibReaction::calculate_rate ( const T T9) const
inlinenodiscardprivate

Template implementation for calculating the reaction rate.

Template Parameters
TThe numeric type (double or CppAD::AD<double>).
Parameters
T9The temperature in units of 10^9 K.
Returns
The calculated reaction rate.

The rate is calculated using the standard REACLIB formula: rate = exp(a0 + a1/T9 + a2/T9^(1/3) + a3*T9^(1/3) + a4*T9 + a5*T9^(5/3) + a6*ln(T9))

◆ calculate_rate() [2/3]

CppAD::AD< double > gridfire::reaction::ReaclibReaction::calculate_rate ( CppAD::AD< double > T9,
CppAD::AD< double > rho,
CppAD::AD< double > Ye,
CppAD::AD< double > mue,
const std::vector< CppAD::AD< double > > & Y,
const std::unordered_map< size_t, fourdst::atomic::Species > & index_to_species_map ) const
nodiscardoverridevirtual

Calculates the reaction rate for a given temperature using CppAD types.

Parameters
T9The temperature in units of 10^9 K, as a CppAD::AD<double>.
rhoDensity, as a CppAD::AD<double> [Not used in this implementation].
Ye
mue
YMolar abundances of species, as a vector of CppAD::AD<double> [Not used in this implementation].
index_to_species_map
Returns
The calculated reaction rate, as a CppAD::AD<double>.

Implements gridfire::reaction::Reaction.

◆ calculate_rate() [3/3]

double gridfire::reaction::ReaclibReaction::calculate_rate ( double T9,
double rho,
double Ye,
double mue,
const std::vector< double > & Y,
const std::unordered_map< size_t, fourdst::atomic::Species > & index_to_species_map ) const
nodiscardoverridevirtual

Calculates the reaction rate for a given temperature.

Parameters
T9The temperature in units of 10^9 K.
rhoDensity [Not used in this implementation].
Ye
mue
Y
index_to_species_map
Returns
The calculated reaction rate.

Implements gridfire::reaction::Reaction.

◆ chapter()

int gridfire::reaction::ReaclibReaction::chapter ( ) const
inlinenodiscard

Gets the REACLIB chapter number.

Returns
The chapter number.

◆ clone()

std::unique_ptr< Reaction > gridfire::reaction::ReaclibReaction::clone ( ) const
nodiscardoverridevirtual

Polymorphic deep copy.

Returns
A std::unique_ptr owning a new Reaction equal to this one.

Implements gridfire::reaction::Reaction.

Reimplemented in gridfire::reaction::WeakReaclibReaction.

◆ contains()

bool gridfire::reaction::ReaclibReaction::contains ( const fourdst::atomic::Species & species) const
nodiscardoverridevirtual

Checks if the reaction involves a given species as a reactant or product.

Parameters
speciesThe species to check for.
Returns
True if the species is involved, false otherwise.

Implements gridfire::reaction::Reaction.

◆ contains_product()

bool gridfire::reaction::ReaclibReaction::contains_product ( const fourdst::atomic::Species & species) const
nodiscardoverridevirtual

Checks if the reaction involves a given species as a product.

Parameters
speciesThe species to check for.
Returns
True if the species is a product, false otherwise.

Implements gridfire::reaction::Reaction.

◆ contains_reactant()

bool gridfire::reaction::ReaclibReaction::contains_reactant ( const fourdst::atomic::Species & species) const
nodiscardoverridevirtual

Checks if the reaction involves a given species as a reactant.

Parameters
speciesThe species to check for.
Returns
True if the species is a reactant, false otherwise.

Implements gridfire::reaction::Reaction.

◆ countProductOccurrences()

size_t gridfire::reaction::ReaclibReaction::countProductOccurrences ( const fourdst::atomic::Species & species) const
inlineoverridevirtual

◆ countReactantOccurrences()

size_t gridfire::reaction::ReaclibReaction::countReactantOccurrences ( const fourdst::atomic::Species & species) const
inlineoverridevirtual

◆ excess_energy()

double gridfire::reaction::ReaclibReaction::excess_energy ( ) const
nodiscard

Calculates the excess energy from the mass difference of reactants and products.

Returns
The excess energy in MeV.

◆ getRateCoefficients()

std::optional< std::vector< RateCoefficientSet > > gridfire::reaction::ReaclibReaction::getRateCoefficients ( ) const
nodiscardoverridevirtual

◆ hash()

uint64_t gridfire::reaction::ReaclibReaction::hash ( uint64_t seed) const
nodiscardoverridevirtual

Computes a hash for the reaction based on its ID.

Parameters
seedThe seed for the hash function.
Returns
A 64-bit hash value.

Uses the XXHash64 algorithm on the reaction's ID string.

Implements gridfire::reaction::Reaction.

◆ id()

std::string_view gridfire::reaction::ReaclibReaction::id ( ) const
inlinenodiscardoverridevirtual

Gets the unique identifier of the reaction.

Returns
The reaction ID.

Implements gridfire::reaction::Reaction.

◆ is_reverse()

bool gridfire::reaction::ReaclibReaction::is_reverse ( ) const
inlinenodiscardoverridevirtual

Checks if this is a reverse reaction rate.

Returns
True if it is a reverse rate, false otherwise.

Implements gridfire::reaction::Reaction.

◆ num_species()

size_t gridfire::reaction::ReaclibReaction::num_species ( ) const
nodiscardoverridevirtual

Gets the number of unique species involved in the reaction.

Returns
The count of unique species.

Implements gridfire::reaction::Reaction.

◆ operator!=()

bool gridfire::reaction::ReaclibReaction::operator!= ( const ReaclibReaction & other) const
inline

Compares this reaction with another for inequality.

Parameters
otherThe other Reaction to compare with.
Returns
True if the reactions are not equal.

◆ operator==()

bool gridfire::reaction::ReaclibReaction::operator== ( const ReaclibReaction & other) const
inline

Compares this reaction with another for equality based on their IDs.

Parameters
otherThe other Reaction to compare with.
Returns
True if the reaction IDs are the same.

◆ peName()

virtual std::string_view gridfire::reaction::ReaclibReaction::peName ( ) const
inlinenodiscardvirtual

Gets the reaction name in (projectile, ejectile) notation.

Returns
The reaction name (e.g., "p(p,g)d").

◆ product_species()

std::unordered_set< Species > gridfire::reaction::ReaclibReaction::product_species ( ) const
nodiscardoverridevirtual

Gets a set of all unique product species.

Returns
An unordered_set of product species.

Implements gridfire::reaction::Reaction.

◆ products()

const std::vector< fourdst::atomic::Species > & gridfire::reaction::ReaclibReaction::products ( ) const
inlinenodiscardoverridevirtual

Gets the vector of product species.

Returns
A const reference to the vector of products.

Implements gridfire::reaction::Reaction.

◆ qValue()

double gridfire::reaction::ReaclibReaction::qValue ( ) const
inlinenodiscardoverridevirtual

Gets the Q-value of the reaction.

Returns
The Q-value in whatever units the reaction was defined in (usually MeV).

Implements gridfire::reaction::Reaction.

◆ rateCoefficients()

const RateCoefficientSet & gridfire::reaction::ReaclibReaction::rateCoefficients ( ) const
inlinenodiscard

Gets the set of rate coefficients.

Returns
A const reference to the RateCoefficientSet.

◆ reactant_species()

std::unordered_set< Species > gridfire::reaction::ReaclibReaction::reactant_species ( ) const
nodiscardoverridevirtual

Gets a set of all unique reactant species.

Returns
An unordered_set of reactant species.

Implements gridfire::reaction::Reaction.

◆ reactants()

const std::vector< fourdst::atomic::Species > & gridfire::reaction::ReaclibReaction::reactants ( ) const
inlinenodiscardoverridevirtual

Gets the vector of reactant species.

Returns
A const reference to the vector of reactants.

Implements gridfire::reaction::Reaction.

◆ sourceLabel()

std::string_view gridfire::reaction::ReaclibReaction::sourceLabel ( ) const
inlinenodiscard

Gets the source label for the rate data.

Returns
The source label (e.g., "wc12w", "st08").

◆ stoichiometry() [1/2]

std::unordered_map< Species, int > gridfire::reaction::ReaclibReaction::stoichiometry ( ) const
nodiscardoverridevirtual

Gets a map of all species to their stoichiometric coefficients.

Returns
An unordered_map from species to their integer coefficients.

Implements gridfire::reaction::Reaction.

◆ stoichiometry() [2/2]

int gridfire::reaction::ReaclibReaction::stoichiometry ( const fourdst::atomic::Species & species) const
nodiscardoverridevirtual

Calculates the stoichiometric coefficient for a given species.

Parameters
speciesThe species for which to find the coefficient.
Returns
The stoichiometric coefficient (negative for reactants, positive for products).

Implements gridfire::reaction::Reaction.

◆ type()

ReactionType gridfire::reaction::ReaclibReaction::type ( ) const
inlinenodiscardoverridevirtual

Category of this reaction (e.g., REACLIB, WEAK, LOGICAL_REACLIB).

Returns
Enumerated reaction type for runtime dispatch and filtering.

Implements gridfire::reaction::Reaction.

Reimplemented in gridfire::reaction::WeakReaclibReaction.

Friends And Related Symbol Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream & os,
const ReaclibReaction & r )
friend

Member Data Documentation

◆ m_chapter

int gridfire::reaction::ReaclibReaction::m_chapter
protected

Chapter number from the REACLIB database, defining the reaction structure.

◆ m_id

std::string gridfire::reaction::ReaclibReaction::m_id
protected

Unique identifier for the reaction (e.g., "h1+h1=>h2+e+nu").

◆ m_logger

quill::Logger* gridfire::reaction::ReaclibReaction::m_logger = fourdst::logging::LogManager::getInstance().getLogger("log")
protected

◆ m_peName

std::string gridfire::reaction::ReaclibReaction::m_peName
protected

Name of the reaction in (projectile, ejectile) notation (e.g. "p(p,g)d").

◆ m_products

std::unordered_map<fourdst::atomic::Species, size_t> gridfire::reaction::ReaclibReaction::m_products
protected

Products of the reaction.

◆ m_productsVec

std::optional<std::vector<fourdst::atomic::Species> > gridfire::reaction::ReaclibReaction::m_productsVec
mutableprotected

◆ m_qValue

double gridfire::reaction::ReaclibReaction::m_qValue = 0.0
protected

Q-value of the reaction in MeV.

◆ m_rateCoefficients

RateCoefficientSet gridfire::reaction::ReaclibReaction::m_rateCoefficients
protected

The seven rate coefficients.

◆ m_reactants

std::unordered_map<fourdst::atomic::Species, size_t> gridfire::reaction::ReaclibReaction::m_reactants
protected

Reactants of the reaction.

◆ m_reactantsVec

std::optional<std::vector<fourdst::atomic::Species> > gridfire::reaction::ReaclibReaction::m_reactantsVec
mutableprotected

◆ m_reverse

bool gridfire::reaction::ReaclibReaction::m_reverse = false
protected

Flag indicating if this is a reverse reaction rate.

◆ m_sourceLabel

std::string gridfire::reaction::ReaclibReaction::m_sourceLabel
protected

Source label for the rate data (e.g., "wc12w", "st08").


The documentation for this class was generated from the following files: