GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
gridfire::screening::WeakScreeningModel Class Referencefinal

Implements the weak screening model based on the Debye-Hückel approximation. More...

#include <screening_weak.h>

Inheritance diagram for gridfire::screening::WeakScreeningModel:
gridfire::screening::ScreeningModel

Public Member Functions

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.
 
std::vector< CppAD::AD< double > > calculateScreeningFactors (const reaction::LogicalReactionSet &reactions, const std::vector< fourdst::atomic::Species > &species, const std::vector< CppAD::AD< double > > &Y, const CppAD::AD< double > T9, const CppAD::AD< double > rho) const override
 Calculates weak screening factors using CppAD types.
 
- Public Member Functions inherited from gridfire::screening::ScreeningModel
virtual ~ScreeningModel ()=default
 Virtual destructor.
 

Private Member Functions

template<typename T>
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.
 

Private Attributes

quill::Logger * m_logger = fourdst::logging::LogManager::getInstance().getLogger("log")
 Logger instance for recording trace and debug information.
 

Additional Inherited Members

- Public Types inherited from gridfire::screening::ScreeningModel
using ADDouble = CppAD::AD<double>
 Alias for CppAD Automatic Differentiation type for double precision.
 

Detailed Description

Implements the weak screening model based on the Debye-Hückel approximation.

This class provides a concrete implementation of the ScreeningModel interface for the weak screening regime, following the formulation of Salpeter (1954). This approach applies the Debye-Hückel theory to model the electrostatic shielding of nuclei in a plasma. It is applicable to non-degenerate, non-relativistic plasmas where thermal energy dominates the electrostatic potential energy.

Definition at line 26 of file screening_weak.h.

Member Function Documentation

◆ calculateFactors_impl()

template<typename T>
std::vector< T > gridfire::screening::WeakScreeningModel::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
nodiscardprivate

Template implementation for calculating weak screening factors.

Core implementation of the weak screening calculation (Debye-Hückel model).

This private helper function contains the core logic for calculating weak screening factors. It is templated to handle both double and CppAD::AD<double> numeric types, avoiding code duplication.

Template Parameters
TThe numeric type, either double or CppAD::AD<double>.
Parameters
reactionsThe set of reactions.
speciesA vector of all species in the network.
YA vector of molar abundances.
T9The temperature in 10^9 K.
rhoThe density in g/cm^3.
Returns
A vector of screening factors of type T.

This function calculates the screening factor exp(H_12) for each reaction, based on the Debye-Hückel approximation as formulated by Salpeter (1954).

Template Parameters
TThe numeric type (double or CppAD::AD<double>).
Parameters
reactionsThe set of reactions to be screened.
speciesThe list of all species in the network.
YThe molar abundances of the species.
T9The temperature in 10^9 K.
rhoThe density in g/cm^3.
Returns
A vector of screening factors, one for each reaction.

Algorithm

  1. Low-Temperature Cutoff: If T9 is below a small threshold (1e-9), screening is effectively turned off to prevent numerical instability.
  2. Zeta Factor (ζ): A composition-dependent term is calculated: ζ = ∑(Z_i² + Z_i) * Y_i, where Z_i is the charge and Y_i is the molar abundance of species i.
  3. Prefactor: A key prefactor is computed: prefactor = 0.188 * sqrt(ρ / T₇³) * sqrt(ζ), where T₇ is the temperature in units of 10^7 K.
  4. Screening Term (H_12): For each reaction, the term H_12 is calculated:
    • For a two-body reaction (reactants Z₁ and Z₂): H_12 = prefactor * Z₁ * Z₂.
    • For the triple-alpha reaction (3 * He4): H_12 = 3 * (prefactor * Z_α * Z_α).
    • For one-body reactions (decays), H_12 is 0, so the factor is 1.
  5. Capping: The value of H_12 is capped at 2.0 to prevent excessively large and unphysical screening factors (exp(2) ≈ 7.4).
  6. Final Factor: The screening factor for the reaction is exp(H_12).

Definition at line 141 of file screening_weak.h.

◆ calculateScreeningFactors() [1/2]

std::vector< ADDouble > gridfire::screening::WeakScreeningModel::calculateScreeningFactors ( const reaction::LogicalReactionSet & reactions,
const std::vector< fourdst::atomic::Species > & species,
const std::vector< CppAD::AD< double > > & Y,
const CppAD::AD< double > T9,
const CppAD::AD< double > rho ) const
nodiscardoverridevirtual

Calculates weak screening factors using CppAD types.

This is the automatic differentiation-compatible version of the method. It allows the derivatives of the screening factors to be computed with respect to plasma conditions.

Parameters
reactionsThe set of logical reactions in the network.
speciesA vector of all atomic species involved in the network.
YA vector of the molar abundances as AD types.
T9The temperature as an AD type.
rhoThe plasma density as an AD type.
Returns
A vector of screening factors as AD types.

Implements gridfire::screening::ScreeningModel.

Definition at line 12 of file screening_weak.cpp.

◆ calculateScreeningFactors() [2/2]

std::vector< double > gridfire::screening::WeakScreeningModel::calculateScreeningFactors ( const reaction::LogicalReactionSet & reactions,
const std::vector< fourdst::atomic::Species > & species,
const std::vector< double > & Y,
const double T9,
const double rho ) const
nodiscardoverridevirtual

Calculates weak screening factors for a set of reactions.

This method computes the screening enhancement factor for each reaction based on the Salpeter (1954) formula.

Parameters
reactionsThe set of logical reactions in the network.
speciesA vector of all atomic species involved in the network.
YA vector of the molar abundances (mol/g) for each species.
T9The temperature in units of 10^9 K.
rhoThe plasma density in g/cm^3.
Returns
A vector of screening factors (dimensionless), one for each reaction.

Usage

WeakScreeningModel weak_model;
// ... (initialize reactions, species, Y, T9, rho)
std::vector<double> factors = weak_model.calculateScreeningFactors(
reactions, species, Y, T9, rho
);
Implements the weak screening model based on the Debye-Hückel approximation.
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.

Implements gridfire::screening::ScreeningModel.

Definition at line 22 of file screening_weak.cpp.

Member Data Documentation

◆ m_logger

quill::Logger* gridfire::screening::WeakScreeningModel::m_logger = fourdst::logging::LogManager::getInstance().getLogger("log")
private

Logger instance for recording trace and debug information.

Definition at line 81 of file screening_weak.h.


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