GridFire 0.6.0
General Purpose Nuclear Network
Loading...
Searching...
No Matches
gridfire Namespace Reference

Namespaces

namespace  approx8
 
namespace  exceptions
 
namespace  expectations
 
namespace  io
 
namespace  partition
 
namespace  reaclib
 
namespace  reaction
 
namespace  screening
 
namespace  solver
 
namespace  utils
 

Classes

class  AdaptiveEngineView
 An engine view that dynamically adapts the reaction network based on runtime conditions. More...
 
class  DefinedEngineView
 
class  DynamicEngine
 Abstract class for engines supporting Jacobian and stoichiometry operations. More...
 
class  Engine
 Abstract base class for a reaction network engine. More...
 
class  EngineView
 Abstract base class for a "view" of a reaction network engine. More...
 
class  FileDefinedEngineView
 
class  GraphEngine
 A reaction network engine that uses a graph-based representation. More...
 
class  MultiscalePartitioningEngineView
 An engine view that partitions the reaction network into multiple groups based on timescales. More...
 
struct  NetIn
 
struct  NetOut
 
class  Network
 
class  NetworkPrimingEngineView
 Provides a view of a DynamicEngine filtered to reactions involving a specified priming species. More...
 
struct  PrimingReport
 Captures the result of a network priming operation. More...
 
struct  QSECacheConfig
 Configuration struct for the QSE cache. More...
 
struct  QSECacheKey
 Key struct for the QSE abundance cache. More...
 
class  Reaction
 Represents a single nuclear reaction from a specific data source. More...
 
struct  StepDerivatives
 Structure holding derivatives and energy generation for a network step. More...
 

Concepts

concept  IsArithmeticOrAD
 Concept for types allowed in engine calculations.
 
concept  EngineType
 Concept for types allowed as engine bases in EngineView.
 

Typedefs

using SparsityPattern = std::vector<std::pair<size_t, size_t>>
 
typedef CppAD::AD< double > ADDouble
 Alias for CppAD AD type for double precision.
 
using BuildDepthType = std::variant<NetworkBuildDepth, int>
 Variant specifying either a predefined NetworkBuildDepth or a custom integer depth.
 
using LogicalReactionSet
 A set of logical reactions.
 
using ReactionSet
 A set of reactions, typically from a single source like REACLIB.
 

Enumerations

enum class  NetworkBuildDepth {
  Full = -1 , Shallow = 1 , SecondOrder = 2 , ThirdOrder = 3 ,
  FourthOrder = 4 , FifthOrder = 5
}
 Specifies supported depths for building the reaction network. More...
 
enum class  PrimingReportStatus {
  FULL_SUCCESS = 0 , NO_SPECIES_TO_PRIME = 1 , MAX_ITERATIONS_REACHED = 2 , FAILED_TO_FINALIZE_COMPOSITION = 3 ,
  FAILED_TO_FIND_CREATION_CHANNEL = 4 , FAILED_TO_FIND_PRIMING_REACTIONS = 5 , BASE_NETWORK_TOO_SHALLOW = 6
}
 Enumerates outcome codes for a network priming operation. More...
 
enum  NetworkFormat { APPROX8 , REACLIB , UNKNOWN }
 

Functions

reaction::LogicalReactionSet build_reaclib_nuclear_network (const fourdst::composition::Composition &composition, BuildDepthType maxLayers=NetworkBuildDepth::Full, bool reverse=false)
 Builds a nuclear reaction network from the Reaclib library based on an initial composition.
 
PrimingReport primeNetwork (const NetIn &, DynamicEngine &engine)
 Primes absent species in the network to their equilibrium abundances.
 
double calculateDestructionRateConstant (const DynamicEngine &engine, const fourdst::atomic::Species &species, const std::vector< double > &Y, double T9, double rho)
 Computes the destruction rate constant for a specific species.
 
