|
GridFire 0.6.0
General Purpose Nuclear Network
|
Namespaces | |
| namespace | approx8 |
| namespace | diagnostics |
| namespace | exceptions |
| namespace | expectations |
| namespace | io |
| namespace | partition |
| namespace | rates |
| namespace | reaclib |
| namespace | reaction |
| namespace | screening |
| namespace | solver |
| namespace | trigger |
| 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... | |
| struct | EnergyDerivatives |
| 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... | |
| class | ReactionSet |
| 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. | |
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::ReactionSet | build_reaclib_nuclear_network (const fourdst::composition::Composition &composition, const rates::weak::WeakRateInterpolator &weakInterpolator, 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 &netIn, DynamicEngine &engine) |
| Primes absent species in the network to their equilibrium abundances. | |
| double | calculateDestructionRateConstant (const DynamicEngine &engine, const fourdst::atomic::Species &species, const fourdst::composition::Composition &composition, double T9, double rho) |
| Computes the destruction rate constant for a specific species. | |
| double | calculateCreationRate (const DynamicEngine &engine, const fourdst::atomic::Species &species, const fourdst::composition::Composition &composition, double T9, double rho) |
| Computes the creation rate for a specific species. | |
| ReactionSet | build_nuclear_network (const Composition &composition, const rates::weak::WeakRateInterpolator &weak_interpolator, BuildDepthType maxLayers, bool reverse_reaclib) |
| const reaction::Reaction * | findDominantCreationChannel (const DynamicEngine &engine, const Species &species, const fourdst::composition::Composition &comp, const double T9, const double rho) |
| double | calculateDestructionRateConstant (const DynamicEngine &engine, const Species &species, const Composition &comp, const double T9, const double rho) |
| double | calculateCreationRate (const DynamicEngine &engine, const Species &species, const Composition &comp, 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 |
| 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.
| using gridfire::BuildDepthType = std::variant<NetworkBuildDepth, int> |
Variant specifying either a predefined NetworkBuildDepth or a custom integer depth.
| using gridfire::SparsityPattern = std::vector<std::pair<size_t, size_t>> |
|
strong |
Specifies supported depths for building the reaction network.
Values:
| Enumerator | |
|---|---|
| Full | |
| Shallow | |
| SecondOrder | |
| ThirdOrder | |
| FourthOrder | |
| FifthOrder | |
|
strong |
Enumerates outcome codes for a network priming operation.
These status codes indicate the reason for success or failure of the priming process:
| 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 | |
| ReactionSet gridfire::build_nuclear_network | ( | const Composition & | composition, |
| const rates::weak::WeakRateInterpolator & | weak_interpolator, | ||
| BuildDepthType | maxLayers, | ||
| bool | reverse_reaclib ) |
| reaction::ReactionSet gridfire::build_reaclib_nuclear_network | ( | const fourdst::composition::Composition & | composition, |
| const rates::weak::WeakRateInterpolator & | weakInterpolator, | ||
| 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.
| composition | Mapping of isotopic species to their mass fractions; species with positive mass fraction seed the network. |
| weakInterpolator | |
| maxLayers | Variant specifying either a predefined NetworkBuildDepth or a custom integer depth; negative depth (Full) collects all reactions, zero is invalid. |
| reverse | If true, collects reverse reactions (decays or back-reactions); if false, uses forward reactions. |
| std::logic_error | If the resolved network depth is zero (no reactions can be collected). |
| double gridfire::calculateCreationRate | ( | const DynamicEngine & | engine, |
| const fourdst::atomic::Species & | species, | ||
| const fourdst::composition::Composition & | composition, | ||
| 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).
| engine | Engine providing the current set of network reactions and flow calculations. |
| species | The atomic species whose creation rate is computed. |
| composition | Composition object containing current abundances. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density of the medium. |
| double gridfire::calculateCreationRate | ( | const DynamicEngine & | engine, |
| const Species & | species, | ||
| const Composition & | comp, | ||
| const double | T9, | ||
| const double | rho ) |
| double gridfire::calculateDestructionRateConstant | ( | const DynamicEngine & | engine, |
| const fourdst::atomic::Species & | species, | ||
| const fourdst::composition::Composition & | composition, | ||
| 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.
| engine | Engine providing the current set of network reactions and flow calculations. |
| species | The atomic species whose destruction rate is computed. |
| composition | Current composition providing abundances for all species. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density of the medium. |
| double gridfire::calculateDestructionRateConstant | ( | const DynamicEngine & | engine, |
| const Species & | species, | ||
| const Composition & | comp, | ||
| const double | T9, | ||
| const double | rho ) |
| const reaction::Reaction * gridfire::findDominantCreationChannel | ( | const DynamicEngine & | engine, |
| const Species & | species, | ||
| const fourdst::composition::Composition & | comp, | ||
| const double | T9, | ||
| const double | rho ) |
| PrimingReport gridfire::primeNetwork | ( | const NetIn & | netIn, |
| DynamicEngine & | engine ) |
Primes absent species in the network to their equilibrium abundances.
Primes absent species in the network to their equilibrium abundances using a robust, two-stage approach.
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.
| netIn | Input network data containing initial composition, temperature, and density. |
| engine | DynamicEngine used to build and evaluate the reaction network. |
This function implements a robust network priming algorithm that avoids the pitfalls of sequential, one-by-one priming. The previous, brittle method could allow an early priming reaction to consume all of a shared reactant, starving later reactions. This new, two-stage method ensures that all priming reactions are considered collectively, competing for the same limited pool of initial reactants in a physically consistent manner.
The algorithm proceeds in three main stages:
scalingFactor. If any reactant is overdrawn, this factor will be less than 1.0, ensuring that no reactant's abundance can go negative.scalingFactor before being applied to the final composition. This ensures that if resources are limited, all primed species are scaled down proportionally.| netIn | Input network data containing initial composition, temperature, and density. |
| engine | DynamicEngine used to build and evaluate the reaction network. |
| std::string gridfire::trim_whitespace | ( | const std::string & | str | ) |
|
inlinestatic |
|
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.
|
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.
|
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.
|
inline |
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.