GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
screening_weak.h
Go to the documentation of this file.
1#pragma once
2
5
6#include "fourdst/logging/logging.h"
7#include "quill/Logger.h"
8#include "quill/LogMacros.h"
9
10#include "cppad/cppad.hpp"
11
12namespace gridfire::screening {
13 class WeakScreeningModel final : public ScreeningModel {
14 public:
15 [[nodiscard]] std::vector<double> calculateScreeningFactors(
16 const reaction::LogicalReactionSet& reactions,
17 const std::vector<fourdst::atomic::Species>& species,
18 const std::vector<double>& Y,
19 const double T9,
20 const double rho
21 ) const override;
22
23 [[nodiscard]] std::vector<CppAD::AD<double>> calculateScreeningFactors(
24 const reaction::LogicalReactionSet& reactions,
25 const std::vector<fourdst::atomic::Species>& species,
26 const std::vector<CppAD::AD<double>>& Y,
27 const CppAD::AD<double> T9,
28 const CppAD::AD<double> rho
29 ) const override;
30 private:
31 quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
32
33 private:
34 template <typename T>
35 [[nodiscard]] std::vector<T> calculateFactors_impl(
36 const reaction::LogicalReactionSet& reactions,
37 const std::vector<fourdst::atomic::Species>& species,
38 const std::vector<T>& Y,
39 const T T9,
40 const T rho
41 ) const;
42 };
43
44 template <typename T>
46 const reaction::LogicalReactionSet& reactions,
47 const std::vector<fourdst::atomic::Species>& species,
48 const std::vector<T>& Y,
49 const T T9,
50 const T rho
51 ) const {
52 LOG_TRACE_L1(
54 "Calculating weak screening factors for {} reactions...",
55 reactions.size()
56 );
57 // --- CppAD Safe low temp checking ---
58 const T zero(0.0);
59 const T one(1.0);
60 const T low_temp_threshold(1e-9);
61
62 const T low_T_flag = CppAD::CondExpLt(T9, low_temp_threshold, zero, one);
63
64 // --- Calculate composition-dependent terms ---
65 // ζ = ∑(Z_i^2 + Z_i) * X_i / A_i
66 // This simplifies somewhat to ζ = ∑ (Z_i^2 + Z_i) * Y_i
67 // Where Y_i is the molar abundance (mol/g)
68 T zeta(0.0);
69 for (size_t i = 0; i < species.size(); ++i) {
70 const T Z = species[i].m_z;
71 zeta += (Z * Z + Z) * Y[i];
72 }
73
74 // --- Constant prefactors ---
75 const T T7 = T9 * static_cast<T>(100.00);
76 const T T7_safe = CppAD::CondExpLe(T7, low_temp_threshold, low_temp_threshold, T7);
77 const T prefactor = static_cast<T>(0.188) * CppAD::sqrt(rho / (T7_safe * T7_safe * T7_safe)) * CppAD::sqrt(zeta);
78
79 // --- Loop through reactions and calculate screening factors for each ---
80 std::vector<T> factors;
81 factors.reserve(reactions.size());
82 for (const auto& reaction : reactions) {
83 T H_12(0.0); // screening abundance term
84 const auto& reactants = reaction.reactants();
85 const bool isTripleAlpha = (
86 reactants.size() == 3 &&
87 reactants[0].m_z == 2 &&
88 reactants[1].m_z == 2 &&
89 reactants[2].m_z == 2 &&
90 reactants[0] == reactants[1] &&
91 reactants[1] == reactants[2]
92 );
93 if (reactants.size() == 2) {
94 LOG_TRACE_L3(m_logger, "Calculating screening factor for reaction: {}", reaction.peName());
95 const T Z1 = static_cast<T>(reactants[0].m_z);
96 const T Z2 = static_cast<T>(reactants[1].m_z);
97 H_12 = prefactor * Z1 * Z2;
98 }
99 else if (isTripleAlpha) {
100 LOG_TRACE_L3(m_logger, "Special case for triple alpha process in reaction: {}", reaction.peName());
101 // Special case for triple alpha process
102 const T Z_alpha = static_cast<T>(2.0);
103 const T H_alpha_alpha = prefactor * Z_alpha * Z_alpha;
104 H_12 = static_cast<T>(3.0) * H_alpha_alpha; // Triple alpha process
105 }
106 // For 1 body reactions H_12 remains 0 so e^H_12 will be 1.0 (screening does not apply)
107 // Aside from triple alpha, all other astrophysically relevant reactions are 2-body in the weak screening regime
108
109 H_12 *= low_T_flag; // Apply low temperature flag to screening factor
110 H_12 = CppAD::CondExpGe(H_12, static_cast<T>(2.0), static_cast<T>(2.0), H_12); // Caps the screening factor at 10 to avoid numerical issues
111 factors.push_back(CppAD::exp(H_12));
112 // std::cout << "Screening factor: " << reaction.peName() << " : " << factors.back() << "(" << H_12 << ")" << std::endl;
113 }
114 return factors;
115 }
116
117}
size_t size() const
Gets the number of reactions in the set.
Definition reaction.h:453
std::vector< T > calculateFactors_impl(const reaction::LogicalReactionSet &reactions, const std::vector< fourdst::atomic::Species > &species, const std::vector< T > &Y, const T T9, const T rho) const
std::vector< double > calculateScreeningFactors(const reaction::LogicalReactionSet &reactions, const std::vector< fourdst::atomic::Species > &species, const std::vector< double > &Y, const double T9, const double rho) const override
TemplatedReactionSet< LogicalReaction > LogicalReactionSet
A set of logical reactions.
Definition reaction.h:557
Defines classes for representing and managing nuclear reactions.