fix(LogicalReaclibReaction): Properly class reverse reactions

Previously there was a bug where some reverse reactions were being
classed as forward reactions. This results in a failure of many
timesteps due to the reverse reactions very high molar flows
This commit is contained in:
2025-11-24 14:57:14 -05:00
parent ce8717b248
commit b335bf7100
5 changed files with 239 additions and 88 deletions

View File

@@ -643,15 +643,6 @@ namespace gridfire::engine {
std::set<fourdst::atomic::Species> seed_species; ///< Dynamic species in this group. std::set<fourdst::atomic::Species> seed_species; ///< Dynamic species in this group.
double mean_timescale; ///< Mean timescale of the group. double mean_timescale; ///< Mean timescale of the group.
// DEBUG METHODS.
// THESE SHOULD NOT BE USED IN PRODUCTION CODE.
[[deprecated("Use for debug only")]] void removeSpecies(const fourdst::atomic::Species& species);
[[deprecated("Use for debug only")]] void addSpeciesToAlgebraic(const fourdst::atomic::Species& species);
[[deprecated("Use for debug only")]] void addSpeciesToSeed(const fourdst::atomic::Species& species);
/** /**
* @brief Less-than operator for QSEGroup, used for sorting. * @brief Less-than operator for QSEGroup, used for sorting.
* @param other The other QSEGroup to compare to. * @param other The other QSEGroup to compare to.
@@ -739,6 +730,8 @@ namespace gridfire::engine {
const std::unordered_map<fourdst::atomic::Species, size_t>& qse_solve_species_index_map; const std::unordered_map<fourdst::atomic::Species, size_t>& qse_solve_species_index_map;
const std::vector<fourdst::atomic::Species>& qse_solve_species; const std::vector<fourdst::atomic::Species>& qse_solve_species;
const QSESolver& instance; const QSESolver& instance;
std::vector<double> row_scaling_factors;
const double initial_group_mass;
}; };
private: private:
@@ -943,6 +936,9 @@ namespace gridfire::engine {
* @brief Builds a connectivity graph from a species pool. * @brief Builds a connectivity graph from a species pool.
* *
* @param species_pool A vector of species indices representing a species pool. * @param species_pool A vector of species indices representing a species pool.
* @param comp
* @param T9
* @param rho
* @return An unordered map representing the adjacency list of the connectivity graph. * @return An unordered map representing the adjacency list of the connectivity graph.
* *
* @par Purpose * @par Purpose
@@ -955,7 +951,8 @@ namespace gridfire::engine {
* that reaction that are also in the pool. * that reaction that are also in the pool.
*/ */
std::unordered_map<fourdst::atomic::Species, std::vector<fourdst::atomic::Species>> buildConnectivityGraph( std::unordered_map<fourdst::atomic::Species, std::vector<fourdst::atomic::Species>> buildConnectivityGraph(
const std::vector<fourdst::atomic::Species>& species_pool const std::vector<fourdst::atomic::Species>& species_pool, const fourdst::composition::Composition &comp, double T9, double
rho
) const; ) const;
/** /**
@@ -989,6 +986,9 @@ namespace gridfire::engine {
* *
* @param timescale_pools A vector of vectors of species indices, where each inner vector * @param timescale_pools A vector of vectors of species indices, where each inner vector
* represents a timescale pool. * represents a timescale pool.
* @param comp
* @param T9
* @param rho
* @return A vector of vectors of species indices, where each inner vector represents a * @return A vector of vectors of species indices, where each inner vector represents a
* single connected component. * single connected component.
* *
@@ -1002,7 +1002,8 @@ namespace gridfire::engine {
* The resulting components from all pools are collected and returned. * The resulting components from all pools are collected and returned.
*/ */
std::vector<std::vector<fourdst::atomic::Species>> analyzeTimescalePoolConnectivity( std::vector<std::vector<fourdst::atomic::Species>> analyzeTimescalePoolConnectivity(
const std::vector<std::vector<fourdst::atomic::Species>> &timescale_pools const std::vector<std::vector<fourdst::atomic::Species>> &timescale_pools, const fourdst::composition::Composition &
comp, double T9, double rho
) const; ) const;
std::vector<QSEGroup> pruneValidatedGroups( std::vector<QSEGroup> pruneValidatedGroups(
@@ -1012,6 +1013,11 @@ namespace gridfire::engine {
double T9, double T9,
double rho double rho
) const; ) const;
static std::vector<QSEGroup> merge_coupled_groups(
const std::vector<QSEGroup> &groups,
const std::vector<reaction::ReactionSet> &groupReactions
);
}; };
} }

