fix(MultiscalePartitioningEngineView): began work to call the QSE solver every timestep

one major issue was that QSE solving was only running at each partition. This was creating effectivley infinite sources of partitioned species. Now we partition when engine updates are triggered, however, solveQSEAbundance is called every timestep. This has major performance implications and so has required a lot of optimization to make it even somewhat viable. For now construction is much slower. Time per iteration is still slower than it was before; however, it is tractable. There is also currently too much stiffness in the network. This is likeley a bug that was introduced in this refactoring which will be addressed soon.
This commit is contained in:
2025-11-12 16:54:12 -05:00
parent 741a11d256
commit 81cca35130
20 changed files with 446 additions and 697 deletions

View File

@@ -98,7 +98,7 @@ namespace gridfire {
// --- The public facing interface can always use the precomputed version since taping is done internally ---
return calculateAllDerivativesUsingPrecomputation(comp, bare_rates, bare_reverse_rates, T9, rho, activeReactions);
} else {
return calculateAllDerivatives<double>(
StepDerivatives<double> result = calculateAllDerivatives<double>(
comp.getMolarAbundanceVector(),
T9,
rho,
@@ -115,6 +115,7 @@ namespace gridfire {
return false;
}
);
return result;
}
}
@@ -238,8 +239,6 @@ namespace gridfire {
throw std::runtime_error("Species not found in global atomic species database: " + std::string(name));
}
}
// TODO: Currently this works. We sort the vector based on mass so that for the same set of species we always get the same ordering and we get the same ordering as a composition with the same set of species
// However, we need some checks so that when we get a composition we confirm that it is the same ordering / contains the same species. This is important for the ODE integrator to work properly.
std::ranges::sort(m_networkSpecies, [](const fourdst::atomic::Species& a, const fourdst::atomic::Species& b) -> bool {
return a.mass() < b.mass(); // Otherwise, sort by mass
});
@@ -581,7 +580,9 @@ namespace gridfire {
}
fourdst::composition::Composition GraphEngine::collectComposition(
fourdst::composition::CompositionAbstract &comp
const fourdst::composition::CompositionAbstract &comp,
double T9,
double rho
) const {
for (const auto &species: comp.getRegisteredSpecies()) {
if (!m_networkSpeciesMap.contains(species.name())) {
@@ -615,6 +616,7 @@ namespace gridfire {
rho
);
// --- Optimized loop ---
std::vector<double> molarReactionFlows;
molarReactionFlows.reserve(m_precomputedReactions.size());
@@ -654,6 +656,7 @@ namespace gridfire {
forwardAbundanceProduct *
std::pow(rho, numReactants > 1 ? static_cast<double>(numReactants) - 1 : 0.0);
// --- Reverse reaction flow ---
// Only do this is the reaction has a non-zero reverse symmetry factor (i.e. is reversible)
double reverseMolarReactionFlow = 0.0;
@@ -697,11 +700,35 @@ namespace gridfire {
const int stoichiometricCoefficient = precomp.stoichiometric_coefficients[i];
// Update the derivative for this species
result.dydt.at(species) += static_cast<double>(stoichiometricCoefficient) * R_j;
double dydt_increment = static_cast<double>(stoichiometricCoefficient) * R_j;
result.dydt.at(species) += dydt_increment;
result.reactionContributions[species][std::string(reaction->id())] = dydt_increment;
}
reactionCounter++;
}
// std::vector<std::string> reactionIDs;
// for (const auto& reaction: activeReactions) {
// reactionIDs.push_back(std::string(reaction->id()));
// }
//
// std::vector<std::unique_ptr<utils::ColumnBase>> columns;
// columns.push_back(std::make_unique<utils::Column<std::string>>("Reaction", reactionIDs));
// for (const auto& [species, contributions] : reactionContributions) {
// std::vector<double> speciesData;
// for (const auto& reactionID : reactionIDs) {
// if (contributions.contains(reactionID)) {
// speciesData.push_back(contributions.at(reactionID));
// } else {
// speciesData.push_back(0.0);
// }
// }
// columns.push_back(std::make_unique<utils::Column<double>>(std::string(species.name()), speciesData));
// }
// utils::print_table("Contributions", columns);
// exit(0);
// --- Calculate the nuclear energy generation rate ---
double massProductionRate = 0.0; // [mol][s^-1]
for (const auto & species : m_networkSpecies) {
@@ -1103,7 +1130,7 @@ namespace gridfire {
) const {
const double Ye = comp.getElectronAbundance();
auto [dydt, _] = calculateAllDerivatives<double>(
auto [dydt, _, __] = calculateAllDerivatives<double>(
comp.getMolarAbundanceVector(),
T9,
rho,
@@ -1155,7 +1182,7 @@ namespace gridfire {
return std::nullopt; // Species not present
};
auto [dydt, _] = calculateAllDerivatives<double>(
auto [dydt, _, __] = calculateAllDerivatives<double>(
Y,
T9,
rho,
@@ -1246,7 +1273,7 @@ namespace gridfire {
// 5. Call the actual templated function
// We let T9 and rho be constant, so we pass them as fixed values.
auto [dydt, nuclearEnergyGenerationRate] = calculateAllDerivatives<CppAD::AD<double>>(
auto [dydt, nuclearEnergyGenerationRate, _] = calculateAllDerivatives<CppAD::AD<double>>(
adY,
adT9,
adRho,

View File

@@ -147,8 +147,8 @@ namespace gridfire {
const double rho
) const {
validateState();
// TODO: Think about if I need to reach in and adjust the composition to zero out inactive species.
auto result = m_baseEngine.calculateRHSAndEnergy(comp, T9, rho);
fourdst::composition::Composition collectedComp = collectComposition(comp, T9, rho);
auto result = m_baseEngine.calculateRHSAndEnergy(collectedComp, T9, rho);
if (!result) {
return std::unexpected{result.error()};
@@ -313,9 +313,11 @@ namespace gridfire {
}
fourdst::composition::Composition AdaptiveEngineView::collectComposition(
fourdst::composition::CompositionAbstract &comp
const fourdst::composition::CompositionAbstract &comp,
const double T9,
const double rho
) const {
fourdst::composition::Composition result = m_baseEngine.collectComposition(comp); // Step one is to bubble the results from lower levels of the engine chain up
fourdst::composition::Composition result = m_baseEngine.collectComposition(comp, T9, rho);
for (const auto& species : m_activeSpecies) {
if (!result.contains(species)) {

View File

@@ -284,9 +284,11 @@ namespace gridfire {
}
fourdst::composition::Composition DefinedEngineView::collectComposition(
fourdst::composition::CompositionAbstract &comp
const fourdst::composition::CompositionAbstract &comp,
const double T9,
const double rho
) const {
fourdst::composition::Composition result = m_baseEngine.collectComposition(comp);
fourdst::composition::Composition result = m_baseEngine.collectComposition(comp, T9, rho);
for (const auto& species : m_activeSpecies) {
if (!result.contains(species)) {

View File

@@ -17,34 +17,11 @@
#include "quill/LogMacros.h"
#include "quill/Logger.h"
#include "xxhash64.h"
#include "fourdst/composition/utils/composition_hash.h"
namespace {
using namespace fourdst::atomic;
//TODO: Replace all calls to this function with composition.getMolarAbundanceVector() so that
// we don't have to keep this function around. (Cant do this right now because there is not a
// guarantee that this function will return the same ordering as the canonical vector representation)
std::vector<double> packCompositionToVector(
const fourdst::composition::Composition& composition,
const gridfire::DynamicEngine& engine
) {
std::vector<double> Y(engine.getNetworkSpecies().size(), 0.0);
const auto& allSpecies = engine.getNetworkSpecies();
for (size_t i = 0; i < allSpecies.size(); ++i) {
const auto& species = allSpecies[i];
if (!composition.contains(species)) {
Y[i] = 0.0; // Species not in the composition, set to zero
} else {
Y[i] = composition.getMolarAbundance(species);
}
}
return Y;
}
template <class T>
void hash_combine(std::size_t& seed, const T& v) {
std::hash<T> hashed;
seed ^= hashed(v) + 0x9e3779b9 + (seed << 6) + (seed >> 2);
}
std::vector<std::vector<Species>> findConnectedComponentsBFS(
const std::unordered_map<Species, std::vector<Species>>& graph,
const std::vector<Species>& nodes
@@ -186,11 +163,14 @@ namespace gridfire {
const double T9,
const double rho
) const {
// Check the cache to see if the network needs to be repartitioned. Note that the QSECacheKey manages binning of T9, rho, and Y_full to ensure that small changes (which would likely not result in a repartitioning) do not trigger a cache miss.
const auto result = m_baseEngine.calculateRHSAndEnergy(comp, T9, rho);
const fourdst::composition::Composition qseComposition = getNormalizedEquilibratedComposition(comp, T9, rho);
const auto result = m_baseEngine.calculateRHSAndEnergy(qseComposition, T9, rho);
if (!result) {
return std::unexpected{result.error()};
}
auto deriv = result.value();
for (const auto& species : m_algebraic_species) {
@@ -204,7 +184,8 @@ namespace gridfire {
const double T9,
const double rho
) const {
return m_baseEngine.calculateEpsDerivatives(comp, T9, rho);
const fourdst::composition::Composition qseComposition = getNormalizedEquilibratedComposition(comp, T9, rho);
return m_baseEngine.calculateEpsDerivatives(qseComposition, T9, rho);
}
void MultiscalePartitioningEngineView::generateJacobianMatrix(
@@ -212,8 +193,8 @@ namespace gridfire {
const double T9,
const double rho
) const {
// We do not need to generate the jacobian for QSE species since those entries are by definition 0
m_baseEngine.generateJacobianMatrix(comp, T9, rho, m_dynamic_species);
const fourdst::composition::Composition qseComposition = getNormalizedEquilibratedComposition(comp, T9, rho);
m_baseEngine.generateJacobianMatrix(qseComposition, T9, rho, m_dynamic_species);
}
void MultiscalePartitioningEngineView::generateJacobianMatrix(
@@ -249,7 +230,9 @@ namespace gridfire {
}
}
m_baseEngine.generateJacobianMatrix(comp, T9, rho, dynamicActiveSpeciesIntersection);
const fourdst::composition::Composition qseComposition = getNormalizedEquilibratedComposition(comp, T9, rho);
m_baseEngine.generateJacobianMatrix(qseComposition, T9, rho, dynamicActiveSpeciesIntersection);
}
void MultiscalePartitioningEngineView::generateJacobianMatrix(
@@ -258,7 +241,8 @@ namespace gridfire {
const double rho,
const SparsityPattern &sparsityPattern
) const {
return m_baseEngine.generateJacobianMatrix(comp, T9, rho, sparsityPattern);
const fourdst::composition::Composition qseComposition = getNormalizedEquilibratedComposition(comp, T9, rho);
return m_baseEngine.generateJacobianMatrix(qseComposition, T9, rho, sparsityPattern);
}
double MultiscalePartitioningEngineView::getJacobianMatrixEntry(
@@ -292,19 +276,9 @@ namespace gridfire {
const double T9,
const double rho
) const {
// Fix the algebraic species to the equilibrium abundances we calculate.
fourdst::composition::Composition comp_mutable;
for (const auto& species : comp.getRegisteredSpecies()) {
comp_mutable.registerSpecies(species);
comp_mutable.setMolarAbundance(species, comp.getMolarAbundance(species));
}
const fourdst::composition::Composition qseComposition = getNormalizedEquilibratedComposition(comp, T9, rho);
for (const auto& species : m_algebraic_species) {
const double Yi = m_algebraic_abundances.at(species);
comp_mutable.setMolarAbundance(species, Yi);
}
return m_baseEngine.calculateMolarReactionFlow(reaction, comp_mutable, T9, rho);
return m_baseEngine.calculateMolarReactionFlow(reaction, qseComposition, T9, rho);
}
const reaction::ReactionSet & MultiscalePartitioningEngineView::getNetworkReactions() const {
@@ -321,7 +295,8 @@ namespace gridfire {
const double T9,
const double rho
) const {
const auto result = m_baseEngine.getSpeciesTimescales(comp, T9, rho);
const fourdst::composition::Composition qseComposition = getNormalizedEquilibratedComposition(comp, T9, rho);
const auto result = m_baseEngine.getSpeciesTimescales(qseComposition, T9, rho);
if (!result) {
return std::unexpected{result.error()};
}
@@ -338,7 +313,8 @@ namespace gridfire {
const double T9,
const double rho
) const {
const auto result = m_baseEngine.getSpeciesDestructionTimescales(comp, T9, rho);
const fourdst::composition::Composition qseComposition = getNormalizedEquilibratedComposition(comp, T9, rho);
const auto result = m_baseEngine.getSpeciesDestructionTimescales(qseComposition, T9, rho);
if (!result) {
return std::unexpected{result.error()};
}
@@ -354,27 +330,13 @@ namespace gridfire {
NetIn baseUpdatedNetIn = netIn;
baseUpdatedNetIn.composition = baseUpdatedComposition;
const fourdst::composition::Composition equilibratedComposition = equilibrateNetwork(baseUpdatedNetIn);
std::unordered_map<Species, double> algebraicAbundances;
for (const auto& species : m_algebraic_species) {
algebraicAbundances.emplace(species, equilibratedComposition.getMolarAbundance(species));
}
m_algebraic_abundances = std::move(algebraicAbundances);
fourdst::composition::Composition equilibratedComposition = partitionNetwork(baseUpdatedNetIn);
return equilibratedComposition;
}
bool MultiscalePartitioningEngineView::isStale(const NetIn &netIn) {
const auto key = QSECacheKey(
netIn.temperature,
netIn.density,
packCompositionToVector(netIn.composition, m_baseEngine)
);
if (m_qse_abundance_cache.contains(key)) {
return m_baseEngine.isStale(netIn); // The cache hit indicates the engine is not stale for the given conditions.
}
return true;
return m_baseEngine.isStale(netIn);
}
void MultiscalePartitioningEngineView::setScreeningModel(
@@ -392,10 +354,7 @@ namespace gridfire {
}
std::vector<std::vector<Species>> MultiscalePartitioningEngineView::analyzeTimescalePoolConnectivity(
const std::vector<std::vector<Species>> &timescale_pools,
const fourdst::composition::Composition &comp,
double T9,
double rho
const std::vector<std::vector<Species>> &timescale_pools
) const {
std::vector<std::vector<Species>> final_connected_pools;
@@ -413,17 +372,22 @@ namespace gridfire {
}
void MultiscalePartitioningEngineView::partitionNetwork(
const fourdst::composition::Composition &comp,
const double T9,
const double rho
fourdst::composition::Composition MultiscalePartitioningEngineView::partitionNetwork(
const NetIn &netIn
) {
// --- Step 0. Clear previous state ---
// --- Step 0. Prime the network ---
const PrimingReport primingReport = m_baseEngine.primeEngine(netIn);
const fourdst::composition::Composition& comp = primingReport.primedComposition;
const double T9 = netIn.temperature / 1e9;
const double rho = netIn.density;
// --- Step 0.5 Clear previous state ---
LOG_TRACE_L1(m_logger, "Partitioning network...");
LOG_TRACE_L1(m_logger, "Clearing previous state...");
m_qse_groups.clear();
m_dynamic_species.clear();
m_algebraic_species.clear();
m_composition_cache.clear(); // We need to clear the cache now cause the same comp, temp, and density may result in a different value
// --- Step 1. Identify distinct timescale regions ---
LOG_TRACE_L1(m_logger, "Identifying fast reactions...");
@@ -449,7 +413,7 @@ namespace gridfire {
}
LOG_TRACE_L1(m_logger, "Preforming connectivity analysis on timescale pools...");
const std::vector<std::vector<Species>> connected_pools = analyzeTimescalePoolConnectivity(candidate_pools, comp, T9, rho);
const std::vector<std::vector<Species>> connected_pools = analyzeTimescalePoolConnectivity(candidate_pools);
LOG_TRACE_L1(m_logger, "Found {} connected pools (compared to {} timescale pools) for QSE analysis.", connected_pools.size(), timescale_pools.size());
@@ -555,15 +519,12 @@ namespace gridfire {
m_qse_groups.size(),
m_qse_groups.size() == 1 ? "" : "s"
);
}
void MultiscalePartitioningEngineView::partitionNetwork(
const NetIn &netIn
) {
const double T9 = netIn.temperature / 1e9; // Convert temperature from Kelvin to T9 (T9 = T / 1e9)
const double rho = netIn.density; // Density in g/cm^3
partitionNetwork(netIn.composition, T9, rho);
// Sort the QSE groups by mean timescale so that fastest groups get equilibrated first (as these may feed slower groups)
std::ranges::sort(m_qse_groups, [](const QSEGroup& a, const QSEGroup& b) {
return a.mean_timescale < b.mean_timescale;
});
return getNormalizedEquilibratedComposition(comp, T9, rho);
}
void MultiscalePartitioningEngineView::exportToDot(
@@ -593,6 +554,7 @@ namespace gridfire {
}
}
const fourdst::composition::Composition qseComposition = getNormalizedEquilibratedComposition(comp, T9, rho);
// Calculate reaction flows and find min/max for logarithmic scaling of transparency
std::vector<double> reaction_flows;
reaction_flows.reserve(all_reactions.size());
@@ -600,7 +562,7 @@ namespace gridfire {
double max_log_flow = std::numeric_limits<double>::lowest();
for (const auto& reaction : all_reactions) {
double flow = std::abs(m_baseEngine.calculateMolarReactionFlow(*reaction, comp, T9, rho));
double flow = std::abs(m_baseEngine.calculateMolarReactionFlow(*reaction, qseComposition, T9, rho));
reaction_flows.push_back(flow);
if (flow > 1e-99) { // Avoid log(0)
double log_flow = std::log10(flow);
@@ -776,7 +738,7 @@ namespace gridfire {
std::vector<double> MultiscalePartitioningEngineView::mapNetInToMolarAbundanceVector(const NetIn &netIn) const {
std::vector<double> Y(m_dynamic_species.size(), 0.0); // Initialize with zeros
std::vector<double> Y(netIn.composition.size(), 0.0); // Initialize with zeros
for (const auto& [sp, y] : netIn.composition) {
Y[getSpeciesIndex(sp)] = y; // Map species to their molar abundance
}
@@ -804,33 +766,6 @@ namespace gridfire {
return m_baseEngine.primeEngine(netIn);
}
fourdst::composition::Composition MultiscalePartitioningEngineView::equilibrateNetwork(
const fourdst::composition::Composition &comp,
const double T9,
const double rho
) {
partitionNetwork(comp, T9, rho);
fourdst::composition::Composition qseComposition = solveQSEAbundances(comp, T9, rho);
for (const auto &[sp, y]: qseComposition) {
if (y < 0.0 && std::abs(y) < 1e-20) {
qseComposition.setMolarAbundance(sp, 0.0); // Avoid negative mass fractions
}
}
return qseComposition;
}
fourdst::composition::Composition MultiscalePartitioningEngineView::equilibrateNetwork(
const NetIn &netIn
) {
const PrimingReport primingReport = m_baseEngine.primeEngine(netIn);
const double T9 = netIn.temperature / 1e9; // Convert temperature from Kelvin to T9 (T9 = T / 1e9)
const double rho = netIn.density; // Density in g/cm^3
return equilibrateNetwork(primingReport.primedComposition, T9, rho);
}
bool MultiscalePartitioningEngineView::involvesSpecies(
const Species &species
@@ -852,30 +787,43 @@ namespace gridfire {
return std::ranges::find(m_dynamic_species, species) != m_dynamic_species.end();
}
fourdst::composition::Composition MultiscalePartitioningEngineView::collectComposition(
fourdst::composition::CompositionAbstract &comp
fourdst::composition::Composition MultiscalePartitioningEngineView::getNormalizedEquilibratedComposition(
const fourdst::composition::CompositionAbstract& comp,
const double T9,
const double rho
) const {
fourdst::composition::Composition result = m_baseEngine.collectComposition(comp);
std::map<Species, double> Ym; // Use an ordered map here so that this is ordered by atomic mass (which is the </> comparator for Species)
for (const auto& [sp, y] : result) {
Ym.emplace(sp, y);
const std::array<uint64_t, 3> hashes = {
fourdst::composition::utils::CompositionHash::hash_exact(comp),
std::hash<double>()(T9),
std::hash<double>()(rho)
};
const uint64_t composite_hash = XXHash64::hash(hashes.begin(), sizeof(uint64_t) * 3, 0);
if (m_composition_cache.contains(composite_hash)) {
LOG_TRACE_L1(m_logger, "Cache Hit in Multiscale Partitioning Engine View for composition at T9 = {}, rho = {}.", T9, rho);
return m_composition_cache.at(composite_hash);
}
for (const auto& [species, Yi] : m_algebraic_abundances) {
if (!Ym.contains(species)) {
throw exceptions::BadCollectionError("MuiltiscalePartitioningEngineView failed to collect composition for species " + std::string(species.name()) + " as the base engine did not report that species present in its composition!");
fourdst::composition::Composition qseComposition = solveQSEAbundances(comp, T9, rho);
for (const auto &[sp, y]: qseComposition) {
if (y < 0.0 && std::abs(y) < 1e-20) {
qseComposition.setMolarAbundance(sp, 0.0); // normalize small negative abundances to zero
}
Ym.at(species) = Yi;
}
std::vector<double> Y;
std::vector<std::string> speciesNames;
Y.reserve(Ym.size());
m_composition_cache[composite_hash] = qseComposition;
for (const auto& [species, Yi] : Ym) {
Y.emplace_back(Yi);
speciesNames.emplace_back(species.name());
}
return qseComposition;
}
return {speciesNames, Y};
fourdst::composition::Composition MultiscalePartitioningEngineView::collectComposition(
const fourdst::composition::CompositionAbstract &comp,
const double T9,
const double rho
) const {
const fourdst::composition::Composition result = m_baseEngine.collectComposition(comp, T9, rho);
fourdst::composition::Composition qseComposition = solveQSEAbundances(result, T9, rho);
return qseComposition;
}
size_t MultiscalePartitioningEngineView::getSpeciesIndex(const Species &species) const {
@@ -1201,33 +1149,21 @@ namespace gridfire {
}
fourdst::composition::Composition MultiscalePartitioningEngineView::solveQSEAbundances(
const fourdst::composition::Composition &comp,
const fourdst::composition::CompositionAbstract &comp,
const double T9,
const double rho
) {
) const {
LOG_TRACE_L1(m_logger, "Solving for QSE abundances...");
// Sort by timescale to solve fastest QSE groups first (these can seed slower groups)
std::ranges::sort(m_qse_groups, [](const QSEGroup& a, const QSEGroup& b) {
return a.mean_timescale < b.mean_timescale;
});
fourdst::composition::Composition outputComposition = comp;
fourdst::composition::Composition outputComposition(comp.getRegisteredSpecies());
for (const auto& [sp, y]: comp) {
outputComposition.setMolarAbundance(sp, y);
}
for (const auto&[is_in_equilibrium, algebraic_species, seed_species, mean_timescale] : m_qse_groups) {
if (!is_in_equilibrium || (algebraic_species.empty() && seed_species.empty())) {
continue;
}
fourdst::composition::Composition normalized_composition = comp;
for (const auto& species: algebraic_species) {
if (!normalized_composition.contains(species)) {
normalized_composition.registerSpecies(species);
}
}
for (const auto& species: seed_species) {
if (!normalized_composition.contains(species)) {
normalized_composition.registerSpecies(species);
}
}
Eigen::VectorXd Y_scale(algebraic_species.size());
Eigen::VectorXd v_initial(algebraic_species.size());
@@ -1235,14 +1171,14 @@ namespace gridfire {
std::unordered_map<Species, size_t> species_to_index_map;
for (const auto& species : algebraic_species) {
constexpr double abundance_floor = 1.0e-100;
const double initial_abundance = normalized_composition.getMolarAbundance(species);
const double initial_abundance = comp.getMolarAbundance(species);
const double Y = std::max(initial_abundance, abundance_floor);
v_initial(i) = std::log(Y);
species_to_index_map.emplace(species, i);
i++;
}
EigenFunctor functor(*this, algebraic_species, normalized_composition, T9, rho, Y_scale, species_to_index_map);
EigenFunctor functor(*this, algebraic_species, comp, T9, rho, Y_scale, species_to_index_map);
Eigen::LevenbergMarquardt lm(functor);
lm.parameters.ftol = 1.0e-10;
lm.parameters.xtol = 1.0e-10;
@@ -1275,7 +1211,7 @@ namespace gridfire {
m_logger,
"During QSE solving species {} started with a molar abundance of {} and ended with an abundance of {}.",
species.name(),
normalized_composition.getMolarAbundance(species),
comp.getMolarAbundance(species),
Y_final_qse(i)
);
// double Xi = Y_final_qse(i) * species.mass(); // Convert from molar abundance to mass fraction
@@ -1488,23 +1424,28 @@ namespace gridfire {
}
//////////////////////////////////////
/// Eigen Functor Member Functions ///
/////////////////////////////////////
int MultiscalePartitioningEngineView::EigenFunctor::operator()(const InputType &v_qse, OutputType &f_qse) const {
fourdst::composition::Composition comp_trial = m_initial_comp;
fourdst::composition::Composition comp_trial(m_initial_comp.getRegisteredSpecies());
for (const auto& [sp, y] : m_initial_comp) {
comp_trial.setMolarAbundance(sp, y);
}
Eigen::VectorXd y_qse = v_qse.array().exp(); // Convert to physical abundances using exponential scaling
for (const auto& species: m_qse_solve_species) {
if (!comp_trial.contains(std::string(species.name()))) {
comp_trial.registerSpecies(species);
}
auto index = static_cast<long>(m_qse_solve_species_index_map.at(species));
comp_trial.setMolarAbundance(species, y_qse[index]);
}
const auto result = m_view->getBaseEngine().calculateRHSAndEnergy(comp_trial, m_T9, m_rho);
const auto result = m_view.getBaseEngine().calculateRHSAndEnergy(comp_trial, m_T9, m_rho);
if (!result) {
throw exceptions::StaleEngineError("Failed to calculate RHS and energy due to stale engine state");
}
const auto&[dydt, nuclearEnergyGenerationRate] = result.value();
const auto&[dydt, nuclearEnergyGenerationRate, _] = result.value();
f_qse.resize(static_cast<long>(m_qse_solve_species.size()));
long i = 0;
// TODO: make sure that just counting up i is a valid approach, this is a possible place an indexing bug may have crept in
@@ -1514,7 +1455,7 @@ namespace gridfire {
i++;
}
LOG_TRACE_L2(
m_view->m_logger,
m_view.m_logger,
"Functor evaluation at T9 = {}, rho = {}, y_qse = <{}> complete. ||f|| = {}",
m_T9,
m_rho,
@@ -1531,7 +1472,7 @@ namespace gridfire {
f_qse.norm()
);
LOG_TRACE_L3(
m_view->m_logger,
m_view.m_logger,
"{}",
[&]() -> std::string {
std::stringstream ss;
@@ -1546,24 +1487,24 @@ namespace gridfire {
}
int MultiscalePartitioningEngineView::EigenFunctor::df(const InputType &v_qse, JacobianType &J_qse) const {
fourdst::composition::Composition comp_trial = m_initial_comp;
fourdst::composition::Composition comp_trial(m_initial_comp.getRegisteredSpecies());
for (const auto& [sp, y] : m_initial_comp) {
comp_trial.setMolarAbundance(sp, y);
}
Eigen::VectorXd y_qse = v_qse.array().exp(); // Convert to physical abundances using exponential scaling
for (const auto& species: m_qse_solve_species) {
if (!comp_trial.contains(std::string(species.name()))) {
comp_trial.registerSpecies(species);
}
const double molarAbundance = y_qse[static_cast<long>(m_qse_solve_species_index_map.at(species))];
comp_trial.setMolarAbundance(species, molarAbundance);
}
std::vector<Species> qse_species_vector(m_qse_solve_species.begin(), m_qse_solve_species.end());
m_view->getBaseEngine().generateJacobianMatrix(comp_trial, m_T9, m_rho, qse_species_vector);
const auto result = m_view->getBaseEngine().calculateRHSAndEnergy(comp_trial, m_T9, m_rho);
m_view.getBaseEngine().generateJacobianMatrix(comp_trial, m_T9, m_rho, qse_species_vector);
const auto result = m_view.getBaseEngine().calculateRHSAndEnergy(comp_trial, m_T9, m_rho);
if (!result) {
throw exceptions::StaleEngineError("Failed to calculate RHS and energy due to stale engine state");
}
const auto&[dydt, nuclearEnergyGenerationRate] = result.value();
const auto&[dydt, nuclearEnergyGenerationRate, _] = result.value();
const long N = static_cast<long>(m_qse_solve_species.size());
J_qse.resize(N, N);
@@ -1571,12 +1512,12 @@ namespace gridfire {
for (const auto& rowSpecies : m_qse_solve_species) {
long colID = 0;
for (const auto& colSpecies: m_qse_solve_species) {
J_qse(rowID, colID) = m_view->getBaseEngine().getJacobianMatrixEntry(
J_qse(rowID, colID) = m_view.getBaseEngine().getJacobianMatrixEntry(
rowSpecies,
colSpecies
);
colID += 1;
LOG_TRACE_L3(m_view->m_logger, "Jacobian[{}, {}] (d(dY({}))/dY({})) = {}", rowID, colID - 1, rowSpecies.name(), colSpecies.name(), J_qse(rowID, colID - 1));
LOG_TRACE_L3(m_view.m_logger, "Jacobian[{}, {}] (d(dY({}))/dY({})) = {}", rowID, colID - 1, rowSpecies.name(), colSpecies.name(), J_qse(rowID, colID - 1));
}
rowID += 1;
}
@@ -1597,46 +1538,12 @@ namespace gridfire {
}
QSECacheKey::QSECacheKey(
const double T9,
const double rho,
const std::vector<double> &Y
) :
m_T9(T9),
m_rho(rho),
m_Y(Y) {
m_hash = hash();
}
size_t QSECacheKey::hash() const {
std::size_t seed = 0;
/////////////////////////////////
/// QSEGroup Member Functions ///
////////////////////////////////
hash_combine(seed, m_Y.size());
hash_combine(seed, bin(m_T9, m_cacheConfig.T9_tol));
hash_combine(seed, bin(m_rho, m_cacheConfig.rho_tol));
double negThresh = 1e-10; // Threshold for considering a value as negative.
for (double Yi : m_Y) {
if (Yi < 0.0 && std::abs(Yi) < negThresh) {
Yi = 0.0; // Avoid negative abundances
} else if (Yi < 0.0 && std::abs(Yi) >= negThresh) {
throw std::invalid_argument("Yi should be positive for valid hashing (expected Yi > 0, received " + std::to_string(Yi) + ")");
}
hash_combine(seed, bin(Yi, m_cacheConfig.Yi_tol));
}
return seed;
}
long QSECacheKey::bin(const double value, const double tol) {
return static_cast<long>(std::floor(value / tol));
}
bool QSECacheKey::operator==(const QSECacheKey &other) const {
return m_hash == other.m_hash;
}
bool MultiscalePartitioningEngineView::QSEGroup::operator==(const QSEGroup &other) const {
return mean_timescale == other.mean_timescale;
@@ -1707,36 +1614,4 @@ namespace gridfire {
return ss.str();
}
void MultiscalePartitioningEngineView::CacheStats::hit(const operators op) {
if (op == operators::All) {
throw std::invalid_argument("Cannot use 'ALL' as an operator for a hit");
}
m_hit ++;
m_operatorHits[op]++;
}
void MultiscalePartitioningEngineView::CacheStats::miss(const operators op) {
if (op == operators::All) {
throw std::invalid_argument("Cannot use 'ALL' as an operator for a miss");
}
m_miss ++;
m_operatorMisses[op]++;
}
size_t MultiscalePartitioningEngineView::CacheStats::hits(const operators op) const {
if (op == operators::All) {
return m_hit;
}
return m_operatorHits.at(op);
}
size_t MultiscalePartitioningEngineView::CacheStats::misses(const operators op) const {
if (op == operators::All) {
return m_miss;
}
return m_operatorMisses.at(op);
}
}

View File

@@ -55,6 +55,7 @@ namespace gridfire::policy {
);
auto& graphRepr = dynamic_cast<GraphEngine&>(*m_network_stack.back().get());
// graphRepr.setPrecomputation(false);
graphRepr.setUseReverseReactions(false);

View File

@@ -4,27 +4,99 @@
#include "gridfire/reaction/reaclib.h"
#include "gridfire/reaction/reactions_data.h"
#include "gridfire/network.h"
#include "gridfire/exceptions/error_reaction.h"
#include <stdexcept>
#include <sstream>
#include <vector>
#include <string>
#include <format>
#include <print>
#include <expected>
std::string trim_whitespace(const std::string& str) {
auto startIt = str.begin();
const auto endIt = str.end();
namespace {
std::string trim_whitespace(const std::string& str) {
auto startIt = str.begin();
const auto endIt = str.end();
while (startIt != endIt && std::isspace(static_cast<unsigned char>(*startIt))) {
++startIt;
while (startIt != endIt && std::isspace(static_cast<unsigned char>(*startIt))) {
++startIt;
}
if (startIt == endIt) {
return "";
}
const auto ritr = std::find_if(str.rbegin(), std::string::const_reverse_iterator(startIt),
[](const unsigned char ch){ return !std::isspace(ch); });
return {startIt, ritr.base()};
}
if (startIt == endIt) {
return "";
enum class ReactionParseError {
MissingOpenParenthesis,
MissingCloseParenthesis,
ParenthesesOutOfOrder,
MissingComma,
BadFormat // Generic error (e.g., from an out_of_range exception)
};
std::string error_to_string(ReactionParseError err) {
switch (err) {
case ReactionParseError::MissingOpenParenthesis:
return "Missing '('";
case ReactionParseError::MissingCloseParenthesis:
return "Missing ')'";
case ReactionParseError::ParenthesesOutOfOrder:
return "')' found before '('";
case ReactionParseError::MissingComma:
return "Missing ',' within parentheses";
case ReactionParseError::BadFormat:
return "Bad format (substr error)";
}
return "Unknown error";
}
std::expected<std::string, ReactionParseError> reverse_pe_name(const std::string& forwardPEName) noexcept {
try {
size_t pos_open_paren = forwardPEName.find('(');
size_t pos_close_paren = forwardPEName.find(')');
if (pos_open_paren == std::string::npos) {
return std::unexpected(ReactionParseError::MissingOpenParenthesis);
}
if (pos_close_paren == std::string::npos) {
return std::unexpected(ReactionParseError::MissingCloseParenthesis);
}
if (pos_close_paren < pos_open_paren) {
return std::unexpected(ReactionParseError::ParenthesesOutOfOrder);
}
std::string target = forwardPEName.substr(0, pos_open_paren);
std::string result = forwardPEName.substr(pos_close_paren + 1);
std::string inner_content = forwardPEName.substr(
pos_open_paren + 1,
pos_close_paren - pos_open_paren - 1
);
size_t pos_comma = inner_content.find(',');
if (pos_comma == std::string::npos) {
return std::unexpected(ReactionParseError::MissingComma);
}
std::string projectiles = inner_content.substr(0, pos_comma);
std::string ejectiles = inner_content.substr(pos_comma + 1);
std::ostringstream oss;
oss << result << "(" << ejectiles << "," << projectiles << ")" << target;
return oss.str();
} catch (const std::out_of_range&) {
return std::unexpected(ReactionParseError::BadFormat);
}
}
const auto ritr = std::find_if(str.rbegin(), std::string::const_reverse_iterator(startIt),
[](const unsigned char ch){ return !std::isspace(ch); });
return {startIt, ritr.base()};
}
namespace gridfire::reaclib {
static std::unique_ptr<reaction::ReactionSet> s_all_reaclib_reactions_ptr = nullptr;
@@ -106,10 +178,20 @@ namespace gridfire::reaclib {
coeffs[6]
};
auto rpName_revNormalized = std::string(rpName_sv);
if (reverse) {
auto result = reverse_pe_name(std::string(rpName_sv));
if (!result) {
std::string msg = std::format("Error reversing stored projectile-ejectile name for marked reverse reaction ({})", error_to_string(result.error()));
throw exceptions::ReactionParsingError(msg, rpName_revNormalized);
}
rpName_revNormalized = result.value();
}
// Construct the Reaction object. We use rpName for both the unique ID and the human-readable name.
reaction_list.emplace_back(std::make_unique<reaction::ReaclibReaction>(
rpName_sv,
rpName_sv,
rpName_revNormalized,
rpName_revNormalized,
chapter,
reactants,
products,

View File

@@ -40,8 +40,8 @@ namespace gridfire::reaction {
m_peName(peName),
m_chapter(chapter),
m_qValue(qValue),
m_reactants(reactants),
m_products(products),
m_reactants(reactants.begin(), reactants.end()),
m_products(products.begin(), products.end()),
m_sourceLabel(label),
m_rateCoefficients(sets),
m_reverse(reverse) {}
@@ -343,10 +343,12 @@ namespace gridfire::reaction {
return; // Case where the reactions will be added later.
}
m_reactionNameMap.reserve(m_reactions.size());
m_reactionHashes.reserve(m_reactions.size());
size_t i = 0;
for (const auto& reaction : m_reactions) {
m_id += reaction->id();
m_reactionNameMap.emplace(std::string(reaction->id()), i);
m_reactionHashes.insert(reaction->hash(0));
i++;
}
}
@@ -354,11 +356,13 @@ namespace gridfire::reaction {
ReactionSet::ReactionSet(const std::vector<Reaction *> &reactions) {
m_reactions.reserve(reactions.size());
m_reactionNameMap.reserve(reactions.size());
m_reactionHashes.reserve(reactions.size());
size_t i = 0;
for (const auto& reaction : reactions) {
m_reactions.push_back(reaction->clone());
m_id += reaction->id();
m_reactionNameMap.emplace(std::string(reaction->id()), i);
m_reactionHashes.insert(reaction->hash(0));
i++;
}
}
@@ -367,14 +371,13 @@ namespace gridfire::reaction {
ReactionSet::ReactionSet(const ReactionSet &other) {
m_reactions.reserve(other.m_reactions.size());
for (const auto& reaction: other.m_reactions) {
m_reactions.push_back(reaction->clone());
}
m_reactionHashes.reserve(other.m_reactions.size());
m_reactionNameMap.reserve(other.m_reactionNameMap.size());
size_t i = 0;
for (const auto& reaction : m_reactions) {
for (const auto& reaction: other.m_reactions) {
m_reactions.push_back(reaction->clone());
m_reactionNameMap.emplace(std::string(reaction->id()), i);
m_reactionHashes.insert(reaction->hash(0));
i++;
}
}
@@ -397,18 +400,23 @@ namespace gridfire::reaction {
m_id += reaction_id;
m_reactionNameMap.emplace(std::move(reaction_id), new_index);
m_reactionHashes.insert(reaction.hash(0));
}
void ReactionSet::add_reaction(std::unique_ptr<Reaction>&& reaction) {
const std::size_t new_index = m_reactionNameMap.size();
auto reaction_id = std::string(reaction->id());
size_t reaction_hash = reaction->hash(0);
m_reactions.emplace_back(std::move(reaction));
m_id += reaction_id;
m_reactionNameMap.emplace(std::move(reaction_id), new_index);
m_reactionHashes.insert(reaction_hash);
}
void ReactionSet::extend(const ReactionSet &other) {
@@ -425,43 +433,32 @@ namespace gridfire::reaction {
}
void ReactionSet::remove_reaction(const Reaction& reaction) {
const auto reaction_id = std::string(reaction.id());
if (!m_reactionNameMap.contains(reaction_id)) {
const size_t rh = reaction.hash(0);
if (!m_reactionHashes.contains(rh)) {
return;
}
std::erase_if(m_reactions, [&reaction_id](const auto& r_ptr) {
return r_ptr->id() == reaction_id;
std::erase_if(m_reactions, [&rh](const auto& r_ptr) {
return r_ptr->hash(0) == rh;
});
m_reactionNameMap.clear();
m_reactionNameMap.reserve(m_reactions.size());
for (size_t i = 0; i < m_reactions.size(); ++i) {
m_reactionNameMap.emplace(std::string(m_reactions[i]->id()), i);
}
m_reactionNameMap.erase(std::string(reaction.id()));
m_id.clear();
for (const auto& r_ptr : m_reactions) {
m_id += r_ptr->id();
}
m_reactionHashes.erase(rh);
}
bool ReactionSet::contains(const std::string_view& id) const {
for (const auto& reaction : m_reactions) {
if (reaction->id() == id) {
return true;
}
}
return false;
return m_reactionNameMap.contains(std::string(id));
}
bool ReactionSet::contains(const Reaction& reaction) const {
for (const auto& r : m_reactions) {
if (r->id() == reaction.id()) {
return true;
}
}
return false;
const size_t rh = reaction.hash(0);
return m_reactionHashes.contains(rh);
}
void ReactionSet::clear() {
@@ -505,11 +502,12 @@ namespace gridfire::reaction {
}
const Reaction& ReactionSet::operator[](const std::string_view& id) const {
if (const auto it = m_reactionNameMap.find(std::string(id)); it != m_reactionNameMap.end()) {
return *m_reactions[it->second];
if (!contains(id)) {
m_logger -> flush_log();
throw std::out_of_range("Reaction " + std::string(id) + " does not exist in ReactionSet.");
}
m_logger -> flush_log();
throw std::out_of_range("Reaction " + std::string(id) + " does not exist in ReactionSet.");
const size_t idx = m_reactionNameMap.at(std::string(id));
return *m_reactions[idx];
}
bool ReactionSet::operator==(const ReactionSet& other) const {
@@ -523,7 +521,7 @@ namespace gridfire::reaction {
return !(*this == other);
}
uint64_t ReactionSet::hash(uint64_t seed) const {
uint64_t ReactionSet::hash(const uint64_t seed) const {
if (m_reactions.empty()) {
return XXHash64::hash(nullptr, 0, seed);
}

View File

@@ -16,10 +16,12 @@
#include <stdexcept>
#include <algorithm>
#include "fourdst/atomic/species.h"
#include "fourdst/composition/exceptions/exceptions_composition.h"
#include "gridfire/engine/engine_graph.h"
#include "gridfire/solver/strategies/triggers/engine_partitioning_trigger.h"
#include "gridfire/trigger/procedures/trigger_pprint.h"
#include "gridfire/exceptions/error_solver.h"
namespace {
std::unordered_map<int, std::string> cvode_ret_code_map {
@@ -59,9 +61,9 @@ namespace {
void check_cvode_flag(const int flag, const std::string& func_name) {
if (flag < 0) {
if (!cvode_ret_code_map.contains(flag)) {
throw std::runtime_error("CVODE error in " + func_name + ": Unknown error code: " + std::to_string(flag));
throw gridfire::exceptions::CVODESolverFailureError("CVODE error in " + func_name + ": Unknown error code: " + std::to_string(flag));
}
throw std::runtime_error("CVODE error in " + func_name + ": " + cvode_ret_code_map.at(flag));
throw gridfire::exceptions::CVODESolverFailureError("CVODE error in " + func_name + ": " + cvode_ret_code_map.at(flag));
}
}
@@ -91,7 +93,8 @@ namespace gridfire::solver {
const DynamicEngine &engine,
const std::vector<fourdst::atomic::Species> &networkSpecies,
const size_t currentConvergenceFailure,
const size_t currentNonlinearIterations
const size_t currentNonlinearIterations,
const std::map<fourdst::atomic::Species, std::unordered_map<std::string, double>> &reactionContributionMap
) :
t(t),
state(state),
@@ -103,7 +106,8 @@ namespace gridfire::solver {
engine(engine),
networkSpecies(networkSpecies),
currentConvergenceFailures(currentConvergenceFailure),
currentNonlinearIterations(currentNonlinearIterations)
currentNonlinearIterations(currentNonlinearIterations),
reactionContributionMap(reactionContributionMap)
{}
std::vector<std::tuple<std::string, std::string>> CVODESolverStrategy::TimestepContext::describe() const {
@@ -241,6 +245,13 @@ namespace gridfire::solver {
<< "\n";
}
static const std::map<fourdst::atomic::Species,
std::unordered_map<std::string,double>> kEmptyMap{};
const auto& rcMap = user_data.reaction_contribution_map.has_value()
? user_data.reaction_contribution_map.value()
: kEmptyMap;
auto ctx = TimestepContext(
current_time,
m_Y,
@@ -252,7 +263,9 @@ namespace gridfire::solver {
m_engine,
m_engine.getNetworkSpecies(),
convFail_diff,
iter_diff);
iter_diff,
rcMap
);
prev_nonlinear_iterations = nliters + total_nonlinear_iterations;
prev_convergence_failures = nlcfails + total_convergence_failures;
@@ -459,7 +472,7 @@ namespace gridfire::solver {
LOG_TRACE_L2(m_logger, "Constructing final composition= with {} species", numSpecies);
fourdst::composition::Composition topLevelComposition(m_engine.getNetworkSpecies(), y_vec);
fourdst::composition::Composition outputComposition = m_engine.collectComposition(topLevelComposition);
fourdst::composition::Composition outputComposition = m_engine.collectComposition(topLevelComposition, netIn.temperature/1e9, netIn.density);
assert(outputComposition.getRegisteredSymbols().size() == equilibratedComposition.getRegisteredSymbols().size());
@@ -514,7 +527,8 @@ namespace gridfire::solver {
const auto* instance = data->solver_instance;
try {
instance->calculate_rhs(t, y, ydot, data);
const CVODERHSOutputData out = instance->calculate_rhs(t, y, ydot, data);
data->reaction_contribution_map = out.reaction_contribution_map;
return 0;
} catch (const exceptions::StaleEngineTrigger& e) {
data->captured_exception = std::make_unique<exceptions::StaleEngineTrigger>(e);
@@ -564,7 +578,7 @@ namespace gridfire::solver {
return 0;
}
void CVODESolverStrategy::calculate_rhs(
CVODESolverStrategy::CVODERHSOutputData CVODESolverStrategy::calculate_rhs(
const sunrealtype t,
N_Vector y,
N_Vector ydot,
@@ -581,14 +595,6 @@ namespace gridfire::solver {
y_data[i] = 0.0;
}
}
// PERF: The trade off of ensured index consistency is some performance here. If this becomes a bottleneck we can revisit.
// The specific trade off is that we have decided to enforce that all interfaces accept composition objects rather
// than raw vectors of molar abundances. This then lets any method lookup the species by name rather than relying on
// the index in the vector being consistent. The trade off is that sometimes we need to construct a composition object
// which, at the moment, requires a somewhat expensive set of operations. Perhaps in the future we could enforce
// some consistent memory layout for the composition object to make this cheeper. That optimization would need to be
// done in the libcomposition library though...
std::vector<double> y_vec(y_data, y_data + numSpecies);
fourdst::composition::Composition composition(m_engine.getNetworkSpecies(), y_vec);
@@ -598,13 +604,15 @@ namespace gridfire::solver {
}
sunrealtype* ydot_data = N_VGetArrayPointer(ydot);
const auto& [dydt, nuclearEnergyGenerationRate] = result.value();
const auto& [dydt, nuclearEnergyGenerationRate, reactionContributions] = result.value();
for (size_t i = 0; i < numSpecies; ++i) {
fourdst::atomic::Species species = m_engine.getNetworkSpecies()[i];
ydot_data[i] = dydt.at(species);
}
ydot_data[numSpecies] = nuclearEnergyGenerationRate; // Set the last element to the specific energy rate
return {reactionContributions};
}
void CVODESolverStrategy::initialize_cvode_integration_resources(