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 {
26 class WeakScreeningModel final : public ScreeningModel {
27 public:
50 [[nodiscard]] std::vector<double> calculateScreeningFactors(
51 const reaction::LogicalReactionSet& reactions,
52 const std::vector<fourdst::atomic::Species>& species,
53 const std::vector<double>& Y,
54 const double T9,
55 const double rho
56 ) const override;
57
72 [[nodiscard]] std::vector<CppAD::AD<double>> calculateScreeningFactors(
73 const reaction::LogicalReactionSet& reactions,
74 const std::vector<fourdst::atomic::Species>& species,
75 const std::vector<CppAD::AD<double>>& Y,
76 const CppAD::AD<double> T9,
77 const CppAD::AD<double> rho
78 ) const override;
79 private:
81 quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
82
83 private:
99 template <typename T>
100 [[nodiscard]] std::vector<T> calculateFactors_impl(
101 const reaction::LogicalReactionSet& reactions,
102 const std::vector<fourdst::atomic::Species>& species,
103 const std::vector<T>& Y,
104 const T T9,
105 const T rho
106 ) const;
107 };
108
140 template <typename T>
142 const reaction::LogicalReactionSet& reactions,
143 const std::vector<fourdst::atomic::Species>& species,
144 const std::vector<T>& Y,
145 const T T9,
146 const T rho
147 ) const {
148 LOG_TRACE_L1(
149 m_logger,
150 "Calculating weak screening factors for {} reactions...",
151 reactions.size()
152 );
153 // --- CppAD Safe low temp checking ---
154 const T zero(0.0);
155 const T one(1.0);
156 const T low_temp_threshold(1e-9);
157
158 const T low_T_flag = CppAD::CondExpLt(T9, low_temp_threshold, zero, one);
159
160 // --- Calculate composition-dependent terms ---
161 // ζ = ∑(Z_i^2 + Z_i) * X_i / A_i
162 // This simplifies somewhat to ζ = ∑ (Z_i^2 + Z_i) * Y_i
163 // Where Y_i is the molar abundance (mol/g)
164 T zeta(0.0);
165 for (size_t i = 0; i < species.size(); ++i) {
166 const T Z = species[i].m_z;
167 zeta += (Z * Z + Z) * Y[i];
168 }
169
170 // --- Constant prefactors ---
171 const T T7 = T9 * static_cast<T>(100.00);
172 const T T7_safe = CppAD::CondExpLe(T7, low_temp_threshold, low_temp_threshold, T7);
173 const T prefactor = static_cast<T>(0.188) * CppAD::sqrt(rho / (T7_safe * T7_safe * T7_safe)) * CppAD::sqrt(zeta);
174
175 // --- Loop through reactions and calculate screening factors for each ---
176 std::vector<T> factors;
177 factors.reserve(reactions.size());
178 for (const auto& reaction : reactions) {
179 T H_12(0.0); // screening abundance term
180 const auto& reactants = reaction.reactants();
181 const bool isTripleAlpha = (
182 reactants.size() == 3 &&
183 reactants[0].m_z == 2 &&
184 reactants[1].m_z == 2 &&
185 reactants[2].m_z == 2 &&
186 reactants[0] == reactants[1] &&
187 reactants[1] == reactants[2]
188 );
189 if (reactants.size() == 2) {
190 LOG_TRACE_L3(m_logger, "Calculating screening factor for reaction: {}", reaction.peName());
191 const T Z1 = static_cast<T>(reactants[0].m_z);
192 const T Z2 = static_cast<T>(reactants[1].m_z);
193 H_12 = prefactor * Z1 * Z2;
194 }
195 else if (isTripleAlpha) {
196 LOG_TRACE_L3(m_logger, "Special case for triple alpha process in reaction: {}", reaction.peName());
197 // Special case for triple alpha process
198 const T Z_alpha = static_cast<T>(2.0);
199 const T H_alpha_alpha = prefactor * Z_alpha * Z_alpha;
200 H_12 = static_cast<T>(3.0) * H_alpha_alpha; // Triple alpha process
201 }
202 // For 1 body reactions H_12 remains 0 so e^H_12 will be 1.0 (screening does not apply)
203 // Aside from triple alpha, all other astrophysically relevant reactions are 2-body in the weak screening regime
204
205 H_12 *= low_T_flag; // Apply low temperature flag to screening factor
206 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
207 factors.push_back(CppAD::exp(H_12));
208 // std::cout << "Screening factor: " << reaction.peName() << " : " << factors.back() << "(" << H_12 << ")" << std::endl;
209 }
210 return factors;
211 }
212
213}
size_t size() const
Gets the number of reactions in the set.
Definition reaction.h:453
An abstract base class for plasma screening models.
Implements the weak screening model based on the Debye-Hückel approximation.
quill::Logger * m_logger
Logger instance for recording trace and debug information.
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
Template implementation for calculating weak screening factors.
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
Calculates weak screening factors for a set of reactions.
TemplatedReactionSet< LogicalReaction > LogicalReactionSet
A set of logical reactions.
Definition reaction.h:557
Defines classes for representing and managing nuclear reactions.