double calculateCreationRate (const DynamicEngine &engine, const fourdst::atomic::Species &species, const std::vector< double > &Y, double T9, double rho)
 Computes the creation rate for a specific species.
 
LogicalReactionSet build_reaclib_nuclear_network (const Composition &composition, BuildDepthType maxLayers, bool reverse)
 
const reaction::ReactionfindDominantCreationChannel (const DynamicEngine &engine, const Species &species, const std::vector< double > &Y, const double T9, const double rho)
 
std::string trim_whitespace (const std::string &str)
 

Variables

static constexpr double MIN_DENSITY_THRESHOLD = 1e-18
 Minimum density threshold below which reactions are ignored.
 
static constexpr double MIN_ABUNDANCE_THRESHOLD = 1e-18
 Minimum abundance threshold below which species are ignored.
 
static constexpr double MIN_JACOBIAN_THRESHOLD = 1e-24
 Minimum value for Jacobian matrix entries.
 
std::map< PrimingReportStatus, std::string > PrimingReportStatusStrings
 Mapping from PrimingReportStatus codes to human-readable strings.
 
static std::unordered_map< NetworkFormat, std::string > FormatStringLookup
 
static int s_operator_parens_called = 0
 

Typedef Documentation

◆ ADDouble

typedef CppAD::AD<double> gridfire::ADDouble

Alias for CppAD AD type for double precision.

This alias simplifies the use of the CppAD automatic differentiation type.

◆ BuildDepthType

using gridfire::BuildDepthType = std::variant<NetworkBuildDepth, int>

Variant specifying either a predefined NetworkBuildDepth or a custom integer depth.

Precondition
If using the integer alternative, the value must be >= 0 or -1 to indicate a full build.
Postcondition
The network builder will interpret and apply the specified depth to control reaction expansion.

◆ LogicalReactionSet

A set of logical reactions.

◆ ReactionSet

A set of reactions, typically from a single source like REACLIB.

◆ SparsityPattern

using gridfire::SparsityPattern = std::vector<std::pair<size_t, size_t>>

Enumeration Type Documentation

◆ NetworkBuildDepth

enum class gridfire::NetworkBuildDepth
strong

Specifies supported depths for building the reaction network.

Values:

  • Full: Build the complete network (infinite depth).
  • Shallow: Build only direct reactions (depth = 1).
  • SecondOrder: Include reactions up to second order (depth = 2).
  • ThirdOrder: Include reactions up to third order (depth = 3).
  • FourthOrder: Include reactions up to fourth order (depth = 4).
  • FifthOrder: Include reactions up to fifth order (depth = 5).
Note
For custom build depths, see BuildDepthType.
Enumerator
Full 
Shallow 
SecondOrder 
ThirdOrder 
FourthOrder 
FifthOrder 

◆ NetworkFormat

Enumerator
APPROX8 

Approx8 nuclear reaction network format.

REACLIB 

General REACLIB nuclear reaction network format.

UNKNOWN 

◆ PrimingReportStatus

enum class gridfire::PrimingReportStatus
strong

Enumerates outcome codes for a network priming operation.

These status codes indicate the reason for success or failure of the priming process:

  • FULL_SUCCESS: Priming completed successfully with all species processed.
  • NO_SPECIES_TO_PRIME: There were no species eligible for priming.
  • MAX_ITERATIONS_REACHED: The algorithm reached its iteration limit without converging.
  • FAILED_TO_FINALIZE_COMPOSITION: Unable to build a valid Composition object at end.
  • FAILED_TO_FIND_CREATION_CHANNEL: No reaction path found to create the priming species.
  • FAILED_TO_FIND_PRIMING_REACTIONS: No reactions containing the priming species were found.
  • BASE_NETWORK_TOO_SHALLOW: The provided base network depth was insufficient for priming.
