fix(engine_multiscale): resolved bug which prevented proper equilibrium abundances from being found
this was done by adjusting the scaling of the QSE operator() residuals from r = dy/dt to r=(dy/dt)/y
This commit is contained in:
@@ -60,6 +60,15 @@ namespace gridfire {
|
||||
const fourdst::composition::Composition &comp,
|
||||
const double T9,
|
||||
const double rho
|
||||
) const {
|
||||
return calculateRHSAndEnergy(comp, T9, rho, m_reactions);
|
||||
}
|
||||
|
||||
std::expected<StepDerivatives<double>, expectations::StaleEngineError> GraphEngine::calculateRHSAndEnergy(
|
||||
const fourdst::composition::Composition &comp,
|
||||
const double T9,
|
||||
const double rho,
|
||||
const reaction::ReactionSet &activeReactions
|
||||
) const {
|
||||
const double Ye = comp.getElectronAbundance();
|
||||
const double mue = 5.0e-3 * std::pow(rho * Ye, 1.0 / 3.0) + 0.5 * T9;
|
||||
@@ -73,7 +82,9 @@ namespace gridfire {
|
||||
|
||||
for (const auto& reaction: m_reactions) {
|
||||
bare_rates.push_back(reaction->calculate_rate(T9, rho, Ye, mue, comp.getMolarAbundanceVector(), m_indexToSpeciesMap));
|
||||
bare_reverse_rates.push_back(calculateReverseRate(*reaction, T9, rho, comp));
|
||||
if (reaction->type() != reaction::ReactionType::WEAK) {
|
||||
bare_reverse_rates.push_back(calculateReverseRate(*reaction, T9, rho, comp));
|
||||
}
|
||||
}
|
||||
|
||||
// --- The public facing interface can always use the precomputed version since taping is done internally ---
|
||||
@@ -90,6 +101,10 @@ namespace gridfire {
|
||||
return comp.getSpeciesIndex(species); // Return the index of the species in the composition
|
||||
}
|
||||
return std::nullopt; // Species not found in the composition
|
||||
},
|
||||
[&activeReactions](const reaction::Reaction& reaction) -> bool {
|
||||
if (activeReactions.contains(reaction)) { return true; }
|
||||
return false;
|
||||
}
|
||||
);
|
||||
}
|
||||
@@ -99,6 +114,15 @@ namespace gridfire {
|
||||
const fourdst::composition::Composition &comp,
|
||||
const double T9,
|
||||
const double rho
|
||||
) const {
|
||||
return calculateEpsDerivatives(comp, T9, rho, m_reactions);
|
||||
}
|
||||
|
||||
EnergyDerivatives GraphEngine::calculateEpsDerivatives(
|
||||
const fourdst::composition::Composition &comp,
|
||||
const double T9,
|
||||
const double rho,
|
||||
const reaction::ReactionSet &activeReactions
|
||||
) const {
|
||||
const size_t numSpecies = m_networkSpecies.size();
|
||||
const size_t numADInputs = numSpecies + 2; // +2 for T9 and rho
|
||||
@@ -192,6 +216,8 @@ 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 teh 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
|
||||
});
|
||||
@@ -434,14 +460,13 @@ namespace gridfire {
|
||||
const fourdst::composition::Composition& comp,
|
||||
const double reverseRate
|
||||
) const {
|
||||
assert(reaction.type() == reaction::ReactionType::LOGICAL_REACLIB || reaction.type() == reaction::ReactionType::REACLIB);
|
||||
|
||||
if (!m_useReverseReactions) {
|
||||
LOG_TRACE_L3_LIMIT_EVERY_N(std::numeric_limits<int>::max(), m_logger, "Reverse reactions are disabled. Returning 0.0 for reverse rate of reaction '{}'.", reaction.id());
|
||||
return 0.0; // If reverse reactions are not used, return 0.0
|
||||
}
|
||||
double Ye = comp.getElectronAbundance();
|
||||
// TODO: This is a dummy value for the electron chemical potential. We eventually need to replace this with an EOS call.
|
||||
double mue = 5.0e-3 * std::pow(rho * Ye, 1.0 / 3.0) + 0.5 * T9;
|
||||
const double d_log_kFwd = reaction.calculate_log_rate_partial_deriv_wrt_T9(T9, rho, Ye, mue, comp);
|
||||
const double d_log_kFwd = reaction.calculate_log_rate_partial_deriv_wrt_T9(T9, rho, {}, {}, {});
|
||||
|
||||
auto log_deriv_pf_op = [&](double acc, const auto& species) {
|
||||
const double g = m_partitionFunction->evaluate(species.z(), species.a(), T9);
|
||||
@@ -947,6 +972,15 @@ namespace gridfire {
|
||||
const fourdst::composition::Composition &comp,
|
||||
const double T9,
|
||||
const double rho
|
||||
) const {
|
||||
return getSpeciesTimescales(comp, T9, rho, m_reactions);
|
||||
}
|
||||
|
||||
std::expected<std::unordered_map<fourdst::atomic::Species, double>, expectations::StaleEngineError> GraphEngine::getSpeciesTimescales(
|
||||
const fourdst::composition::Composition &comp,
|
||||
const double T9,
|
||||
const double rho,
|
||||
const reaction::ReactionSet &activeReactions
|
||||
) const {
|
||||
const double Ye = comp.getElectronAbundance();
|
||||
|
||||
@@ -961,6 +995,9 @@ namespace gridfire {
|
||||
return comp.getSpeciesIndex(species);
|
||||
}
|
||||
return std::nullopt; // Species not present
|
||||
},
|
||||
[&activeReactions](const reaction::Reaction& reaction) -> bool {
|
||||
return activeReactions.contains(reaction);
|
||||
}
|
||||
);
|
||||
std::unordered_map<fourdst::atomic::Species, double> speciesTimescales;
|
||||
@@ -979,6 +1016,15 @@ namespace gridfire {
|
||||
const fourdst::composition::Composition &comp,
|
||||
const double T9,
|
||||
const double rho
|
||||
) const {
|
||||
return getSpeciesDestructionTimescales(comp, T9, rho, m_reactions);
|
||||
}
|
||||
|
||||
std::expected<std::unordered_map<fourdst::atomic::Species, double>, expectations::StaleEngineError> GraphEngine::getSpeciesDestructionTimescales(
|
||||
const fourdst::composition::Composition &comp,
|
||||
const double T9,
|
||||
const double rho,
|
||||
const reaction::ReactionSet &activeReactions
|
||||
) const {
|
||||
const double Ye = comp.getElectronAbundance();
|
||||
const std::vector<double>& Y = comp.getMolarAbundanceVector();
|
||||
@@ -996,7 +1042,10 @@ namespace gridfire {
|
||||
rho,
|
||||
Ye,
|
||||
0.0,
|
||||
speciesLookup
|
||||
speciesLookup,
|
||||
[&activeReactions](const reaction::Reaction& reaction) -> bool {
|
||||
return activeReactions.contains(reaction);
|
||||
}
|
||||
);
|
||||
|
||||
std::unordered_map<fourdst::atomic::Species, double> speciesDestructionTimescales;
|
||||
@@ -1095,6 +1144,9 @@ namespace gridfire {
|
||||
adMue,
|
||||
[&](const fourdst::atomic::Species& querySpecies) -> size_t {
|
||||
return m_speciesToIndexMap.at(querySpecies);
|
||||
},
|
||||
[](const reaction::Reaction& reaction) -> bool {
|
||||
return true; // Use all reactions
|
||||
}
|
||||
);
|
||||
|
||||
@@ -1102,6 +1154,9 @@ namespace gridfire {
|
||||
// Extract the raw vector from the associative map
|
||||
std::vector<CppAD::AD<double>> dydt_vec;
|
||||
dydt_vec.reserve(dydt.size());
|
||||
// TODO: There is a possibility for a bug here if the map ordering is not consistent with the ordering of the species in m_networkSpecies.
|
||||
// right now this works but that's because I am careful to build the map in the right order. This should be made less fragile
|
||||
// so that if map construction order changes this still works.
|
||||
std::ranges::transform(dydt, std::back_inserter(dydt_vec),[](const auto& kv) { return kv.second; });
|
||||
|
||||
m_rhsADFun.Dependent(adInput, dydt_vec);
|
||||
@@ -1147,6 +1202,9 @@ namespace gridfire {
|
||||
adMue,
|
||||
[&](const fourdst::atomic::Species& querySpecies) -> size_t {
|
||||
return m_speciesToIndexMap.at(querySpecies); // TODO: This is bad, needs to be fixed
|
||||
},
|
||||
[](const reaction::Reaction& reaction) -> bool {
|
||||
return true; // Use all reactions
|
||||
}
|
||||
);
|
||||
|
||||
@@ -1188,6 +1246,7 @@ namespace gridfire {
|
||||
const auto& reaction = m_reactions[i];
|
||||
PrecomputedReaction precomp;
|
||||
precomp.reaction_index = i;
|
||||
precomp.reaction_type = reaction.type();
|
||||
|
||||
// --- Precompute forward reaction information ---
|
||||
// Count occurrences for each reactant to determine powers and symmetry
|
||||
@@ -1208,7 +1267,7 @@ namespace gridfire {
|
||||
precomp.symmetry_factor = 1.0/symmetryDenominator;
|
||||
|
||||
// --- Precompute reverse reaction information ---
|
||||
if (reaction.qValue() != 0.0) {
|
||||
if (reaction.qValue() != 0.0 && reaction.type() != reaction::ReactionType::WEAK) {
|
||||
std::unordered_map<size_t, int> productCounts;
|
||||
for (const auto& product : reaction.products()) {
|
||||
productCounts[speciesIndexMap.at(product)]++;
|
||||
@@ -1224,7 +1283,7 @@ namespace gridfire {
|
||||
} else {
|
||||
precomp.unique_product_indices.clear();
|
||||
precomp.product_powers.clear();
|
||||
precomp.reverse_symmetry_factor = 0.0; // No reverse reaction for Q = 0 reactions
|
||||
precomp.reverse_symmetry_factor = 0.0; // No reverse reaction for weak reactions
|
||||
}
|
||||
|
||||
// --- Precompute stoichiometry information ---
|
||||
|
||||
@@ -78,6 +78,8 @@ namespace gridfire {
|
||||
return dominateReaction;
|
||||
}
|
||||
|
||||
|
||||
|
||||
/**
|
||||
* @brief Primes absent species in the network to their equilibrium abundances using a robust, two-stage approach.
|
||||
*
|
||||
@@ -112,7 +114,7 @@ namespace gridfire {
|
||||
*/
|
||||
PrimingReport primeNetwork(
|
||||
const NetIn& netIn,
|
||||
DynamicEngine& engine,
|
||||
GraphEngine& engine,
|
||||
const std::optional<std::vector<reaction::ReactionType>>& ignoredReactionTypes
|
||||
) {
|
||||
auto logger = fourdst::logging::LogManager::getInstance().getLogger("log");
|
||||
|
||||
@@ -91,7 +91,11 @@ namespace gridfire {
|
||||
|
||||
updatedNetIn.composition = baseUpdatedComposition;
|
||||
|
||||
updatedNetIn.composition.finalize(false);
|
||||
bool didFinalize = updatedNetIn.composition.finalize(false);
|
||||
if (!didFinalize) {
|
||||
LOG_ERROR(m_logger, "Failed to finalize composition during adaptive engine view update. Check input mass fractions for validity.");
|
||||
throw std::runtime_error("Failed to finalize composition during adaptive engine view update.");
|
||||
}
|
||||
|
||||
LOG_TRACE_L1(m_logger, "Updating AdaptiveEngineView with new network input...");
|
||||
|
||||
|
||||
@@ -1,4 +1,5 @@
|
||||
#include "gridfire/engine/views/engine_defined.h"
|
||||
#include "gridfire/engine/engine_graph.h"
|
||||
|
||||
#include <ranges>
|
||||
|
||||
@@ -7,14 +8,33 @@
|
||||
#include <string>
|
||||
#include <vector>
|
||||
#include <unordered_set>
|
||||
#include <set>
|
||||
#include <stdexcept>
|
||||
#include <unordered_map>
|
||||
#include <utility>
|
||||
|
||||
namespace {
|
||||
class MaskedComposition final : public fourdst::composition::Composition {
|
||||
private:
|
||||
std::set<fourdst::atomic::Species> m_activeSpecies;
|
||||
public:
|
||||
MaskedComposition(const Composition& baseComposition, const std::set<fourdst::atomic::Species>& activeSpecies) :
|
||||
Composition(baseComposition),
|
||||
m_activeSpecies(activeSpecies) {}
|
||||
|
||||
bool contains(const fourdst::atomic::Species& species) const override {
|
||||
return Composition::contains(species) && m_activeSpecies.contains(species);
|
||||
}
|
||||
};
|
||||
}
|
||||
|
||||
namespace gridfire {
|
||||
using fourdst::atomic::Species;
|
||||
|
||||
DefinedEngineView::DefinedEngineView(const std::vector<std::string>& peNames, DynamicEngine& baseEngine) :
|
||||
DefinedEngineView::DefinedEngineView(
|
||||
const std::vector<std::string>& peNames,
|
||||
GraphEngine& baseEngine
|
||||
) :
|
||||
m_baseEngine(baseEngine) {
|
||||
collect(peNames);
|
||||
}
|
||||
@@ -24,7 +44,11 @@ namespace gridfire {
|
||||
}
|
||||
|
||||
const std::vector<Species> & DefinedEngineView::getNetworkSpecies() const {
|
||||
return m_activeSpecies;
|
||||
if (m_activeSpeciesVectorCache.has_value()) {
|
||||
return m_activeSpeciesVectorCache.value();
|
||||
}
|
||||
m_activeSpeciesVectorCache = std::vector<Species>(m_activeSpecies.begin(), m_activeSpecies.end());
|
||||
return m_activeSpeciesVectorCache.value();
|
||||
}
|
||||
|
||||
std::expected<StepDerivatives<double>, expectations::StaleEngineError> DefinedEngineView::calculateRHSAndEnergy(
|
||||
@@ -34,7 +58,8 @@ namespace gridfire {
|
||||
) const {
|
||||
validateNetworkState();
|
||||
|
||||
const auto result = m_baseEngine.calculateRHSAndEnergy(comp, T9, rho);
|
||||
const MaskedComposition masked(comp, m_activeSpecies);
|
||||
const auto result = m_baseEngine.calculateRHSAndEnergy(masked, T9, rho, m_activeReactions);
|
||||
|
||||
if (!result) {
|
||||
return std::unexpected{result.error()};
|
||||
@@ -50,7 +75,9 @@ namespace gridfire {
|
||||
) const {
|
||||
validateNetworkState();
|
||||
|
||||
return m_baseEngine.calculateEpsDerivatives(comp, T9, rho);
|
||||
const MaskedComposition masked(comp, m_activeSpecies);
|
||||
|
||||
return m_baseEngine.calculateEpsDerivatives(masked, T9, rho, m_activeReactions);
|
||||
}
|
||||
|
||||
void DefinedEngineView::generateJacobianMatrix(
|
||||
@@ -60,7 +87,10 @@ namespace gridfire {
|
||||
) const {
|
||||
validateNetworkState();
|
||||
|
||||
m_baseEngine.generateJacobianMatrix(comp, T9, rho);
|
||||
const MaskedComposition masked(comp, m_activeSpecies);
|
||||
|
||||
// TODO: We likely want to be able to think more carefully about this so that the jacobian matches the active species/reactions
|
||||
m_baseEngine.generateJacobianMatrix(masked, T9, rho);
|
||||
}
|
||||
|
||||
double DefinedEngineView::getJacobianMatrixEntry(
|
||||
@@ -69,6 +99,17 @@ namespace gridfire {
|
||||
) const {
|
||||
validateNetworkState();
|
||||
|
||||
if (!m_activeSpecies.contains(rowSpecies)) {
|
||||
LOG_ERROR(m_logger, "Row species '{}' is not part of the active species in the DefinedEngineView.", rowSpecies.name());
|
||||
m_logger -> flush_log();
|
||||
throw std::runtime_error("Row species not found in active species: " + std::string(rowSpecies.name()));
|
||||
}
|
||||
if (!m_activeSpecies.contains(colSpecies)) {
|
||||
LOG_ERROR(m_logger, "Column species '{}' is not part of the active species in the DefinedEngineView.", colSpecies.name());
|
||||
m_logger -> flush_log();
|
||||
throw std::runtime_error("Column species not found in active species: " + std::string(colSpecies.name()));
|
||||
}
|
||||
|
||||
return m_baseEngine.getJacobianMatrixEntry(rowSpecies, colSpecies);
|
||||
}
|
||||
|
||||
@@ -84,6 +125,18 @@ namespace gridfire {
|
||||
) const {
|
||||
validateNetworkState();
|
||||
|
||||
if (!m_activeSpecies.contains(species)) {
|
||||
LOG_ERROR(m_logger, "Species '{}' is not part of the active species in the DefinedEngineView.", species.name());
|
||||
m_logger -> flush_log();
|
||||
throw std::runtime_error("Species not found in active species: " + std::string(species.name()));
|
||||
}
|
||||
|
||||
if (!m_activeReactions.contains(reaction)) {
|
||||
LOG_ERROR(m_logger, "Reaction '{}' is not part of the active reactions in the DefinedEngineView.", reaction.id());
|
||||
m_logger -> flush_log();
|
||||
throw std::runtime_error("Reaction not found in active reactions: " + std::string(reaction.id()));
|
||||
}
|
||||
|
||||
return m_baseEngine.getStoichiometryMatrixEntry(species, reaction);
|
||||
}
|
||||
|
||||
@@ -100,7 +153,9 @@ namespace gridfire {
|
||||
m_logger -> flush_log();
|
||||
throw std::runtime_error("Reaction not found in active reactions: " + std::string(reaction.id()));
|
||||
}
|
||||
return m_baseEngine.calculateMolarReactionFlow(reaction, comp, T9, rho);
|
||||
|
||||
const MaskedComposition masked(comp, m_activeSpecies);
|
||||
return m_baseEngine.calculateMolarReactionFlow(reaction, masked, T9, rho);
|
||||
}
|
||||
|
||||
const reaction::ReactionSet & DefinedEngineView::getNetworkReactions() const {
|
||||
@@ -115,6 +170,7 @@ namespace gridfire {
|
||||
peNames.emplace_back(reaction->id());
|
||||
}
|
||||
collect(peNames);
|
||||
m_activeSpeciesVectorCache = std::nullopt; // Invalidate species vector cache
|
||||
}
|
||||
|
||||
std::expected<std::unordered_map<Species, double>, expectations::StaleEngineError> DefinedEngineView::getSpeciesTimescales(
|
||||
@@ -123,8 +179,9 @@ namespace gridfire {
|
||||
const double rho
|
||||
) const {
|
||||
validateNetworkState();
|
||||
const MaskedComposition masked(comp, m_activeSpecies);
|
||||
|
||||
const auto result = m_baseEngine.getSpeciesTimescales(comp, T9, rho);
|
||||
const auto result = m_baseEngine.getSpeciesTimescales(masked, T9, rho, m_activeReactions);
|
||||
if (!result) {
|
||||
return std::unexpected{result.error()};
|
||||
}
|
||||
@@ -139,15 +196,15 @@ namespace gridfire {
|
||||
return definedTimescales;
|
||||
}
|
||||
|
||||
std::expected<std::unordered_map<fourdst::atomic::Species, double>, expectations::StaleEngineError>
|
||||
DefinedEngineView::getSpeciesDestructionTimescales(
|
||||
std::expected<std::unordered_map<Species, double>, expectations::StaleEngineError> DefinedEngineView::getSpeciesDestructionTimescales(
|
||||
const fourdst::composition::Composition &comp,
|
||||
const double T9,
|
||||
const double rho
|
||||
) const {
|
||||
validateNetworkState();
|
||||
const MaskedComposition masked(comp, m_activeSpecies);
|
||||
|
||||
const auto result = m_baseEngine.getSpeciesDestructionTimescales(comp, T9, rho);
|
||||
const auto result = m_baseEngine.getSpeciesDestructionTimescales(masked, T9, rho, m_activeReactions);
|
||||
|
||||
if (!result) {
|
||||
return std::unexpected{result.error()};
|
||||
@@ -182,6 +239,7 @@ namespace gridfire {
|
||||
}
|
||||
|
||||
size_t DefinedEngineView::getSpeciesIndex(const Species &species) const {
|
||||
// TODO: We are working to phase out all of these methods, its probably broken but it also should no longer be used and will be removed soon
|
||||
validateNetworkState();
|
||||
|
||||
const auto it = std::ranges::find(m_activeSpecies, species);
|
||||
@@ -328,13 +386,13 @@ namespace gridfire {
|
||||
for (const auto& reactant : reaction->reactants()) {
|
||||
if (!seenSpecies.contains(reactant)) {
|
||||
seenSpecies.insert(reactant);
|
||||
m_activeSpecies.push_back(reactant);
|
||||
m_activeSpecies.emplace(reactant);
|
||||
}
|
||||
}
|
||||
for (const auto& product : reaction->products()) {
|
||||
if (!seenSpecies.contains(product)) {
|
||||
seenSpecies.insert(product);
|
||||
m_activeSpecies.push_back(product);
|
||||
m_activeSpecies.emplace(product);
|
||||
}
|
||||
}
|
||||
m_activeReactions.add_reaction(*reaction);
|
||||
@@ -373,7 +431,7 @@ namespace gridfire {
|
||||
/////////////////////////////////////////////
|
||||
|
||||
FileDefinedEngineView::FileDefinedEngineView(
|
||||
DynamicEngine &baseEngine,
|
||||
GraphEngine &baseEngine,
|
||||
const std::string &fileName,
|
||||
const io::NetworkFileParser &parser
|
||||
):
|
||||
|
||||
@@ -1,12 +1,14 @@
|
||||
#include "gridfire/engine/views/engine_multiscale.h"
|
||||
#include "gridfire/exceptions/error_engine.h"
|
||||
#include "gridfire/engine/procedures/priming.h"
|
||||
#include "gridfire/utils/general_composition.h"
|
||||
|
||||
#include <stdexcept>
|
||||
#include <vector>
|
||||
#include <ranges>
|
||||
#include <unordered_map>
|
||||
#include <unordered_set>
|
||||
#include <fstream>
|
||||
|
||||
#include <queue>
|
||||
|
||||
@@ -15,6 +17,8 @@
|
||||
#include "quill/LogMacros.h"
|
||||
#include "quill/Logger.h"
|
||||
|
||||
static std::ofstream debug_multiscale_log("debug_multiscale_log.txt");
|
||||
|
||||
namespace {
|
||||
using namespace fourdst::atomic;
|
||||
//TODO: Replace all calls to this function with composition.getMolarAbundanceVector() so that
|
||||
@@ -247,16 +251,30 @@ namespace gridfire {
|
||||
) const {
|
||||
// Fix the algebraic species to the equilibrium abundances we calculate.
|
||||
fourdst::composition::Composition comp_mutable = comp;
|
||||
for (const auto& species : m_algebraic_species) {
|
||||
// TODO: Check this conversion to mass fraction (also consider adding the ability to set molar abundance directly)
|
||||
const double Yi = m_algebraic_abundances.at(species);
|
||||
comp_mutable.setMassFraction(species, Yi * species.a() / (rho * 1e-3)); // Convert Yi (mol/g) to Xi (mass fraction)
|
||||
const bool didFinalize = comp_mutable.finalize(false);
|
||||
if (!didFinalize) {
|
||||
LOG_ERROR(m_logger, "Failed to finalize composition before setting algebraic species abundances.");
|
||||
m_logger->flush_log();
|
||||
throw std::runtime_error("Failed to finalize composition before setting algebraic species abundances.");
|
||||
}
|
||||
|
||||
for (const auto& species : m_algebraic_species) {
|
||||
const double Yi = m_algebraic_abundances.at(species);
|
||||
double Xi = utils::massFractionFromMolarAbundance(comp_mutable, species, Yi);
|
||||
comp_mutable.setMassFraction(species, Xi); // Convert Yi (mol/g) to Xi (mass fraction)
|
||||
if (!comp_mutable.finalize(false)) {
|
||||
LOG_ERROR(m_logger, "Failed to finalize composition after setting algebraic species abundance for species '{}'.", species.name());
|
||||
m_logger->flush_log();
|
||||
throw std::runtime_error("Failed to finalize composition after setting algebraic species abundance for species: " + std::string(species.name()));
|
||||
}
|
||||
}
|
||||
|
||||
if (!comp_mutable.finalize()) {
|
||||
LOG_ERROR(m_logger, "Failed to finalize composition after setting algebraic species abundances.");
|
||||
m_logger->flush_log();
|
||||
throw std::runtime_error("Failed to finalize composition after setting algebraic species abundances.");
|
||||
}
|
||||
|
||||
return m_baseEngine.calculateMolarReactionFlow(reaction, comp_mutable, T9, rho);
|
||||
}
|
||||
|
||||
@@ -285,7 +303,7 @@ namespace gridfire {
|
||||
return speciesTimescales;
|
||||
}
|
||||
|
||||
std::expected<std::unordered_map<fourdst::atomic::Species, double>, expectations::StaleEngineError>
|
||||
std::expected<std::unordered_map<Species, double>, expectations::StaleEngineError>
|
||||
MultiscalePartitioningEngineView::getSpeciesDestructionTimescales(
|
||||
const fourdst::composition::Composition &comp,
|
||||
const double T9,
|
||||
@@ -854,13 +872,15 @@ namespace gridfire {
|
||||
|
||||
constexpr double ABSOLUTE_QSE_TIMESCALE_THRESHOLD = 3.156e7; // Absolute threshold for QSE timescale (1 yr)
|
||||
constexpr double MIN_GAP_THRESHOLD = 2.0; // Require a 2 order of magnitude gap
|
||||
constexpr double MIN_MOLAR_ABUNDANCE_THRESHOLD = 1e-10; // Minimum abundance threshold to consider a species for QSE. Any species above this will always be considered dynamic.
|
||||
constexpr double MAX_MOLAR_ABUNDANCE_THRESHOLD = 1e-10; // Maximum molar abundance which a fast species is allowed to have (anything more abundant is always considered dynamic)
|
||||
constexpr double MIN_MOLAR_ABUNDANCE_THRESHOLD = 1e-50; // Minimum molar abundance to consider a species at all (anything less abundance will be classed as dynamic but with the intent that some latter view will deal with it)
|
||||
|
||||
LOG_TRACE_L1(m_logger, "Found {} species with finite timescales.", sorted_destruction_timescales.size());
|
||||
LOG_TRACE_L1(m_logger, "Absolute QSE timescale threshold: {} seconds ({} years).",
|
||||
ABSOLUTE_QSE_TIMESCALE_THRESHOLD, ABSOLUTE_QSE_TIMESCALE_THRESHOLD / 3.156e7);
|
||||
LOG_TRACE_L1(m_logger, "Minimum gap threshold: {} orders of magnitude.", MIN_GAP_THRESHOLD);
|
||||
LOG_TRACE_L1(m_logger, "Minimum molar abundance threshold: {}.", MIN_MOLAR_ABUNDANCE_THRESHOLD);
|
||||
LOG_TRACE_L1(m_logger, "Maximum molar abundance threshold for fast species consideration : {}.", MAX_MOLAR_ABUNDANCE_THRESHOLD);
|
||||
LOG_TRACE_L1(m_logger, "Minimum molar abundance threshold for species consideration : {}.", MIN_MOLAR_ABUNDANCE_THRESHOLD);
|
||||
|
||||
std::vector<Species> dynamic_pool_species;
|
||||
std::vector<std::pair<double, Species>> fast_candidates;
|
||||
@@ -878,8 +898,14 @@ namespace gridfire {
|
||||
dynamic_pool_species.push_back(species);
|
||||
} else {
|
||||
const double Yi = comp.getMolarAbundance(species);
|
||||
if (Yi > MIN_MOLAR_ABUNDANCE_THRESHOLD) {
|
||||
if (Yi > MAX_MOLAR_ABUNDANCE_THRESHOLD) {
|
||||
LOG_TRACE_L3(m_logger, "Species {} with abundance {} is considered dynamic (above minimum abundance threshold of {}).",
|
||||
species.name(), Yi, MAX_MOLAR_ABUNDANCE_THRESHOLD);
|
||||
dynamic_pool_species.push_back(species);
|
||||
continue;
|
||||
}
|
||||
if (Yi < MIN_MOLAR_ABUNDANCE_THRESHOLD) {
|
||||
LOG_TRACE_L3(m_logger, "Species {} with abundance {} is considered dynamic (below minimum abundance threshold of {}). Likely another network view (such as adaptive engine view) will be needed to deal with this species",
|
||||
species.name(), Yi, MIN_MOLAR_ABUNDANCE_THRESHOLD);
|
||||
dynamic_pool_species.push_back(species);
|
||||
continue;
|
||||
@@ -1057,10 +1083,7 @@ namespace gridfire {
|
||||
|
||||
}
|
||||
|
||||
bool group_is_coupled = (coupling_flux / leakage_flux) > FLUX_RATIO_THRESHOLD;
|
||||
bool group_is_balanced = std::abs(std::log(creationFlux) - std::log(destructionFlux)) < LOG_FLOW_RATIO_THRESHOLD;
|
||||
|
||||
if (group_is_coupled) {
|
||||
if (bool group_is_coupled = (coupling_flux / leakage_flux) > FLUX_RATIO_THRESHOLD) {
|
||||
LOG_TRACE_L1(
|
||||
m_logger,
|
||||
"Group containing {} is in equilibrium due to high coupling flux and balanced creation and destruction rate: <coupling: leakage flux = {}, coupling flux = {}, ratio = {} (Threshold: {})>, <creation: creation flux = {}, destruction flux = {}, ratio = {} order of mag (Threshold: {} order of mag)>",
|
||||
@@ -1160,10 +1183,10 @@ namespace gridfire {
|
||||
long i = 0;
|
||||
std::unordered_map<Species, size_t> species_to_index_map;
|
||||
for (const auto& species : algebraic_species) {
|
||||
constexpr double abundance_floor = 1.0e-15;
|
||||
constexpr double abundance_floor = 1.0e-100;
|
||||
const double initial_abundance = normalized_composition.getMolarAbundance(species);
|
||||
Y_scale(i) = std::max(initial_abundance, abundance_floor);
|
||||
v_initial(i) = std::asinh(initial_abundance / Y_scale(i)); // Scale the initial abundances using asinh
|
||||
const double Y = std::max(initial_abundance, abundance_floor);
|
||||
v_initial(i) = std::log(Y);
|
||||
species_to_index_map.emplace(species, i);
|
||||
i++;
|
||||
}
|
||||
@@ -1194,7 +1217,7 @@ namespace gridfire {
|
||||
throw std::runtime_error(msg.str());
|
||||
}
|
||||
LOG_TRACE_L1(m_logger, "QSE Group minimization succeeded with status: {}", lm_status_map.at(status));
|
||||
Eigen::VectorXd Y_final_qse = Y_scale.array() * v_initial.array().sinh(); // Convert back to physical abundances using asinh scaling
|
||||
Eigen::VectorXd Y_final_qse = v_initial.array().exp(); // Convert back to physical abundances using exponential scaling
|
||||
i = 0;
|
||||
for (const auto& species: algebraic_species) {
|
||||
LOG_TRACE_L1(
|
||||
@@ -1204,8 +1227,8 @@ namespace gridfire {
|
||||
normalized_composition.getMolarAbundance(species),
|
||||
Y_final_qse(i)
|
||||
);
|
||||
//TODO: Check this conversion
|
||||
double Xi = Y_final_qse(i) * species.mass(); // Convert from molar abundance to mass fraction
|
||||
// double Xi = Y_final_qse(i) * species.mass(); // Convert from molar abundance to mass fraction
|
||||
double Xi = utils::massFractionFromMolarAbundance(normalized_composition, species, Y_final_qse(i));
|
||||
if (!outputComposition.hasSpecies(species)) {
|
||||
outputComposition.registerSpecies(species);
|
||||
}
|
||||
@@ -1423,24 +1446,25 @@ namespace gridfire {
|
||||
|
||||
int MultiscalePartitioningEngineView::EigenFunctor::operator()(const InputType &v_qse, OutputType &f_qse) const {
|
||||
fourdst::composition::Composition comp_trial = m_initial_comp;
|
||||
Eigen::VectorXd y_qse = m_Y_scale.array() * v_qse.array().sinh(); // Convert to physical abundances using asinh scaling
|
||||
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.hasSymbol(std::string(species.name()))) {
|
||||
comp_trial.registerSpecies(species);
|
||||
}
|
||||
const double molarAbundance = y_qse[static_cast<long>(m_qse_solve_species_index_map.at(species))];
|
||||
double massFraction = molarAbundance * species.mass();
|
||||
if (massFraction < 0 && std::abs(massFraction) < 1e-20) { // if there is a larger negative mass fraction, let the composition module throw an error
|
||||
massFraction = 0.0; // Avoid setting minuscule negative mass fractions due to numerical noise
|
||||
}
|
||||
auto index = static_cast<long>(m_qse_solve_species_index_map.at(species));
|
||||
const double molarAbundance = y_qse[index];
|
||||
double massFraction = utils::massFractionFromMolarAbundance(m_initial_comp, species, molarAbundance);
|
||||
comp_trial.setMassFraction(species, massFraction);
|
||||
}
|
||||
|
||||
|
||||
const bool didFinalize = comp_trial.finalize(false);
|
||||
if (!didFinalize) {
|
||||
std::string msg = std::format("Failed to finalize composition (comp_trial) in {} at line {}", __FILE__, __LINE__);
|
||||
throw std::runtime_error(msg);
|
||||
LOG_TRACE_L1(m_view->m_logger, "While evaluating the functor, failed to finalize composition. This is likely because the solver took a step outside of physical abundances. This is not an error; rather, the solver will be told to take a different step.");
|
||||
f_qse.resize(static_cast<long>(m_qse_solve_species.size()));
|
||||
f_qse.setConstant(1.0e20); // Return a large residual to indicate failure
|
||||
return 0;
|
||||
}
|
||||
|
||||
const auto result = m_view->getBaseEngine().calculateRHSAndEnergy(comp_trial, m_T9, m_rho);
|
||||
@@ -1452,33 +1476,69 @@ namespace gridfire {
|
||||
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
|
||||
for (const auto& species: m_qse_solve_species) {
|
||||
f_qse(i) = dydt.at(species);
|
||||
const double dydt_i = dydt.at(species);
|
||||
f_qse(i) = dydt_i/y_qse(i); // We square the residuals to improve numerical stability in the solver
|
||||
i++;
|
||||
}
|
||||
|
||||
return 0; // Success
|
||||
LOG_TRACE_L2(
|
||||
m_view->m_logger,
|
||||
"Functor evaluation at T9 = {}, rho = {}, y_qse = <{}> complete. ||f|| = {}",
|
||||
m_T9,
|
||||
m_rho,
|
||||
[&]() -> std::string {
|
||||
std::stringstream ss;
|
||||
for (long j = 0; j < y_qse.size(); ++j) {
|
||||
ss << y_qse(j);
|
||||
if (j < y_qse.size() - 1) {
|
||||
ss << ", ";
|
||||
}
|
||||
}
|
||||
return ss.str();
|
||||
}(),
|
||||
f_qse.norm()
|
||||
);
|
||||
LOG_TRACE_L3(
|
||||
m_view->m_logger,
|
||||
"{}",
|
||||
[&]() -> std::string {
|
||||
std::stringstream ss;
|
||||
const std::vector species(m_qse_solve_species.begin(), m_qse_solve_species.end());
|
||||
for (long j = 0; j < f_qse.size(); ++j) {
|
||||
ss << "Residual for species " << species.at(j).name() << " f(" << j << ") = " << f_qse(j) << "\n";
|
||||
}
|
||||
return ss.str();
|
||||
}()
|
||||
);
|
||||
return 0;
|
||||
}
|
||||
|
||||
int MultiscalePartitioningEngineView::EigenFunctor::df(const InputType &v_qse, JacobianType &J_qse) const {
|
||||
fourdst::composition::Composition comp_trial = m_initial_comp;
|
||||
Eigen::VectorXd y_qse = m_Y_scale.array() * v_qse.array().sinh(); // Convert to physical abundances using asinh scaling
|
||||
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.hasSymbol(std::string(species.name()))) {
|
||||
comp_trial.registerSpecies(species);
|
||||
}
|
||||
const double molarAbundance = y_qse[static_cast<long>(m_qse_solve_species_index_map.at(species))];
|
||||
const double massFraction = molarAbundance * species.mass();
|
||||
double massFraction = utils::massFractionFromMolarAbundance(m_initial_comp, species, molarAbundance);
|
||||
comp_trial.setMassFraction(species, massFraction);
|
||||
}
|
||||
|
||||
const bool didFinalize = comp_trial.finalize(false);
|
||||
if (!didFinalize) {
|
||||
const std::string msg = std::format("Failed to finalize composition (comp_trial) in {} at line {}", __FILE__, __LINE__);
|
||||
throw std::runtime_error(msg);
|
||||
LOG_TRACE_L1(m_view->m_logger, "While evaluating the Jacobian, failed to finalize composition. This is likely because the solver took a step outside of physical abundances. This is not an error; rather, the solver will be told to take a different step. Returning Identity");
|
||||
J_qse.resize(static_cast<long>(m_qse_solve_species.size()), static_cast<long>(m_qse_solve_species.size()));
|
||||
J_qse.setIdentity();
|
||||
return 0;
|
||||
}
|
||||
|
||||
m_view->getBaseEngine().generateJacobianMatrix(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 long N = static_cast<long>(m_qse_solve_species.size());
|
||||
J_qse.resize(N, N);
|
||||
@@ -1496,11 +1556,18 @@ namespace gridfire {
|
||||
rowID += 1;
|
||||
}
|
||||
|
||||
// Chain rule for asinh scaling:
|
||||
for (long j = 0; j < J_qse.cols(); ++j) {
|
||||
const double dY_dv = m_Y_scale(j) * std::cosh(v_qse(j));
|
||||
J_qse.col(j) *= dY_dv; // Scale the column by the derivative of the asinh scaling
|
||||
for (long i = 0; i < J_qse.rows(); ++i) {
|
||||
for (long j = 0; j < J_qse.cols(); ++j) {
|
||||
double on_diag_correction = 0.0;
|
||||
if (i == j) {
|
||||
auto rowSpecies = *(std::next(m_qse_solve_species.begin(), i));
|
||||
const double Fi = dydt.at(rowSpecies);
|
||||
on_diag_correction = Fi / y_qse(i);
|
||||
}
|
||||
J_qse(i, j) = y_qse(j) * (J_qse(i, j) - on_diag_correction) / y_qse(i); // Apply chain rule J'(i,j) = y_j * (J(i,j) - δ_ij(F_i/Y_i)) / Y_i
|
||||
}
|
||||
}
|
||||
|
||||
return 0; // Success
|
||||
}
|
||||
|
||||
@@ -1550,6 +1617,35 @@ namespace gridfire {
|
||||
return mean_timescale == other.mean_timescale;
|
||||
}
|
||||
|
||||
void MultiscalePartitioningEngineView::QSEGroup::removeSpecies(const Species &species) {
|
||||
if (algebraic_species.contains(species)) {
|
||||
algebraic_species.erase(species);
|
||||
}
|
||||
if (seed_species.contains(species)) {
|
||||
seed_species.erase(species);
|
||||
}
|
||||
}
|
||||
|
||||
void MultiscalePartitioningEngineView::QSEGroup::addSpeciesToAlgebraic(const Species &species) {
|
||||
if (seed_species.contains(species)) {
|
||||
const std::string msg = std::format("Cannot add species {} to algebraic set as it is already in the seed set.", species.name());
|
||||
throw std::invalid_argument(msg);
|
||||
}
|
||||
if (!algebraic_species.contains(species)) {
|
||||
algebraic_species.insert(species);
|
||||
}
|
||||
}
|
||||
|
||||
void MultiscalePartitioningEngineView::QSEGroup::addSpeciesToSeed(const Species &species) {
|
||||
if (algebraic_species.contains(species)) {
|
||||
const std::string msg = std::format("Cannot add species {} to seed set as it is already in the algebraic set.", species.name());
|
||||
throw std::invalid_argument(msg);
|
||||
}
|
||||
if (!seed_species.contains(species)) {
|
||||
seed_species.insert(species);
|
||||
}
|
||||
}
|
||||
|
||||
bool MultiscalePartitioningEngineView::QSEGroup::operator<(const QSEGroup &other) const {
|
||||
return mean_timescale < other.mean_timescale;
|
||||
}
|
||||
|
||||
@@ -18,7 +18,7 @@ namespace gridfire {
|
||||
|
||||
NetworkPrimingEngineView::NetworkPrimingEngineView(
|
||||
const std::string &primingSymbol,
|
||||
DynamicEngine &baseEngine
|
||||
GraphEngine &baseEngine
|
||||
) :
|
||||
DefinedEngineView(
|
||||
constructPrimingReactionSet(
|
||||
@@ -31,7 +31,7 @@ namespace gridfire {
|
||||
|
||||
NetworkPrimingEngineView::NetworkPrimingEngineView(
|
||||
const fourdst::atomic::Species &primingSpecies,
|
||||
DynamicEngine &baseEngine
|
||||
GraphEngine &baseEngine
|
||||
) :
|
||||
DefinedEngineView(
|
||||
constructPrimingReactionSet(
|
||||
@@ -46,7 +46,7 @@ namespace gridfire {
|
||||
|
||||
std::vector<std::string> NetworkPrimingEngineView::constructPrimingReactionSet(
|
||||
const fourdst::atomic::Species &primingSpecies,
|
||||
const DynamicEngine &baseEngine
|
||||
const GraphEngine &baseEngine
|
||||
) const {
|
||||
std::unordered_set<std::string> primeReactions;
|
||||
for (const auto &reaction : baseEngine.getNetworkReactions()) {
|
||||
|
||||
@@ -268,6 +268,7 @@ namespace gridfire::reaction {
|
||||
const double T9_p23 = std::pow(T9, r_p23);
|
||||
|
||||
|
||||
// ReSharper disable once CppUseStructuredBinding
|
||||
for (const auto& coeffs : m_rates) {
|
||||
const double exponent = coeffs.a0 +
|
||||
coeffs.a1 * T9_m1 +
|
||||
@@ -318,7 +319,7 @@ namespace gridfire::reaction {
|
||||
if (m_reactions.empty()) {
|
||||
return; // Case where the reactions will be added later.
|
||||
}
|
||||
m_reactionNameMap.reserve(reactions.size());
|
||||
m_reactionNameMap.reserve(m_reactions.size());
|
||||
size_t i = 0;
|
||||
for (const auto& reaction : m_reactions) {
|
||||
m_id += reaction->id();
|
||||
@@ -468,11 +469,11 @@ namespace gridfire::reaction {
|
||||
}
|
||||
|
||||
const Reaction& ReactionSet::operator[](const std::string_view& id) const {
|
||||
if (auto it = m_reactionNameMap.find(std::string(id)); it != m_reactionNameMap.end()) {
|
||||
if (const auto it = m_reactionNameMap.find(std::string(id)); it != m_reactionNameMap.end()) {
|
||||
return *m_reactions[it->second];
|
||||
}
|
||||
m_logger -> flush_log();
|
||||
throw std::out_of_range("Species " + std::string(id) + " does not exist in ReactionSet.");
|
||||
throw std::out_of_range("Reaction " + std::string(id) + " does not exist in ReactionSet.");
|
||||
}
|
||||
|
||||
bool ReactionSet::operator==(const ReactionSet& other) const {
|
||||
|
||||
@@ -144,6 +144,8 @@ namespace gridfire::solver {
|
||||
const auto relTol = m_config.get<double>("gridfire:solver:CVODESolverStrategy:relTol", 1.0e-8);
|
||||
|
||||
fourdst::composition::Composition equilibratedComposition = m_engine.update(netIn);
|
||||
std::cout << "Equilibrium d molar abundances: " << equilibratedComposition.getMolarAbundance(fourdst::atomic::H_2) << std::endl;
|
||||
std::cout << "Equilibrium d mass fraction: " << equilibratedComposition.getMassFraction(fourdst::atomic::H_2) << std::endl;
|
||||
std::cout << "EXITED AT EXPECTED TESTING POINT" << std::endl;
|
||||
exit(0);
|
||||
|
||||
@@ -251,7 +253,11 @@ namespace gridfire::solver {
|
||||
mass_fractions.push_back(y_data[i] * species.mass()); // Convert from molar abundance to mass fraction
|
||||
}
|
||||
temp_comp.setMassFraction(m_engine.getNetworkSpecies(), mass_fractions);
|
||||
temp_comp.finalize(true);
|
||||
bool didFinalize = temp_comp.finalize(true);
|
||||
if (!didFinalize) {
|
||||
LOG_ERROR(m_logger, "Failed to finalize composition during engine update. Check input mass fractions for validity.");
|
||||
throw std::runtime_error("Failed to finalize composition during engine update.");
|
||||
}
|
||||
|
||||
NetIn netInTemp = netIn;
|
||||
netInTemp.temperature = T9 * 1e9; // Convert back to Kelvin
|
||||
@@ -301,7 +307,11 @@ namespace gridfire::solver {
|
||||
|
||||
fourdst::composition::Composition outputComposition(speciesNames);
|
||||
outputComposition.setMassFraction(speciesNames, finalMassFractions);
|
||||
outputComposition.finalize(true);
|
||||
bool didFinalize = outputComposition.finalize(true);
|
||||
if (!didFinalize) {
|
||||
LOG_ERROR(m_logger, "Failed to finalize output composition after CVODE integration. Check output mass fractions for validity.");
|
||||
throw std::runtime_error("Failed to finalize output composition after CVODE integration.");
|
||||
}
|
||||
|
||||
NetOut netOut;
|
||||
netOut.composition = outputComposition;
|
||||
|
||||
Reference in New Issue
Block a user