GridFire v0.7.1_rc2
General Purpose Nuclear Network
Loading...
Searching...
No Matches
gridfire::screening::BareScreeningModel Class Referencefinal

A screening model that applies no screening effect. More...

#include <screening_bare.h>

Inheritance diagram for gridfire::screening::BareScreeningModel:
[legend]
Collaboration diagram for gridfire::screening::BareScreeningModel:
[legend]

Public Member Functions

std::vector< double > calculateScreeningFactors (const reaction::ReactionSet &reactions, const std::vector< fourdst::atomic::Species > &species, const std::vector< double > &Y, double T9, double rho) const override
 Calculates screening factors, which are always 1.0.
 
std::vector< ADDoublecalculateScreeningFactors (const reaction::ReactionSet &reactions, const std::vector< fourdst::atomic::Species > &species, const std::vector< CppAD::AD< double > > &Y, ADDouble T9, ADDouble rho) const override
 Calculates screening factors for AD types, which are always 1.0.
 
- Public Member Functions inherited from gridfire::screening::ScreeningModel
virtual ~ScreeningModel ()=default
 Virtual destructor.
 

Private Types

using ADDouble = CppAD::AD< double >
 Alias for CppAD Automatic Differentiation type for double precision.
 

Private Member Functions

template<typename T >
std::vector< T > calculateFactors_impl (const reaction::ReactionSet &reactions, const std::vector< fourdst::atomic::Species > &species, const std::vector< T > &Y, const T &T9, const T &rho) const
 Template implementation for calculating screening factors.
 

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

A screening model that applies no screening effect.

This class implements the ScreeningModel interface but returns a screening factor of 1.0 for all reactions, regardless of the plasma conditions. It represents the case of bare, unscreened nuclei and serves as a baseline or can be used when screening effects are negligible or intentionally ignored.

Member Typedef Documentation

◆ ADDouble

using gridfire::screening::BareScreeningModel::ADDouble = CppAD::AD<double>
private

Alias for CppAD Automatic Differentiation type for double precision.

Member Function Documentation

◆ calculateFactors_impl()

template<typename T >
std::vector< T > gridfire::screening::BareScreeningModel::calculateFactors_impl ( const reaction::ReactionSet reactions,
const std::vector< fourdst::atomic::Species > &  species,
const std::vector< T > &  Y,
const T &  T9,
const T &  rho 
) const
private

Template implementation for calculating screening factors.

Template implementation for the bare screening model.

This private helper function contains the core logic for both the double and ADDouble versions of calculateScreeningFactors. It is templated to handle both numeric types seamlessly.

Template Parameters
TThe numeric type, either double or CppAD::AD<double>.
Parameters
reactionsThe set of reactions for which to calculate factors.
speciesA vector of all atomic species (unused).
YThe current molar composition (unused).
T9The temperature (unused).
rhoThe density (unused).
Returns
A vector of type T with all elements initialized to 1.0.

This function provides the actual implementation for calculateFactors_impl. It creates a vector of the appropriate numeric type (T) and size, and initializes all its elements to 1.0, representing no screening.

Template Parameters
TThe numeric type, either double or CppAD::AD<double>.
Parameters
reactionsThe set of reactions, used to determine the size of the output vector.
speciesUnused parameter.
YUnused parameter. @param T9 Unused parameter. @param rho Unused parameter. @return Astd::vector<T>of the same size asreactions`, with all elements set to 1.0.

◆ calculateScreeningFactors() [1/2]

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

Calculates screening factors for AD types, which are always 1.0.

This implementation returns a vector of AD-typed screening factors where every element is 1.0. This is the automatic differentiation-compatible version.

Parameters
reactionsThe set of logical reactions in the network.
speciesA vector of all atomic species (unused).
YThe current composition, providing molar abundances (mol/g) for each species (unused).
T9The temperature as an AD type (unused).
rhoThe plasma density as an AD type (unused).
Returns
A vector of ADDouble, with each element being 1.0, of the same size as the reactions set.

Implements gridfire::screening::ScreeningModel.

◆ calculateScreeningFactors() [2/2]

std::vector< double > gridfire::screening::BareScreeningModel::calculateScreeningFactors ( const reaction::ReactionSet reactions,
const std::vector< fourdst::atomic::Species > &  species,
const std::vector< double > &  Y,
double  T9,
double  rho 
) const
overridevirtual

Calculates screening factors, which are always 1.0.

This implementation returns a vector of screening factors where every element is 1.0, effectively applying no screening correction to the reaction rates.

Parameters
reactionsThe set of logical reactions in the network.
speciesA vector of all atomic species (unused).
YA vector of the molar abundances.
T9The temperature (unused).
rhoThe plasma density (unused).
Returns
A vector of doubles, with each element being 1.0, of the same size as the reactions set.

Algorithm The function simply creates and returns a std::vector<double> of the same size as the input reactions set, with all elements initialized to 1.0.

Usage

BareScreeningModel bare_model;
// ... (initialize reactions, species, Y, T9, rho)
std::vector<double> factors = bare_model.calculateScreeningFactors(
reactions, species, Y, T9, rho
);
// 'factors' will contain [1.0, 1.0, ...]
A screening model that applies no screening effect.
Definition screening_bare.h:21
std::vector< double > calculateScreeningFactors(const reaction::ReactionSet &reactions, const std::vector< fourdst::atomic::Species > &species, const std::vector< double > &Y, double T9, double rho) const override
Calculates screening factors, which are always 1.0.
Definition screening_bare.cpp:22

Implements gridfire::screening::ScreeningModel.


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