See also
PrimingReport for data associated with each status.
Enumerator
FULL_SUCCESS 
NO_SPECIES_TO_PRIME 
MAX_ITERATIONS_REACHED 
FAILED_TO_FINALIZE_COMPOSITION 
FAILED_TO_FIND_CREATION_CHANNEL 
FAILED_TO_FIND_PRIMING_REACTIONS 
BASE_NETWORK_TOO_SHALLOW 

Function Documentation

◆ build_reaclib_nuclear_network() [1/2]

LogicalReactionSet gridfire::build_reaclib_nuclear_network ( const Composition & composition,
BuildDepthType maxLayers,
bool reverse )

◆ build_reaclib_nuclear_network() [2/2]

reaction::LogicalReactionSet gridfire::build_reaclib_nuclear_network ( const fourdst::composition::Composition & composition,
BuildDepthType maxLayers = NetworkBuildDepth::Full,
bool reverse = false )

Builds a nuclear reaction network from the Reaclib library based on an initial composition.

Constructs a layered reaction network by collecting reactions up to the specified depth from the Reaclib dataset. Starting species are those with non-zero mass fractions in the input composition. Layers expand by including products of collected reactions until the depth limit. Optionally selects reverse reactions instead of forward.

See implementation in construction.cpp for details on the layering algorithm, logging, and performance.

Parameters
compositionMapping of isotopic species to their mass fractions; species with positive mass fraction seed the network.
maxLayersVariant specifying either a predefined NetworkBuildDepth or a custom integer depth; negative depth (Full) collects all reactions, zero is invalid.
reverseIf true, collects reverse reactions (decays or back-reactions); if false, uses forward reactions.
Precondition
composition must have at least one species with positive mass fraction.
Resolved integer depth from maxLayers must not be zero.
Postcondition
Returned network includes only reactions satisfying the depth and reverse criteria.
Returns
A LogicalReactionSet encapsulating the collected reactions for graph-based engines.
Exceptions
std::logic_errorIf the resolved network depth is zero (no reactions can be collected).

◆ calculateCreationRate()

double gridfire::calculateCreationRate ( const DynamicEngine & engine,
const fourdst::atomic::Species & species,
const std::vector< double > & Y,
double T9,
double rho )

Computes the creation rate for a specific species.

Sums molar reaction flows for all reactions where the species appears as a product (positive stoichiometry).

Parameters
engineEngine providing the current set of network reactions and flow calculations.
speciesThe atomic species whose creation rate is computed.
YVector of molar abundances for all species in the engine.
T9Temperature in units of 10^9 K.
rhoDensity of the medium.
Precondition
Y.size() matches engine.getNetworkReactions().size() mapping species order.
Postcondition
Returned creation rate is non-negative.
Returns
Sum of stoichiometry-weighted creation flows for the species.

◆ calculateDestructionRateConstant()

double gridfire::calculateDestructionRateConstant ( const DynamicEngine & engine,
const fourdst::atomic::Species & species,
const std::vector< double > & Y,
double T9,
double rho )

Computes the destruction rate constant for a specific species.

Calculates the sum of molar reaction flows for all reactions where the species is a reactant (negative stoichiometry) after scaling its abundance to unity.

Parameters
engineEngine providing the current set of network reactions and flow calculations.
speciesThe atomic species whose destruction rate is computed.
YVector of molar abundances for all species in the engine.
T9Temperature in units of 10^9 K.
rhoDensity of the medium.
Precondition
Y.size() matches engine.getNetworkReactions().size() mapping species order.
Postcondition
Returned rate constant is non-negative.
Returns
Sum of absolute stoichiometry-weighted destruction flows for the species.

◆ findDominantCreationChannel()

const reaction::Reaction * gridfire::findDominantCreationChannel ( const DynamicEngine & engine,
const Species & species,
const std::vector< double > & Y,
const double T9,
const double rho )

◆ primeNetwork()

PrimingReport gridfire::primeNetwork ( const NetIn & netIn,
DynamicEngine & engine )

Primes absent species in the network to their equilibrium abundances.