View File

@@ -46,20 +46,28 @@ namespace gridfire::utils {
}; };
static inline std::unordered_map<int, std::string> kinsol_ret_code_map { static inline std::unordered_map<int, std::string> kinsol_ret_code_map {
{0, "KIN_SUCCESS: The solver succeeded."}, {0, "KIN_SUCCESS: The solver succeeded."},
{1, "KIN_STEP_LT_STPTOL: The solver step size became less than the stopping tolerance."}, {1, "KIN_INITIAL_GUESS_OKAY: The guess, u=u0, satisfied the system F(u) = 0 within the tolerances specified"},
{2, "KIN_RES_REPTD_ERR: The residual function repeatedly failed recoverably."}, {2, "KIN_STEP_LT_STPTOL: KINSOL stopped based on scaled step length. This means that the current iterate may be an approximate solution of the given nonlinear system, but it is also quite possible that the algorithm is “stalled” (making insufficient progress) near an invalid solution, or that the scalar, scsteptol, is too large."},
{-1, "KIN_MEM_NULL: The KINSOL memory structure is NULL."}, {99, "KIN_WARNING: KINSOL succeeded but in an unusual way"},
{-2, "KIN_ILL_INPUT: An illegal input was detected."}, {-1, "KIN_MEM_NULL: The KINSOL memory pointer is NULL."},
{-2, "KIN_ILL_INPUT: An illegal value was specified for an input argument."},
{-3, "KIN_NO_MALLOC: The KINSOL memory structure has not been allocated."}, {-3, "KIN_NO_MALLOC: The KINSOL memory structure has not been allocated."},
{-4, "KIN_MEM_FAIL: Memory allocation failed."}, {-4, "KIN_MEM_FAIL: A memory allocation failed."},
{-5, "KIN_LINIT_FAIL: The linear solver's initialization function failed."}, {-5, "KIN_LINESEARCH_NONCONV: The line search algorithm was unable to find an iterate sufficiently distinct from the current iterate."},
{-6, "KIN_LSETUP_FAIL: The linear solver's setup function failed."}, {-6, "KIN_MAXITER_REACHED: The maximum number of iterations was reached before convergence."},
{-7, "KIN_LSOLVE_FAIL: The linear solver's solve function failed."}, {-7, "KIN_MXNEWT_5X_EXCEEDED: Five consecutive steps have been taken that satisfy a scaled step length test."},
{-8, "KIN_RESFUNC_FAIL: The residual function failed in an unrecoverable manner."}, {-8, "KIN_LINESEARCH_BCFAIL: The line search algorithm was unable to satisfy the beta-condition for nbcf fails iterations."},
{-9, "KIN_CONSTR_FAIL: The inequality constraint was violated and the solver was unable to recover."}, {-9, "KIN_LINSOLV_NO_RECOVERY: The linear solver's solve function failed recoverably, but the Jacobian data is already current."},
{-10, "KIN_NLS_INIT_FAIL: The nonlinear solver's initialization function failed."}, {-10, "KIN_LINIT_FAIL: The linear solver's init routine failed."},
{-11, "KIN_NLS_SETUP_FAIL: The nonlinear solver's setup function failed."}, {-11, "KIN_LSETUP_FAIL: The linear solver's setup function failed in an unrecoverable manner."},
{-12, "KIN_NLS_FAIL: The nonlinear solver's solve function failed."} {-12, "KIN_LSOLVE_FAIL: The linear solver's solve function failed in an unrecoverable manner."},
{-13, "KIN_SYSFUNC_FAIL: The system function failed in an unrecoverable manner."},
{-14, "KIN_FIRST_SYSFUNC_ERR: The system function failed at the first call."},
{-15, "KIN_REPTD_SYSFUNC_ERR: Unable to correct repeated recoverable system function errors."},
{-16, "KIN_VECTOROP_ERR: A vector operation failed."},
{-17, "KIN_CONTEXT_ERR: A context error occurred."},
{-18, "KIN_DAMPING_FN_ERR: The damping function failed."},
{-19, "KIN_DEPTH_FN_ERR: The depth function failed."}
}; };
inline const std::unordered_map<int, std::string>& sundials_retcode_map(const SUNDIALS_RET_CODE_TYPES type) { inline const std::unordered_map<int, std::string>& sundials_retcode_map(const SUNDIALS_RET_CODE_TYPES type) {
@@ -96,4 +104,27 @@ namespace gridfire::utils {
return vec; return vec;
} }
inline void check_sundials_flag(const int flag, const std::string& func_name, const SUNDIALS_RET_CODE_TYPES type) {
if (flag < 0) {
const auto& ret_code_map = sundials_retcode_map(type);
if (!ret_code_map.contains(flag)) {
switch (type) {
case (SUNDIALS_RET_CODE_TYPES::CVODE):
throw exceptions::CVODESolverFailureError("CVODE error in " + func_name + ": Unknown error code: " + std::to_string(flag));
case (SUNDIALS_RET_CODE_TYPES::KINSOL):
throw exceptions::KINSolSolverFailureError("KINSol error in " + func_name + ": Unknown error code: " + std::to_string(flag));
default:
throw exceptions::CVODESolverFailureError("SUNDIALS error in " + func_name + ": Unknown error code: " + std::to_string(flag));
}
}
switch (type) {
case (SUNDIALS_RET_CODE_TYPES::CVODE):
throw exceptions::CVODESolverFailureError("CVODE error in " + func_name + ": " + ret_code_map.at(flag));
case (SUNDIALS_RET_CODE_TYPES::KINSOL):
throw exceptions::KINSolSolverFailureError("KINSol error in " + func_name + ": " + ret_code_map.at(flag));
default:
throw exceptions::CVODESolverFailureError("SUNDIALS error in " + func_name + ": " + ret_code_map.at(flag));
}
}
}
} }

