fix(weakRates): major progress in resolving bugs

bigs were introduced by the interface change from accepting raw molar abundance vectors to using the composition vector. This commit resolves many of these, including preformant ways to report that a species is not present in the composition and unified index lookups using composition object tooling.

BREAKING CHANGE:
This commit is contained in:
2025-10-10 09:12:40 -04:00
parent 13e2ea9ffa
commit 2f1077c02d
21 changed files with 17953 additions and 375 deletions

View File

@@ -18,6 +18,7 @@
#include <vector>
#include <memory>
#include <ranges>
#include <functional>
#include <boost/numeric/ublas/matrix_sparse.hpp>
@@ -707,6 +708,7 @@ namespace gridfire {
* @param rho Density in g/cm^3.
* @param Ye
* @param mue
* @param speciesIDLookup
* @return Molar flow rate for the reaction (e.g., mol/g/s).
*
* This method computes the net rate at which the given reaction proceeds
@@ -716,14 +718,17 @@ namespace gridfire {
T calculateMolarReactionFlow(
const reaction::Reaction &reaction,
const std::vector<T>& Y,
const T T9,
const T rho, T Ye, T mue
T T9,
T rho,
T Ye,
T mue,
const std::function<std::optional<size_t>(const fourdst::atomic::Species &)>&speciesIDLookup
) const;
template<IsArithmeticOrAD T>
T calculateReverseMolarReactionFlow(
const T T9,
const T rho,
T T9,
T rho,
std::vector<T> screeningFactors,
const std::vector<T>& Y,
size_t reactionIndex,
@@ -739,6 +744,7 @@ namespace gridfire {
* @param rho Density in g/cm^3.
* @param Ye
* @param mue
* @param speciesLookup
* @return StepDerivatives<T> containing dY/dt and energy generation rate.
*
* This method calculates the time derivatives of all species and the
@@ -748,7 +754,10 @@ namespace gridfire {
[[nodiscard]] StepDerivatives<T> calculateAllDerivatives(
const std::vector<T>& Y_in,
T T9,
T rho, T Ye, T mue
T rho,
T Ye,
T mue,
std::function<std::optional<size_t>(const fourdst::atomic::Species &)> speciesLookup
) const;
// /**
@@ -869,7 +878,8 @@ namespace gridfire {
const T T9,
const T rho,
const T Ye,
const T mue
const T mue,
const std::function<std::optional<size_t>(const fourdst::atomic::Species &)> speciesLookup
) const {
std::vector<T> screeningFactors = m_screeningModel->calculateScreeningFactors(
m_reactions,
@@ -913,16 +923,21 @@ namespace gridfire {
const T N_A = static_cast<T>(m_constants.Na); // Avogadro's number in mol^-1
const T c = static_cast<T>(m_constants.c); // Speed of light in cm/s
// TODO: It may be prudent to introduce assertions here which validate units but will be removed in release builds (to ensure that unit inconsistencies do not creep in during future development)
// libconstants already has units built in so this should be straightforward.
// --- SINGLE LOOP OVER ALL REACTIONS ---
for (size_t reactionIndex = 0; reactionIndex < m_reactions.size(); ++reactionIndex) {
const auto& reaction = m_reactions[reactionIndex];
// 1. Calculate forward reaction rate
const T forwardMolarReactionFlow = screeningFactors[reactionIndex] *
calculateMolarReactionFlow<T>(reaction, Y, T9, rho, Ye, mue);
calculateMolarReactionFlow<T>(
reaction,
Y,
T9,
rho,
Ye,
mue,
speciesLookup
);
// 2. Calculate reverse reaction rate
T reverseMolarFlow = static_cast<T>(0.0);
@@ -965,7 +980,8 @@ namespace gridfire {
const T T9,
const T rho,
const T Ye,
const T mue
const T mue,
const std::function<std::optional<size_t>(const fourdst::atomic::Species &)>& speciesIDLookup
) const {
// --- Pre-setup (flags to control conditionals in an AD safe / branch aware manner) ---
@@ -989,10 +1005,12 @@ namespace gridfire {
// --- Loop through each unique reactant species and calculate the molar concentration for that species then multiply that into the accumulator ---
for (const auto& [species_name, count] : reactant_counts) {
// --- Resolve species to molar abundance ---
// PERF: Could probably optimize out this lookup
const auto species_it = m_speciesToIndexMap.find(m_networkSpeciesMap.at(species_name));
const size_t species_index = species_it->second;
const T Yi = Y[species_index];
// TODO: We need some way to handle the case when a species in the reaction is not part of the composition being tracked
const std::optional<size_t> species_index = speciesIDLookup(m_networkSpeciesMap.at(species_name));
if (!species_index.has_value()) {
return static_cast<T>(0.0); // If any reactant is not present, the reaction cannot proceed
}
const T Yi = Y[species_index.value()];
// --- If count is > 1 , we need to raise the molar concentration to the power of count since there are really count bodies in that reaction ---
molar_concentration_product *= CppAD::pow(Yi, static_cast<T>(count)); // ni^count

View File

@@ -2,6 +2,7 @@
#include <exception>
#include <string>
#include <utility>
#include <vector>
namespace gridfire::exceptions {
@@ -17,41 +18,41 @@ namespace gridfire::exceptions {
int m_total_steps;
double m_eps_nuc;
};
explicit StaleEngineTrigger(const state &s)
: m_state(s) {}
explicit StaleEngineTrigger(state s)
: m_state(std::move(s)) {}
const char* what() const noexcept override{
[[nodiscard]] const char* what() const noexcept override{
return "Engine reports stale state. This means that the caller should trigger a update of the engine state before continuing with the integration. If you as an end user are seeing this error, it is likely a bug in the code that should be reported. Please provide the input parameters and the context in which this error occurred. Thank you for your help!";
}
state getState() const {
[[nodiscard]] state getState() const {
return m_state;
}
size_t numSpecies() const {
[[nodiscard]] size_t numSpecies() const {
return m_state.m_Y.size();
}
size_t totalSteps() const {
[[nodiscard]] size_t totalSteps() const {
return m_state.m_total_steps;
}
double energy() const {
[[nodiscard]] double energy() const {
return m_state.m_eps_nuc;
}
double getMolarAbundance(const size_t index) const {
[[nodiscard]] double getMolarAbundance(const size_t index) const {
if (index > m_state.m_Y.size() - 1) {
throw std::out_of_range("Index out of bounds for molar abundance vector.");
}
return m_state.m_Y[index];
}
double temperature() const {
[[nodiscard]] double temperature() const {
return m_state.m_T9 * 1e9; // Convert T9 back to Kelvin
}
double density() const {
[[nodiscard]] double density() const {
return m_state.m_rho;
}
private:
@@ -61,10 +62,10 @@ namespace gridfire::exceptions {
class StaleEngineError final : public EngineError {
public:
explicit StaleEngineError(const std::string& message)
: m_message(message) {}
explicit StaleEngineError(std::string message)
: m_message(std::move(message)) {}
const char* what() const noexcept override {
[[nodiscard]] const char* what() const noexcept override {
return m_message.c_str();
}
@@ -74,10 +75,10 @@ namespace gridfire::exceptions {
class FailedToPartitionEngineError final : public EngineError {
public:
explicit FailedToPartitionEngineError(const std::string& message)
: m_message(message) {}
explicit FailedToPartitionEngineError(std::string message)
: m_message(std::move(message)) {}
const char* what() const noexcept override {
[[nodiscard]] const char* what() const noexcept override {
return m_message.c_str();
}
private:
@@ -86,10 +87,10 @@ namespace gridfire::exceptions {
class NetworkResizedError final : public EngineError {
public:
explicit NetworkResizedError(const std::string& message)
: m_message(message) {}
explicit NetworkResizedError(std::string message)
: m_message(std::move(message)) {}
const char* what() const noexcept override {
[[nodiscard]] const char* what() const noexcept override {
return m_message.c_str();
}
private:
@@ -98,10 +99,10 @@ namespace gridfire::exceptions {
class UnableToSetNetworkReactionsError final : public EngineError {
public:
explicit UnableToSetNetworkReactionsError(const std::string& message)
: m_message(message) {}
explicit UnableToSetNetworkReactionsError(std::string message)
: m_message(std::move(message)) {}
const char* what() const noexcept override {
[[nodiscard]] const char* what() const noexcept override {
return m_message.c_str();
}

View File

@@ -42,7 +42,7 @@ namespace gridfire::expectations {
}
};
struct StaleEngineError : EngineError {
struct StaleEngineError final : EngineError {
StaleEngineErrorTypes staleType;
explicit StaleEngineError(const StaleEngineErrorTypes sType)

View File

@@ -6,7 +6,6 @@
#include "fourdst/logging/logging.h"
#include <string>
#include <unordered_map>
#include <vector>
#include <memory>
@@ -99,6 +98,6 @@ namespace gridfire::partition {
* @return Unique pointer to a new PartitionFunction instance of the given type.
* @throws std::runtime_error If the given type is not recognized.
+ */
[[nodiscard]] std::unique_ptr<PartitionFunction> selectPartitionFunction(const BasePartitionType type) const;
[[nodiscard]] std::unique_ptr<PartitionFunction> selectPartitionFunction(BasePartitionType type) const;
};
}

View File

@@ -1,6 +1,6 @@
#pragma once
#define GRIDFIRE_WEAK_REACTION_LIB_SENTINEL -60.0
#define GRIDFIRE_WEAK_REACTION_LIB_SENTINEL (-60.0)
#include "gridfire/reaction/reaction.h"
#include "gridfire/reaction/weak/weak_types.h"
@@ -541,21 +541,74 @@ namespace gridfire::rates::weak {
m_atomic(ax, ay);
rateConstant = static_cast<T>(ay[0]);
} else { // The case where T is of type double
const std::expected<WeakRatePayload, InterpolationError> result = m_interpolator.get_rates(
std::expected<WeakRatePayload, InterpolationError> result = m_interpolator.get_rates(
static_cast<uint16_t>(m_reactant_a),
static_cast<uint8_t>(m_reactant_z),
T9,
log_rhoYe,
mue
log_rhoYe
);
// TODO: Clean this up. When a bit of code needs this many comments to make it clear it is bad code
if (!result.has_value()) {
const InterpolationErrorType type = result.error().type;
const std::string msg = std::format(
"Failed to interpolate weak rate for (A={}, Z={}) at T9={}, log10(rho*Ye)={}, mu_e={} with error: {}",
m_reactant.name(), m_reactant_a, m_reactant_z, T9, log_rhoYe, mue, InterpolationErrorTypeMap.at(type)
bool okayToClamp = true;
const auto&[errorType, boundsErrorInfo] = result.error();
// The logic here is
// 1. If there is any bounds error in T9 then we do not allow clamping. T9 should be a large enough grid
// that the user should not be asking for values outside the grid.
// 2. If there is no bounds error in T9, but there is a bounds error in log_rhoYe, then we only allow
// clamping if the query value is below the minimum of the grid. If it is above the maximum
// of the grid, then we do not allow clamping. The reason for this is that at high density,
// screening and other effects can make a significant difference to the rates, and
// the user should be aware that they are asking for a value outside the grid.
// There are a couple of safety asserts in here that are only active in debug builds. These are to
// ensure that our assumptions about the error information are correct. These should really never
// be triggered, but if they are, they will help us to identify any issues.
if (errorType == InterpolationErrorType::BOUNDS_ERROR) {
assert(boundsErrorInfo.has_value()); // must be true if type is BOUNDS_ERROR, removed in release builds
if (boundsErrorInfo->contains(TableAxes::T9)) {
okayToClamp = false;
} else {
assert(boundsErrorInfo->contains(TableAxes::LOG_RHOYE)); // must be true if T9 is not, removed in release builds
const BoundsErrorInfo& boundsError = boundsErrorInfo->at(TableAxes::LOG_RHOYE);
if (boundsError.queryValue > boundsError.axisMaxValue) {
okayToClamp = false;
}
assert(boundsError.queryValue < boundsError.axisMinValue); // Given the above logic, this must be true, removed in release builds
}
}
if (!okayToClamp) {
const InterpolationErrorType type = result.error().type;
const std::string msg = std::format(
"Failed to interpolate weak rate for {} (A={}, Z={}) at T9={}, log10(rho*Ye)={}, with error: {}. Clamping disallowed due to either query value being out of bounds in T9 or being above the maximum in log10(rho*Ye).",
m_reactant.name(), m_reactant_a, m_reactant_z, T9, log_rhoYe, InterpolationErrorTypeMap.at(type)
);
throw std::runtime_error(msg);
}
// In the case we get here the error was a bounds error in log_rhoYe and the query value was below the minimum of the grid
// so the solution is to clamp the query value to the minimum of the grid and try again.
result = m_interpolator.get_rates(
static_cast<uint16_t>(m_reactant_a),
static_cast<uint8_t>(m_reactant_z),
T9,
boundsErrorInfo->at(TableAxes::LOG_RHOYE).axisMinValue
);
throw std::runtime_error(msg);
// Check the result again. If it fails this time then we have a real problem and we throw.
if (!result.has_value()) {
const InterpolationErrorType type = result.error().type;
const std::string msg = std::format(
"After clamping, failed to interpolate weak rate for {} (A={}, Z={}) at T9={}, log10(rho*Ye)={}, with error: {}",
m_reactant.name(), m_reactant_a, m_reactant_z, T9, log_rhoYe, InterpolationErrorTypeMap.at(type)
);
throw std::runtime_error(msg);
}
}
const WeakRatePayload payload = result.value();
const double logRate = get_log_rate_from_payload(payload);

View File

@@ -2,6 +2,7 @@
#include "gridfire/reaction/weak/weak_types.h"
#include "fourdst/composition/atomicSpecies.h"
#include "fourdst/logging/logging.h"
#include <unordered_map>
#include <cstdint>
@@ -66,7 +67,6 @@ namespace gridfire::rates::weak {
* @param Z Proton number of the isotope.
* @param t9 Temperature in GK (10^9 K).
* @param log_rhoYe Log10 of rho*Ye (cgs density times electron fraction).
* @param mu_e Electron chemical potential (MeV).
* @return expected<WeakRatePayload, InterpolationError>: payload on success;
* InterpolationError::UNKNOWN_SPECIES_ERROR if (A,Z) not present; or
* InterpolationError::BOUNDS_ERROR if any coordinate is outside the table
@@ -84,8 +84,7 @@ namespace gridfire::rates::weak {
uint16_t A,
uint8_t Z,
double t9,
double log_rhoYe,
double mu_e
double log_rhoYe
) const;
/**
@@ -100,7 +99,6 @@ namespace gridfire::rates::weak {
* @param Z Proton number of the isotope.
* @param t9 Temperature in GK (10^9 K).
* @param log_rhoYe Log10 of rho*Ye (cgs density times electron fraction).
* @param mu_e Electron chemical potential (MeV).
* @return expected<WeakRateDerivatives, InterpolationError>: derivative payload on success;
* otherwise an InterpolationError as described above.
* @par Example
@@ -114,10 +112,10 @@ namespace gridfire::rates::weak {
uint16_t A,
uint8_t Z,
double t9,
double log_rhoYe,
double mu_e
double log_rhoYe
) const;
private:
quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
/**
* @brief Pack (A,Z) into a 32-bit key used for the internal map.
*

View File

@@ -93,18 +93,18 @@ namespace gridfire::rates::weak {
};
/**
* @brief Partial derivatives of the log10() fields w.r.t. (T9, log10(rho*Ye), mu_e).
* @brief Partial derivatives of the log10() fields w.r.t. (T9, log10(rho*Ye)).
*
* Array ordering is [d/dT9, d/dlogRhoYe, d/dMuE] for each corresponding field.
* Array ordering is [d/dT9, d/dlogRhoYe] for each corresponding field.
*/
struct WeakRateDerivatives {
// Each array holds [d/dT9, d/dlogRhoYe, d/dMuE]
std::array<double, 3> d_log_beta_plus;
std::array<double, 3> d_log_electron_capture;
std::array<double, 3> d_log_neutrino_loss_ec;
std::array<double, 3> d_log_beta_minus;
std::array<double, 3> d_log_positron_capture;
std::array<double, 3> d_log_antineutrino_loss_bd;
std::array<double, 2> d_log_beta_plus;
std::array<double, 2> d_log_electron_capture;
std::array<double, 2> d_log_neutrino_loss_ec;
std::array<double, 2> d_log_beta_minus;
std::array<double, 2> d_log_positron_capture;
std::array<double, 2> d_log_antineutrino_loss_bd;
};
/**
@@ -133,6 +133,19 @@ namespace gridfire::rates::weak {
LOG_RHOYE, ///< log10(rho*Ye).
MUE ///< Electron chemical potential (MeV).
};
}
// This need to be here to avoid compiler issues related to the order of specialization
namespace std {
template <>
struct hash<gridfire::rates::weak::TableAxes> {
std::size_t operator()(gridfire::rates::weak::TableAxes t) const noexcept {
return std::hash<int>()(static_cast<int>(t));
}
};
}
namespace gridfire::rates::weak {
/**
* @brief Detailed bounds information for a BOUNDS_ERROR.
@@ -154,20 +167,20 @@ namespace gridfire::rates::weak {
std::optional<std::unordered_map<TableAxes, BoundsErrorInfo>> boundsErrorInfo = std::nullopt;
};
/**
* @brief Regular 3D grid and payloads for a single isotope (A,Z).
* @brief Regular 2D grid and payloads for a single isotope (A,Z).
*
* Axes are monotonically increasing per dimension. Data vector is laid out in
* row-major order with index computed as:
* index = ((i_t9 * rhoYe_axis.size() + j_rhoYe) * mue_axis.size()) + k_mue
*
* index = i_t9 * N_rhoYe + j_rhoYe
*
*/
struct IsotopeGrid {
std::vector<double> t9_axis; ///< Unique sorted T9 grid.
std::vector<double> rhoYe_axis;///< Unique sorted log10(rho*Ye) grid.
std::vector<double> mue_axis; ///< Unique sorted mu_e grid.
// index = ((i_t9 * rhoYe_axis.size() + j_rhoYe) * mue_axis.size()) + k_mue
std::vector<WeakRatePayload> data; ///< Payloads at each grid node.
std::vector<double> t9_axis; ///< Unique sorted T9 grid.
std::vector<double> rhoYe_axis; ///< Unique sorted log10(rho*Ye) grid.
std::vector<WeakRatePayload> data; ///< MuE axis for each (T9, log_rhoYe) pair (the table is ragged in mu_e). This is also where the payloads are stored.
};
/**
@@ -216,4 +229,5 @@ namespace gridfire::rates::weak {
return os;
}
};
}
}

View File

@@ -21,11 +21,15 @@
// Include headers for linear solvers and N_Vectors
// We will use preprocessor directives to select the correct ones
#include <cvode/cvode.h> // For CVDls (serial dense linear solver)
#include <sundials/sundials_context.h>
#include <sunmatrix/sunmatrix_dense.h>
#include <sunlinsol/sunlinsol_dense.h>
// These are the possible N_Vector implementations. We use the compiler defines to select the most appropriate one for the build.
// If none are defined, we default to serial.
// For precompiled binaries we will need to ensure that we have versions built for all three types (ideally with some runtime
// checks that will fail gracefully if the user tries to use an unsupported type).
#ifdef SUNDIALS_HAVE_OPENMP
#include <nvector/nvector_openmp.h>
#endif
@@ -103,7 +107,7 @@ namespace gridfire::solver {
* if present after a step, it is rethrown for upstream handling.
* - Prints/collects diagnostics per step (step size, energy, solver iterations).
* - On trigger activation, rebuilds CVODE resources to reflect a changed network and
* reinitializes the state using the latest engine composition, preserving energy.
* reinitialized the state using the latest engine composition, preserving energy.
* - At the end, converts molar abundances to mass fractions and assembles NetOut,
* including derivatives of energy w.r.t. T and rho from the engine.
*
@@ -250,7 +254,7 @@ namespace gridfire::solver {
* @brief Compute and print per-component error ratios; run diagnostic helpers.
*
* Gathers CVODE's estimated local errors, converts the state to a Composition, and prints a
* sorted table of species with highest error ratios; then invokes diagnostic routines to
* sorted table of species with the highest error ratios; then invokes diagnostic routines to
* inspect Jacobian stiffness and species balance.
*/
void log_step_diagnostics(const CVODEUserData& user_data) const;

View File

@@ -1,10 +1,9 @@
#pragma once
#include "gridfire/engine/engine_abstract.h"
#include "fourdst/composition/composition.h"
#include <string>
#include <vector>
namespace gridfire::utils {
/**
@@ -17,7 +16,7 @@ namespace gridfire::utils {
*
* @param engine A constant reference to a `DynamicEngine` object, used to
* calculate the species timescales.
* @param Y A vector of the molar abundances (mol/g) for each species.
* @param composition The current composition of the plasma
* @param T9 The temperature in units of 10^9 K.
* @param rho The plasma density in g/cm^3.
* @return A std::string containing the formatted table of species and their