Executes a network priming algorithm that iteratively rebuilds the reaction network, calculates equilibrium mass fractions for species with zero initial abundance, and applies mass transfers based on reaction flows.

Refer to priming.cpp for implementation details on logging, algorithmic steps, and error handling.

Parameters
netInInput network data containing initial composition, temperature, and density.
engineDynamicEngine used to build and evaluate the reaction network.
Precondition
netIn.composition defines species and their mass fractions; engine is constructed with a valid network.
Postcondition
engine.networkReactions restored to its initial state; returned report contains primedComposition, massFractionChanges for each species, success flag, and status code.
Returns
PrimingReport encapsulating the results of the priming operation.

◆ trim_whitespace()

std::string gridfire::trim_whitespace ( const std::string & str)

Variable Documentation

◆ FormatStringLookup

std::unordered_map<NetworkFormat, std::string> gridfire::FormatStringLookup
inlinestatic
Initial value:
= {
{APPROX8, "Approx8"},
{REACLIB, "REACLIB"},
{UNKNOWN, "Unknown"}
}
@ APPROX8
Approx8 nuclear reaction network format.
Definition network.h:42
@ REACLIB
General REACLIB nuclear reaction network format.
Definition network.h:43
@ UNKNOWN
Definition network.h:44

◆ MIN_ABUNDANCE_THRESHOLD

double gridfire::MIN_ABUNDANCE_THRESHOLD = 1e-18
staticconstexpr

Minimum abundance threshold below which species are ignored.

Species with abundances below this threshold are treated as zero in reaction rate calculations. This helps to improve performance by avoiding unnecessary calculations for trace species.

◆ MIN_DENSITY_THRESHOLD

double gridfire::MIN_DENSITY_THRESHOLD = 1e-18
staticconstexpr

Minimum density threshold below which reactions are ignored.

Reactions are not calculated if the density falls below this threshold. This helps to improve performance by avoiding unnecessary calculations in very low-density regimes.

◆ MIN_JACOBIAN_THRESHOLD

double gridfire::MIN_JACOBIAN_THRESHOLD = 1e-24
staticconstexpr

Minimum value for Jacobian matrix entries.

Jacobian matrix entries with absolute values below this threshold are treated as zero to maintain sparsity and improve performance.

◆ PrimingReportStatusStrings

std::map<PrimingReportStatus, std::string> gridfire::PrimingReportStatusStrings
inline
Initial value:
= {
{PrimingReportStatus::NO_SPECIES_TO_PRIME, "No Species to Prime"},
{PrimingReportStatus::MAX_ITERATIONS_REACHED, "Max Iterations Reached"},
{PrimingReportStatus::FAILED_TO_FINALIZE_COMPOSITION, "Failed to Finalize Composition"},
{PrimingReportStatus::FAILED_TO_FIND_CREATION_CHANNEL, "Failed to Find Creation Channel"},
{PrimingReportStatus::FAILED_TO_FIND_PRIMING_REACTIONS, "Failed to Find Priming Reactions"},
{PrimingReportStatus::BASE_NETWORK_TOO_SHALLOW, "Base Network Too Shallow"}
}
@ FAILED_TO_FIND_PRIMING_REACTIONS
Definition reporting.h:37
@ MAX_ITERATIONS_REACHED
Definition reporting.h:34
@ FULL_SUCCESS
Definition reporting.h:32
@ NO_SPECIES_TO_PRIME
Definition reporting.h:33
@ FAILED_TO_FIND_CREATION_CHANNEL
Definition reporting.h:36
@ BASE_NETWORK_TOO_SHALLOW
Definition reporting.h:38
@ FAILED_TO_FINALIZE_COMPOSITION
Definition reporting.h:35

Mapping from PrimingReportStatus codes to human-readable strings.

Used when formatting or logging the priming status. No preconditions. The map contains entries for all PrimingReportStatus values.

◆ s_operator_parens_called

int gridfire::s_operator_parens_called = 0
static