View File

@@ -14,6 +14,7 @@
#include <queue> #include <queue>
#include <algorithm> #include <algorithm>
#include <numeric>
#include "fourdst/atomic/species.h" #include "fourdst/atomic/species.h"
#include "quill/LogMacros.h" #include "quill/LogMacros.h"
@@ -161,6 +162,28 @@ namespace {
// For everything else, use the default SUNDIALS logger (or your own) // For everything else, use the default SUNDIALS logger (or your own)
SUNLogErrHandlerFn(line, func, file, msg, err_code, err_user_data, sunctx); SUNLogErrHandlerFn(line, func, file, msg, err_code, err_user_data, sunctx);
} }
struct DisjointSet {
std::vector<size_t> parent;
explicit DisjointSet(const size_t n) {
parent.resize(n);
std::iota(parent.begin(), parent.end(), 0);
}
size_t find(const size_t i) {
if (parent.at(i) == i) return i;
return parent.at(i) = find(parent.at(i)); // Path compression
}
void unite(const size_t i, const size_t j) {
const size_t root_i = find(i);
const size_t root_j = find(j);
if (root_i != root_j) {
parent.at(root_i) = root_j;
}
}
};
} }
namespace gridfire::engine { namespace gridfire::engine {
@@ -405,7 +428,8 @@ namespace gridfire::engine {
} }
std::vector<std::vector<Species>> MultiscalePartitioningEngineView::analyzeTimescalePoolConnectivity( std::vector<std::vector<Species>> MultiscalePartitioningEngineView::analyzeTimescalePoolConnectivity(
const std::vector<std::vector<Species>> &timescale_pools const std::vector<std::vector<Species>> &timescale_pools, const fourdst::composition::Composition &comp, double T9, double
rho
) const { ) const {
std::vector<std::vector<Species>> final_connected_pools; std::vector<std::vector<Species>> final_connected_pools;
@@ -415,7 +439,7 @@ namespace gridfire::engine {
} }
// For each timescale pool, we need to analyze connectivity. // For each timescale pool, we need to analyze connectivity.
auto connectivity_graph = buildConnectivityGraph(pool); auto connectivity_graph = buildConnectivityGraph(pool, comp, T9, rho);
auto components = findConnectedComponentsBFS(connectivity_graph, pool); auto components = findConnectedComponentsBFS(connectivity_graph, pool);
final_connected_pools.insert(final_connected_pools.end(), components.begin(), components.end()); final_connected_pools.insert(final_connected_pools.end(), components.begin(), components.end());
} }
@@ -608,6 +632,60 @@ namespace gridfire::engine {
return newGroups; return newGroups;
} }
std::vector<MultiscalePartitioningEngineView::QSEGroup> MultiscalePartitioningEngineView::merge_coupled_groups(
const std::vector<QSEGroup> &groups,
const std::vector<reaction::ReactionSet> &groupReactions
) {
if (groups.empty()) { return {}; }
DisjointSet dsu(groups.size());
std::unordered_map<std::string_view, size_t> reaction_owner;
for (size_t group_idx = 0; group_idx < groups.size(); ++group_idx) {
for (const auto& reaction : groupReactions[group_idx]) {
const std::string_view rxn_id = reaction->id();
if (reaction_owner.contains(rxn_id)) {
dsu.unite(group_idx, reaction_owner[rxn_id]);
} else {
reaction_owner[rxn_id] = group_idx;
}
}
}
std::unordered_map<size_t, std::pair<QSEGroup, size_t>> merged_groups_map;
for (size_t i = 0; i < groups.size(); ++i) {
size_t root = dsu.find(i);
if (!merged_groups_map.contains(root)) {
merged_groups_map[root].first = QSEGroup();
merged_groups_map.at(root).first.mean_timescale = 0.0;
merged_groups_map.at(root).first.is_in_equilibrium = true;
}
QSEGroup& target_group = merged_groups_map.at(root).first;
const QSEGroup& source_group = groups[i];
target_group.seed_species.insert(source_group.seed_species.begin(), source_group.seed_species.end());
target_group.algebraic_species.insert(source_group.algebraic_species.begin(), source_group.algebraic_species.end());
target_group.mean_timescale += source_group.mean_timescale;
merged_groups_map.at(root).second += 1;
}
std::vector<QSEGroup> final_merged_groups;
for (auto &groupInfo: merged_groups_map | std::views::values) {
QSEGroup& group = groupInfo.first;
for (const auto& alg_species : group.algebraic_species) {
group.seed_species.erase(alg_species); // A species cannot be both seed and algebraic. Algebraic takes precedence.
}
group.mean_timescale /= static_cast<double>(groupInfo.second);
final_merged_groups.push_back(std::move(group));
}
return final_merged_groups;
}
fourdst::composition::Composition MultiscalePartitioningEngineView::partitionNetwork( fourdst::composition::Composition MultiscalePartitioningEngineView::partitionNetwork(
const NetIn &netIn const NetIn &netIn
) { ) {
@@ -650,7 +728,7 @@ namespace gridfire::engine {
} }
LOG_TRACE_L1(m_logger, "Preforming connectivity analysis on timescale pools..."); LOG_TRACE_L1(m_logger, "Preforming connectivity analysis on timescale pools...");
const std::vector<std::vector<Species>> connected_pools = analyzeTimescalePoolConnectivity(candidate_pools); const std::vector<std::vector<Species>> connected_pools = analyzeTimescalePoolConnectivity(candidate_pools, comp, T9, rho);
LOG_TRACE_L1(m_logger, "Found {} connected pools (compared to {} timescale pools) for QSE analysis.", connected_pools.size(), timescale_pools.size()); LOG_TRACE_L1(m_logger, "Found {} connected pools (compared to {} timescale pools) for QSE analysis.", connected_pools.size(), timescale_pools.size());
@@ -663,14 +741,19 @@ namespace gridfire::engine {
const auto [validated_groups, invalidate_groups, validated_group_reactions] = validateGroupsWithFluxAnalysis(candidate_groups, comp, T9, rho); const auto [validated_groups, invalidate_groups, validated_group_reactions] = validateGroupsWithFluxAnalysis(candidate_groups, comp, T9, rho);
LOG_TRACE_L1(m_logger, "Validated {} group(s) QSE groups. {}", validated_groups.size(), utils::iterable_to_delimited_string(validated_groups)); LOG_TRACE_L1(m_logger, "Validated {} group(s) QSE groups. {}", validated_groups.size(), utils::iterable_to_delimited_string(validated_groups));
LOG_TRACE_L1(m_logger, "Pruning groups based on log abundance-normalized flux analysis..."); LOG_TRACE_L1(m_logger, "Pruning groups based on log abundance-normalized flux analysis...");
const std::vector<QSEGroup> prunedGroups = pruneValidatedGroups(validated_groups, validated_group_reactions, comp, T9, rho); const std::vector<QSEGroup> prunedGroups = pruneValidatedGroups(validated_groups, validated_group_reactions, comp, T9, rho);
LOG_TRACE_L1(m_logger, "After Pruning remaining groups are: {}", utils::iterable_to_delimited_string(prunedGroups)); LOG_TRACE_L1(m_logger, "After Pruning remaining groups are: {}", utils::iterable_to_delimited_string(prunedGroups));
LOG_TRACE_L1(m_logger, "Re-validating pruned groups with flux analysis..."); LOG_TRACE_L1(m_logger, "Re-validating pruned groups with flux analysis...");
auto [pruned_validated_groups, _, __] = validateGroupsWithFluxAnalysis(prunedGroups, comp, T9, rho); auto [pruned_validated_groups, _, pruned_validated_reactions] = validateGroupsWithFluxAnalysis(prunedGroups, comp, T9, rho);
LOG_TRACE_L1(m_logger, "After re-validation, {} QSE groups remain. ({})",pruned_validated_groups.size(), utils::iterable_to_delimited_string(pruned_validated_groups)); LOG_TRACE_L1(m_logger, "After re-validation, {} QSE groups remain. ({})",pruned_validated_groups.size(), utils::iterable_to_delimited_string(pruned_validated_groups));
LOG_TRACE_L1(m_logger, "Merging coupled QSE groups...");
const std::vector<QSEGroup> merged_groups = merge_coupled_groups(pruned_validated_groups, pruned_validated_reactions);
LOG_TRACE_L1(m_logger, "After merging, {} QSE groups remain. ({})", merged_groups.size(), utils::iterable_to_delimited_string(merged_groups));
m_qse_groups = pruned_validated_groups; m_qse_groups = pruned_validated_groups;
LOG_TRACE_L1(m_logger, "Pushing all identified algebraic species into algebraic set..."); LOG_TRACE_L1(m_logger, "Pushing all identified algebraic species into algebraic set...");
@@ -1509,7 +1592,10 @@ namespace gridfire::engine {
} }
std::unordered_map<Species, std::vector<Species>> MultiscalePartitioningEngineView::buildConnectivityGraph( std::unordered_map<Species, std::vector<Species>> MultiscalePartitioningEngineView::buildConnectivityGraph(
const std::vector<Species> &species_pool const std::vector<Species> &species_pool,
const fourdst::composition::Composition &comp,
double T9,
double rho
) const { ) const {
std::unordered_map<Species, std::vector<Species>> connectivity_graph; std::unordered_map<Species, std::vector<Species>> connectivity_graph;
const std::set<Species> pool_set(species_pool.begin(), species_pool.end()); const std::set<Species> pool_set(species_pool.begin(), species_pool.end());
@@ -1714,26 +1800,28 @@ namespace gridfire::engine {
N_VConst(1.0, m_constraints); N_VConst(1.0, m_constraints);
m_kinsol_mem = KINCreate(m_sun_ctx); m_kinsol_mem = KINCreate(m_sun_ctx);
utils::check_cvode_flag(m_kinsol_mem ? 0 : -1, "KINCreate"); utils::check_sundials_flag(m_kinsol_mem ? 0 : -1, "KINCreate", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
utils::check_cvode_flag(KINInit(m_kinsol_mem, sys_func, m_func_tmpl), "KINInit"); utils::check_sundials_flag(KINInit(m_kinsol_mem, sys_func, m_func_tmpl), "KINInit", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
utils::check_cvode_flag(KINSetConstraints(m_kinsol_mem, m_constraints), "KINSetConstraints"); utils::check_sundials_flag(KINSetConstraints(m_kinsol_mem, m_constraints), "KINSetConstraints", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
m_J = SUNDenseMatrix(static_cast<sunindextype>(m_N), static_cast<sunindextype>(m_N), m_sun_ctx); m_J = SUNDenseMatrix(static_cast<sunindextype>(m_N), static_cast<sunindextype>(m_N), m_sun_ctx);
utils::check_cvode_flag(m_J ? 0 : -1, "SUNDenseMatrix"); utils::check_sundials_flag(m_J ? 0 : -1, "SUNDenseMatrix", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
m_LS = SUNLinSol_Dense(m_Y, m_J, m_sun_ctx); m_LS = SUNLinSol_Dense(m_Y, m_J, m_sun_ctx);
utils::check_cvode_flag(m_LS ? 0 : -1, "SUNLinSol_Dense"); utils::check_sundials_flag(m_LS ? 0 : -1, "SUNLinSol_Dense", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
utils::check_cvode_flag(KINSetLinearSolver(m_kinsol_mem, m_LS, m_J), "KINSetLinearSolver"); utils::check_sundials_flag(KINSetLinearSolver(m_kinsol_mem, m_LS, m_J), "KINSetLinearSolver", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
utils::check_cvode_flag(KINSetJacFn(m_kinsol_mem, sys_jac), "KINSetJacFn"); utils::check_sundials_flag(KINSetJacFn(m_kinsol_mem, sys_jac), "KINSetJacFn", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
utils::check_cvode_flag(KINSetMaxSetupCalls(m_kinsol_mem, 20), "KINSetMaxSetupCalls"); utils::check_sundials_flag(KINSetMaxSetupCalls(m_kinsol_mem, 20), "KINSetMaxSetupCalls", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
utils::check_cvode_flag(KINSetFuncNormTol(m_kinsol_mem, 1e-6), "KINSetFuncNormTol"); utils::check_sundials_flag(KINSetFuncNormTol(m_kinsol_mem, 1e-6), "KINSetFuncNormTol", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
utils::check_cvode_flag(KINSetNumMaxIters(m_kinsol_mem, 200), "KINSetNumMaxIters"); utils::check_sundials_flag(KINSetNumMaxIters(m_kinsol_mem, 200), "KINSetNumMaxIters", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
utils::check_sundials_flag(KINSetScaledStepTol(m_kinsol_mem, 1e-10), "KINSetScaledStepTol", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
// We want to effectively disable this since enormous changes in order of magnitude are realistic for this problem. // We want to effectively disable this since enormous changes in order of magnitude are realistic for this problem.
utils::check_cvode_flag(KINSetMaxNewtonStep(m_kinsol_mem, std::numeric_limits<double>::infinity()), "KINSetMaxNewtonStep"); utils::check_sundials_flag(KINSetMaxNewtonStep(m_kinsol_mem, std::numeric_limits<double>::infinity()), "KINSetMaxNewtonStep", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
LOG_TRACE_L1(getLogger(), "Finished initializing QSE Solver."); LOG_TRACE_L1(getLogger(), "Finished initializing QSE Solver.");
} }
@@ -1789,7 +1877,7 @@ namespace gridfire::engine {
*this *this
}; };
utils::check_cvode_flag(KINSetUserData(m_kinsol_mem, &data), "KINSetUserData"); utils::check_sundials_flag(KINSetUserData(m_kinsol_mem, &data), "KINSetUserData", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
sunrealtype* y_data = N_VGetArrayPointer(m_Y); sunrealtype* y_data = N_VGetArrayPointer(m_Y);
sunrealtype* scale_data = N_VGetArrayPointer(m_scale); sunrealtype* scale_data = N_VGetArrayPointer(m_scale);
@@ -1819,9 +1907,9 @@ namespace gridfire::engine {
if (m_solves > 0 && m_has_jacobian) { if (m_solves > 0 && m_has_jacobian) {
// After the initial solution we want to allow kinsol to reuse its state // After the initial solution we want to allow kinsol to reuse its state
utils::check_cvode_flag(KINSetNoInitSetup(m_kinsol_mem, SUNTRUE), "KINSetNoInitSetup"); utils::check_sundials_flag(KINSetNoInitSetup(m_kinsol_mem, SUNTRUE), "KINSetNoInitSetup", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
} else { } else {
utils::check_cvode_flag(KINSetNoInitSetup(m_kinsol_mem, SUNFALSE), "KINSetNoInitSetup"); utils::check_sundials_flag(KINSetNoInitSetup(m_kinsol_mem, SUNFALSE), "KINSetNoInitSetup", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
} }
LOG_TRACE_L2( LOG_TRACE_L2(
@@ -1874,10 +1962,21 @@ namespace gridfire::engine {
); );
const int flag = KINSol(m_kinsol_mem, m_Y, KIN_LINESEARCH, m_scale, m_f_scale); const int flag = KINSol(m_kinsol_mem, m_Y, KIN_LINESEARCH, m_scale, m_f_scale);
if (flag < 0) { if (flag < 0) {
LOG_WARNING(getLogger(), "KINSol failed to converge while solving QSE abundances with flag {}.", utils::cvode_ret_code_map.at(flag)); double fnorm = 0.0;
throw exceptions::InvalidQSESolutionError("KINSol failed to converge while solving QSE abundances. Check the log file for more details regarding the specific failure mode."); constexpr double ACCEPTABLE_FTOL = 1e-4; // Due to row scaling this is a relative threshold. So this is saying that the imbalance is 0.01% of the abundance scale.
utils::check_sundials_flag(KINGetFuncNorm(m_kinsol_mem, &fnorm), "KINGetFuncNorm", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
if (fnorm < ACCEPTABLE_FTOL) {
LOG_INFO(getLogger(), "KINSol failed to converge within the maximum number of iterations, but achieved acceptable accuracy with function norm {} < {}. Proceeding with solution.",
fnorm, ACCEPTABLE_FTOL);
} else {
LOG_WARNING(getLogger(), "KINSol failed to converge while solving QSE abundances with flag {}. Error {}", utils::kinsol_ret_code_map.at(flag), fnorm);
throw exceptions::InvalidQSESolutionError("KINSol failed to converge while solving QSE abundances. " + utils::kinsol_ret_code_map.at(flag));
} }
}
for (size_t i = 0; i < m_N; ++i) { for (size_t i = 0; i < m_N; ++i) {
const auto& species = m_species[i]; const auto& species = m_species[i];
@@ -1895,11 +1994,11 @@ namespace gridfire::engine {
void MultiscalePartitioningEngineView::QSESolver::log_diagnostics(const QSEGroup &group, const fourdst::composition::Composition &comp) const { void MultiscalePartitioningEngineView::QSESolver::log_diagnostics(const QSEGroup &group, const fourdst::composition::Composition &comp) const {
long int nni, nfe, nje; long int nni, nfe, nje;
utils::check_cvode_flag(KINGetNumNonlinSolvIters(m_kinsol_mem, &nni), "KINGetNumNonlinSolvIters"); utils::check_sundials_flag(KINGetNumNonlinSolvIters(m_kinsol_mem, &nni), "KINGetNumNonlinSolvIters", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
utils::check_cvode_flag(KINGetNumFuncEvals(m_kinsol_mem, &nfe), "KINGetNumFuncEvals"); utils::check_sundials_flag(KINGetNumFuncEvals(m_kinsol_mem, &nfe), "KINGetNumFuncEvals", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
utils::check_cvode_flag(KINGetNumJacEvals(m_kinsol_mem, &nje), "KINGetNumJacEvals"); utils::check_sundials_flag(KINGetNumJacEvals(m_kinsol_mem, &nje), "KINGetNumJacEvals", utils::SUNDIALS_RET_CODE_TYPES::KINSOL);
LOG_TRACE_L1(getLogger(), LOG_TRACE_L1(getLogger(),
"QSE Stats | Iters: {} | RHS Evals: {} | Jac Evals: {} | Ratio (J/I): {:.2f} | Algebraic Species: {}", "QSE Stats | Iters: {} | RHS Evals: {} | Jac Evals: {} | Ratio (J/I): {:.2f} | Algebraic Species: {}",
@@ -2008,7 +2107,8 @@ namespace gridfire::engine {
#endif #endif
for (const auto& [species, index] : map) { for (const auto& [species, index] : map) {
f_data[index] = dydt.at(species); const double val = dydt.at(species);
f_data[index] = val;
} }
return 0; // Success return 0; // Success
@@ -2018,12 +2118,12 @@ namespace gridfire::engine {
int MultiscalePartitioningEngineView::QSESolver::sys_jac( int MultiscalePartitioningEngineView::QSESolver::sys_jac(
const N_Vector y, const N_Vector y,
N_Vector fy, N_Vector fy,
SUNMatrix J, const SUNMatrix J,
void *user_data, void *user_data,
N_Vector tmp1, N_Vector tmp1,
N_Vector tmp2 N_Vector tmp2
) { ) {
const auto* data = static_cast<UserData*>(user_data); auto* data = static_cast<UserData*>(user_data);
const sunrealtype* y_data = N_VGetArrayPointer(y); const sunrealtype* y_data = N_VGetArrayPointer(y);
const auto& map = data->qse_solve_species_index_map; const auto& map = data->qse_solve_species_index_map;
@@ -2041,12 +2141,48 @@ namespace gridfire::engine {
sunrealtype* J_data = SUNDenseMatrix_Data(J); sunrealtype* J_data = SUNDenseMatrix_Data(J);
const sunindextype N = SUNDenseMatrix_Columns(J); const sunindextype N = SUNDenseMatrix_Columns(J);
if (data->row_scaling_factors.size() != static_cast<size_t>(N)) {
data->row_scaling_factors.resize(N, 0.0);
}
for (const auto& [row_species, row_idx]: map) {
double max_value = std::numeric_limits<double>::lowest();
for (const auto &col_species: map | std::views::keys) {
const double val = std::abs(jac(row_species, col_species));
if (val > max_value) {
max_value = val;
}
}
for (const auto& [col_species, col_idx] : map) { for (const auto& [col_species, col_idx] : map) {
for (const auto& [row_species, row_idx] : map) {
J_data[col_idx * N + row_idx] = jac(row_species, col_species); J_data[col_idx * N + row_idx] = jac(row_species, col_species);
} }
} }
LOG_TRACE_L2( getLogger(), "After filling, Jacobian in sys_jac is: {}",
[&]() -> std::string {
std::ostringstream oss;
const sunrealtype* J_log_data = SUNDenseMatrix_Data(J);
const sunindextype N_log = SUNDenseMatrix_Columns(J);
oss << "[";
for (size_t i = 0; i < data->qse_solve_species.size(); ++i) {
oss << "[";
for (size_t j = 0; j < data->qse_solve_species.size(); ++j) {
oss << J_log_data[i * N_log + j];
if (j < data->qse_solve_species.size() - 1) {
oss << ", ";
}
}
oss << "]";
if (i < data->qse_solve_species.size() - 1) {
oss << ", ";
}
}
oss << "]";
return oss.str();
}()
);
data->instance.m_has_jacobian = true; data->instance.m_has_jacobian = true;
return 0; return 0;
@@ -2062,35 +2198,6 @@ namespace gridfire::engine {
return mean_timescale == other.mean_timescale; 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 {:8} 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 {:8} 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 { bool MultiscalePartitioningEngineView::QSEGroup::operator<(const QSEGroup &other) const {
return mean_timescale < other.mean_timescale; return mean_timescale < other.mean_timescale;
} }

View File

@@ -573,7 +573,7 @@ namespace gridfire::reaction {
} }
// Now, process the grouped REACLIB reactions // Now, process the grouped REACLIB reactions
for (const auto &reactionsGroup: groupedReaclibReactions | std::views::values) { for (const auto &[key, reactionsGroup]: groupedReaclibReactions) {
if (reactionsGroup.empty()) { if (reactionsGroup.empty()) {
continue; continue;
} }
@@ -581,7 +581,14 @@ namespace gridfire::reaction {
finalReactionSet.add_reaction(reactionsGroup.front()); finalReactionSet.add_reaction(reactionsGroup.front());
} }
else { else {
const auto logicalReaction = std::make_unique<LogicalReaclibReaction>(reactionsGroup); // Check that is_reverse is consistent across the group
assert(std::ranges::all_of(
reactionsGroup,
[&reactionsGroup](const ReaclibReaction& r) {
return r.is_reverse() == reactionsGroup.front().is_reverse();
}
) && "Inconsistent is_reverse values in grouped REACLIB reactions.");
const auto logicalReaction = std::make_unique<LogicalReaclibReaction>(reactionsGroup, reactionsGroup.front().is_reverse());
finalReactionSet.add_reaction(logicalReaction->clone()); finalReactionSet.add_reaction(logicalReaction->clone());
} }
} }

View File

@@ -566,11 +566,11 @@ namespace gridfire::solver {
LOG_TRACE_L2(instance->m_logger, "CVODE RHS wrapper completed successfully at time {}", t); LOG_TRACE_L2(instance->m_logger, "CVODE RHS wrapper completed successfully at time {}", t);
return 0; return 0;
} catch (const exceptions::EngineError& e) { } catch (const exceptions::EngineError& e) {
LOG_ERROR(instance->m_logger, "EngineError caught in CVODE RHS wrapper at time {}: {}", t, e.what()); LOG_ERROR(instance->m_logger, "EngineError caught in CVODE RHS wrapper at time {}: {}. Will attempt to recover...", t, e.what());
data->captured_exception = std::make_unique<exceptions::EngineError>(e);
return 1; // 1 Indicates a recoverable error, CVODE will retry the step return 1; // 1 Indicates a recoverable error, CVODE will retry the step
} catch (...) { } catch (const std::exception& e) {
LOG_CRITICAL(instance->m_logger, "Unrecoverable and Unknown exception caught in CVODE RHS wrapper at time {}", t); LOG_CRITICAL(instance->m_logger, "Unrecoverable and Unknown exception caught in CVODE RHS wrapper at time {} ({})", t, e.what());
// data->captured_exception = std::make_unique<exceptions::GridFireError>(e.what());
return -1; // Some unrecoverable error return -1; // Some unrecoverable error
} }
} }