GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
engine_multiscale.cpp
Go to the documentation of this file.
3
4#include <stdexcept>
5#include <vector>
6#include <unordered_map>
7#include <unordered_set>
8
9#include <queue>
10
11#include <ranges>
12#include <algorithm>
13
14#include "quill/LogMacros.h"
15#include "quill/Logger.h"
16
17namespace {
18 using namespace fourdst::atomic;
19 std::vector<double> packCompositionToVector(const fourdst::composition::Composition& composition, const gridfire::GraphEngine& engine) {
20 std::vector<double> Y(engine.getNetworkSpecies().size(), 0.0);
21 const auto& allSpecies = engine.getNetworkSpecies();
22 for (size_t i = 0; i < allSpecies.size(); ++i) {
23 const auto& species = allSpecies[i];
24 if (!composition.contains(species)) {
25 Y[i] = 0.0; // Species not in the composition, set to zero
26 } else {
27 Y[i] = composition.getMolarAbundance(species);
28 }
29 }
30 return Y;
31 }
32
33 template <class T>
34 void hash_combine(std::size_t& seed, const T& v) {
35 std::hash<T> hashed;
36 seed ^= hashed(v) + 0x9e3779b9 + (seed << 6) + (seed >> 2);
37 }
38
39 std::vector<std::vector<size_t>> findConnectedComponentsBFS(
40 const std::unordered_map<size_t, std::vector<size_t>>& graph,
41 const std::vector<size_t>& nodes
42 ) {
43 std::vector<std::vector<size_t>> components;
44 std::unordered_set<size_t> visited;
45
46 for (const size_t& start_node : nodes) {
47 if (visited.find(start_node) == visited.end()) {
48 std::vector<size_t> current_component;
49 std::queue<size_t> q;
50
51 q.push(start_node);
52 visited.insert(start_node);
53
54 while (!q.empty()) {
55 size_t u = q.front();
56 q.pop();
57 current_component.push_back(u);
58
59 if (graph.count(u)) {
60 for (const auto& v : graph.at(u)) {
61 if (visited.find(v) == visited.end()) {
62 visited.insert(v);
63 q.push(v);
64 }
65 }
66 }
67 }
68 components.push_back(current_component);
69 }
70 }
71 return components;
72 }
73
74 struct SpeciesSetIntersection {
75 const Species species;
76 std::size_t count;
77 };
78
79 std::expected<SpeciesSetIntersection, std::string> get_intersection_info (
80 const std::unordered_set<Species>& setA,
81 const std::unordered_set<Species>& setB
82 ) {
83 // Iterate over the smaller of the two
84 auto* outerSet = &setA;
85 auto* innerSet = &setB;
86 if (setA.size() > setB.size()) {
87 outerSet = &setB;
88 innerSet = &setA;
89 }
90
91 std::size_t matchCount = 0;
92 const Species* firstMatch = nullptr;
93
94 for (const Species& sp : *outerSet) {
95 if (innerSet->contains(sp)) {
96 if (matchCount == 0) {
97 firstMatch = &sp;
98 }
99 ++matchCount;
100 if (matchCount > 1) {
101 break;
102 }
103 }
104 }
105 if (!firstMatch) {
106 // No matches found
107 return std::unexpected{"Intersection is empty"};
108 }
109 if (matchCount == 0) {
110 // No matches found
111 return std::unexpected{"No intersection found"};
112 }
113
114 // Return the first match and the count of matches
115 return SpeciesSetIntersection{*firstMatch, matchCount};
116
117 }
118
119 bool has_distinct_reactant_and_product_species (
120 const std::unordered_set<Species>& poolSpecies,
121 const std::unordered_set<Species>& reactants,
122 const std::unordered_set<Species>& products
123 ) {
124 const auto reactant_result = get_intersection_info(poolSpecies, reactants);
125 if (!reactant_result) {
126 return false; // No reactants found
127 }
128 const auto [reactantSample, reactantCount] = reactant_result.value();
129
130 const auto product_result = get_intersection_info(poolSpecies, products);
131 if (!product_result) {
132 return false; // No products found
133 }
134
135 const auto [productSample, productCount] = product_result.value();
136
137 // If either side has ≥2 distinct matches, we can always pick
138 // one from each that differ.
139 if (reactantCount > 1 || productCount > 1) {
140 return true;
141 }
142
143 // Exactly one match on each side → they must differ
144 return reactantSample != productSample;
145 }
146
147}
148
149namespace gridfire {
151 using fourdst::atomic::Species;
152
156
157 const std::vector<Species> & MultiscalePartitioningEngineView::getNetworkSpecies() const {
158 return m_baseEngine.getNetworkSpecies();
159 }
160
162 const std::vector<double> &Y_full,
163 const double T9,
164 const double rho
165 ) const {
166 if (Y_full.size() != getNetworkSpecies().size()) {
167 LOG_ERROR(
168 m_logger,
169 "Y_full size ({}) does not match the number of species in the network ({})",
170 Y_full.size(),
171 getNetworkSpecies().size()
172 );
173 throw std::runtime_error("Y_full size does not match the number of species in the network. See logs for more details...");
174 }
175
176 // 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.
177 const QSECacheKey key(T9, rho, Y_full);
178 if (! m_qse_abundance_cache.contains(key)) {
180 LOG_ERROR(
181 m_logger,
182 "QSE abundance cache miss for T9 = {}, rho = {} (misses: {}, hits: {}). calculateRHSAndEnergy does not receive sufficient context to partition and stabilize the network. Throwing an error which should be caught by the caller and trigger a re-partition stage.",
183 T9,
184 rho,
185 m_cacheStats.misses(),
186 m_cacheStats.hits()
187 );
189 }
191 const auto result = m_baseEngine.calculateRHSAndEnergy(Y_full, T9, rho);
192 if (!result) {
193 return std::unexpected{result.error()};
194 }
195 auto deriv = result.value();
196
197 for (size_t i = 0; i < m_algebraic_species_indices.size(); ++i) {
198 const size_t species_index = m_algebraic_species_indices[i];
199 deriv.dydt[species_index] = 0.0; // Fix the algebraic species to the equilibrium abundances we calculate.
200 }
201 return deriv;
202 }
203
205 const std::vector<double> &Y_full,
206 const double T9,
207 const double rho
208 ) const {
209 const QSECacheKey key(T9, rho, Y_full);
210 if (!m_qse_abundance_cache.contains(key)) {
212 LOG_ERROR(
213 m_logger,
214 "QSE abundance cache miss for T9 = {}, rho = {} (misses: {}, hits: {}). generateJacobianMatrix does not receive sufficient context to partition and stabilize the network. Throwing an error which should be caught by the caller and trigger a re-partition stage.",
215 T9,
216 rho,
217 m_cacheStats.misses(),
218 m_cacheStats.hits()
219 );
220 throw exceptions::StaleEngineError("QSE Cache Miss while lacking context for partitioning. This should be caught by the caller and trigger a re-partition stage.");
221 }
223
224 // TODO: Add sparsity pattern to this to prevent base engine from doing unnecessary work.
225 m_baseEngine.generateJacobianMatrix(Y_full, T9, rho);
226 }
227
229 const int i_full,
230 const int j_full
231 ) const {
232 // // Check if the species we are differentiating with respect to is algebraic or dynamic. If it is algebraic we can reduce the work significantly...
233 // if (std::ranges::contains(m_algebraic_species_indices, j_full)) {
234 // const auto& species = m_baseEngine.getNetworkSpecies()[j_full];
235 // // If j is algebraic, we can return 0.0 since the Jacobian entry for algebraic species is always zero.
236 // return 0.0;
237 // }
238 // Otherwise we need to query the full jacobian
239 return m_baseEngine.getJacobianMatrixEntry(i_full, j_full);
240 }
241
243 m_baseEngine.generateStoichiometryMatrix();
244 }
245
247 const int speciesIndex,
248 const int reactionIndex
249 ) const {
250 return m_baseEngine.getStoichiometryMatrixEntry(speciesIndex, reactionIndex);
251 }
252
255 const std::vector<double> &Y_full,
256 const double T9,
257 const double rho
258 ) const {
259 const auto key = QSECacheKey(T9, rho, Y_full);
260 if (!m_qse_abundance_cache.contains(key)) {
262 LOG_ERROR(
263 m_logger,
264 "QSE abundance cache miss for T9 = {}, rho = {} (misses: {}, hits: {}). calculateMolarReactionFlow does not receive sufficient context to partition and stabilize the network. Throwing an error which should be caught by the caller and trigger a re-partition stage.",
265 T9,
266 rho,
267 m_cacheStats.misses(),
268 m_cacheStats.hits()
269 );
270 throw exceptions::StaleEngineError("QSE Cache Miss while lacking context for partitioning. This should be caught by the caller and trigger a re-partition stage.");
271 }
273 std::vector<double> Y_algebraic = m_qse_abundance_cache.at(key);
274
275 assert(Y_algebraic.size() == m_algebraic_species_indices.size());
276
277 // Fix the algebraic species to the equilibrium abundances we calculate.
278 std::vector<double> Y_mutable = Y_full;
279 for (const auto& [index, Yi] : std::views::zip(m_algebraic_species_indices, Y_algebraic)) {
280 Y_mutable[index] = Yi;
281
282 }
283 return m_baseEngine.calculateMolarReactionFlow(reaction, Y_mutable, T9, rho);
284 }
285
289
291 LOG_CRITICAL(m_logger, "setNetworkReactions is not supported in MultiscalePartitioningEngineView. Did you mean to call this on the base engine?");
292 throw exceptions::UnableToSetNetworkReactionsError("setNetworkReactions is not supported in MultiscalePartitioningEngineView. Did you mean to call this on the base engine?");
293 }
294
296 const std::vector<double> &Y,
297 const double T9,
298 const double rho
299 ) const {
300 const auto key = QSECacheKey(T9, rho, Y);
301 if (!m_qse_abundance_cache.contains(key)) {
303 LOG_ERROR(
304 m_logger,
305 "QSE abundance cache miss for T9 = {}, rho = {} (misses: {}, hits: {}). getSpeciesTimescales does not receive sufficient context to partition and stabilize the network. Throwing an error which should be caught by the caller and trigger a re-partition stage.",
306 T9,
307 rho,
308 m_cacheStats.misses(),
309 m_cacheStats.hits()
310 );
311 throw exceptions::StaleEngineError("QSE Cache Miss while lacking context for partitioning. This should be caught by the caller and trigger a re-partition stage.");
312 }
314 const auto result = m_baseEngine.getSpeciesTimescales(Y, T9, rho);
315 if (!result) {
316 return std::unexpected{result.error()};
317 }
318 std::unordered_map<Species, double> speciesTimescales = result.value();
319 for (const auto& algebraicSpecies : m_algebraic_species) {
320 speciesTimescales[algebraicSpecies] = std::numeric_limits<double>::infinity(); // Algebraic species have infinite timescales.
321 }
322 return speciesTimescales;
323 }
324
325 std::expected<std::unordered_map<fourdst::atomic::Species, double>, expectations::StaleEngineError>
327 const std::vector<double> &Y,
328 double T9,
329 double rho
330 ) const {
331 const auto key = QSECacheKey(T9, rho, Y);
332 if (!m_qse_abundance_cache.contains(key)) {
334 LOG_ERROR(
335 m_logger,
336 "QSE abundance cache miss for T9 = {}, rho = {} (misses: {}, hits: {}). getSpeciesDestructionTimescales does not receive sufficient context to partition and stabilize the network. Throwing an error which should be caught by the caller and trigger a re-partition stage.",
337 T9,
338 rho,
339 m_cacheStats.misses(),
340 m_cacheStats.hits()
341 );
342 throw exceptions::StaleEngineError("QSE Cache Miss while lacking context for partitioning. This should be caught by the caller and trigger a re-partition stage.");
343 }
345 const auto result = m_baseEngine.getSpeciesDestructionTimescales(Y, T9, rho);
346 if (!result) {
347 return std::unexpected{result.error()};
348 }
349 std::unordered_map<Species, double> speciesDestructionTimescales = result.value();
350 for (const auto& algebraicSpecies : m_algebraic_species) {
351 speciesDestructionTimescales[algebraicSpecies] = std::numeric_limits<double>::infinity(); // Algebraic species have infinite destruction timescales.
352 }
353 return speciesDestructionTimescales;
354 }
355
356 fourdst::composition::Composition MultiscalePartitioningEngineView::update(const NetIn &netIn) {
357 const fourdst::composition::Composition baseUpdatedComposition = m_baseEngine.update(netIn);
358 double T9 = netIn.temperature / 1.0e9; // Convert temperature from Kelvin to T9 (T9 = T / 1e9)
359
360 const auto preKey = QSECacheKey(
361 T9,
362 netIn.density,
363 packCompositionToVector(baseUpdatedComposition, m_baseEngine)
364 );
365 if (m_qse_abundance_cache.contains(preKey)) {
366 return baseUpdatedComposition;
367 }
368 NetIn baseUpdatedNetIn = netIn;
369 baseUpdatedNetIn.composition = baseUpdatedComposition;
370 const fourdst::composition::Composition equilibratedComposition = equilibrateNetwork(baseUpdatedNetIn);
371 std::vector<double> Y_algebraic(m_algebraic_species_indices.size(), 0.0);
372 for (size_t i = 0; i < m_algebraic_species_indices.size(); ++i) {
373 const size_t species_index = m_algebraic_species_indices[i];
374 Y_algebraic[i] = equilibratedComposition.getMolarAbundance(m_baseEngine.getNetworkSpecies()[species_index]);
375 }
376
377 // We store the algebraic abundances in the cache for both pre- and post-conditions to avoid recalculating them.
378 m_qse_abundance_cache[preKey] = Y_algebraic;
379
380 const auto postKey = QSECacheKey(
381 T9,
382 netIn.density,
383 packCompositionToVector(equilibratedComposition, m_baseEngine)
384 );
385 m_qse_abundance_cache[postKey] = Y_algebraic;
386
387 return equilibratedComposition;
388 }
389
391 const auto key = QSECacheKey(
392 netIn.temperature,
393 netIn.density,
394 packCompositionToVector(netIn.composition, m_baseEngine)
395 );
396 if (m_qse_abundance_cache.contains(key)) {
397 return m_baseEngine.isStale(netIn); // The cache hit indicates the engine is not stale for the given conditions.
398 }
399 return true;
400 }
401
403 const screening::ScreeningType model
404 ) {
405 m_baseEngine.setScreeningModel(model);
406 }
407
411
415
417 const std::vector<std::vector<size_t>> &timescale_pools,
418 const std::vector<double> &Y,
419 double T9,
420 double rho
421 ) const {
422 std::vector<std::vector<size_t>> final_connected_pools;
423
424 for (const auto& pool : timescale_pools) {
425 if (pool.empty()) {
426 continue; // Skip empty pools
427 }
428
429 // For each timescale pool, we need to analyze connectivity.
430 auto connectivity_graph = buildConnectivityGraph(pool);
431 auto components = findConnectedComponentsBFS(connectivity_graph, pool);
432 final_connected_pools.insert(final_connected_pools.end(), components.begin(), components.end());
433 }
434 return final_connected_pools;
435
436 }
437
439 const std::vector<double> &Y,
440 const double T9,
441 const double rho
442 ) {
443 // --- Step 0. Clear previous state ---
444 LOG_TRACE_L1(m_logger, "Partitioning network...");
445 LOG_TRACE_L1(m_logger, "Clearing previous state...");
446 m_qse_groups.clear();
447 m_dynamic_species.clear();
449 m_algebraic_species.clear();
451
452 // --- Step 1. Identify distinct timescale regions ---
453 LOG_TRACE_L1(m_logger, "Identifying fast reactions...");
454 const std::vector<std::vector<size_t>> timescale_pools = partitionByTimescale(Y, T9, rho);
455 LOG_TRACE_L1(m_logger, "Found {} timescale pools.", timescale_pools.size());
456
457 // --- Step 2. Select the mean slowest pool as the base dynamical group ---
458 LOG_TRACE_L1(m_logger, "Identifying mean slowest pool...");
459 const size_t mean_slowest_pool_index = identifyMeanSlowestPool(timescale_pools, Y, T9, rho);
460 LOG_TRACE_L1(m_logger, "Mean slowest pool index: {}", mean_slowest_pool_index);
461
462 // --- Step 3. Push the slowest pool into the dynamic species list ---
463 m_dynamic_species_indices = timescale_pools[mean_slowest_pool_index];
464 for (const auto& index : m_dynamic_species_indices) {
465 m_dynamic_species.push_back(m_baseEngine.getNetworkSpecies()[index]);
466 }
467
468 // --- Step 4. Pack Candidate QSE Groups ---
469 std::vector<std::vector<size_t>> candidate_pools;
470 for (size_t i = 0; i < timescale_pools.size(); ++i) {
471 if (i == mean_slowest_pool_index) continue; // Skip the slowest pool
472 LOG_TRACE_L1(m_logger, "Group {} with {} species identified for potential QSE.", i, timescale_pools[i].size());
473 candidate_pools.push_back(timescale_pools[i]);
474 }
475
476 LOG_TRACE_L1(m_logger, "Preforming connectivity analysis on timescale pools...");
477 const std::vector<std::vector<size_t>> connected_pools = analyzeTimescalePoolConnectivity(candidate_pools, Y, T9, rho);
478 LOG_TRACE_L1(m_logger, "Found {} connected pools (compared to {} timescale pools) for QSE analysis.", connected_pools.size(), timescale_pools.size());
479
480
481 // --- Step 5. Identify potential seed species for each candidate pool ---
482 LOG_TRACE_L1(m_logger, "Identifying potential seed species for candidate pools...");
483 const std::vector<QSEGroup> candidate_groups = constructCandidateGroups(connected_pools, Y, T9, rho);
484 LOG_TRACE_L1(m_logger, "Found {} candidate QSE groups for further analysis", candidate_groups.size());
485
486 LOG_TRACE_L1(m_logger, "Validating candidate groups with flux analysis...");
487 const std::vector<QSEGroup> validated_groups = validateGroupsWithFluxAnalysis(candidate_groups, Y, T9, rho);
488 LOG_TRACE_L1(
489 m_logger,
490 "Validated {} group(s) QSE groups. {}",
491 validated_groups.size(),
492 [&]() -> std::string {
493 std::stringstream ss;
494 int count = 0;
495 for (const auto& group : validated_groups) {
496 ss << "Group " << count + 1;
497 if (group.is_in_equilibrium) {
498 ss << " is in equilibrium";
499 } else {
500 ss << " is not in equilibrium";
501 }
502 if (count < validated_groups.size() - 1) {
503 ss << ", ";
504 }
505 count++;
506 }
507 return ss.str();
508 }()
509 );
510
511 m_qse_groups = std::move(validated_groups);
512 LOG_TRACE_L1(m_logger, "Identified {} QSE groups.", m_qse_groups.size());
513
514 for (const auto& group : m_qse_groups) {
515 // Add algebraic species to the algebraic set
516 for (const auto& index : group.algebraic_indices) {
517 if (std::ranges::find(m_algebraic_species_indices, index) == m_algebraic_species_indices.end()) {
518 m_algebraic_species.push_back(m_baseEngine.getNetworkSpecies()[index]);
519 m_algebraic_species_indices.push_back(index);
520
521 }
522 }
523 }
524
525 LOG_INFO(
526 m_logger,
527 "Partitioning complete. Found {} dynamic species, {} algebraic (QSE) species ({}) spread over {} QSE group{}.",
528 m_dynamic_species.size(),
529 m_algebraic_species.size(),
530 [&]() -> std::string {
531 std::stringstream ss;
532 size_t count = 0;
533 for (const auto& species : m_algebraic_species) {
534 ss << species.name();
535 if (m_algebraic_species.size() > 1 && count < m_algebraic_species.size() - 1) {
536 ss << ", ";
537 }
538 count++;
539 }
540 return ss.str();
541 }(),
542 m_qse_groups.size(),
543 m_qse_groups.size() == 1 ? "" : "s"
544 );
545
546 }
547
549 const NetIn &netIn
550 ) {
551 const std::vector<double> Y = packCompositionToVector(netIn.composition, m_baseEngine);
552 const double T9 = netIn.temperature / 1e9; // Convert temperature from Kelvin to T9 (T9 = T / 1e9)
553 const double rho = netIn.density; // Density in g/cm^3
554
555 partitionNetwork(Y, T9, rho);
556 }
557
559 const std::string &filename,
560 const std::vector<double>& Y,
561 const double T9,
562 const double rho
563 ) const {
564 std::ofstream dotFile(filename);
565 if (!dotFile.is_open()) {
566 LOG_ERROR(m_logger, "Failed to open file for writing: {}", filename);
567 throw std::runtime_error("Failed to open file for writing: " + filename);
568 }
569
570 const auto& all_species = m_baseEngine.getNetworkSpecies();
571 const auto& all_reactions = m_baseEngine.getNetworkReactions();
572
573 // --- 1. Pre-computation and Categorization ---
574
575 // Categorize species into algebraic, seed, and core dynamic
576 std::unordered_set<size_t> algebraic_indices;
577 std::unordered_set<size_t> seed_indices;
578 for (const auto& group : m_qse_groups) {
579 if (group.is_in_equilibrium) {
580 algebraic_indices.insert(group.algebraic_indices.begin(), group.algebraic_indices.end());
581 seed_indices.insert(group.seed_indices.begin(), group.seed_indices.end());
582 }
583 }
584
585 // Calculate reaction flows and find min/max for logarithmic scaling of transparency
586 std::vector<double> reaction_flows;
587 reaction_flows.reserve(all_reactions.size());
588 double min_log_flow = std::numeric_limits<double>::max();
589 double max_log_flow = std::numeric_limits<double>::lowest();
590
591 for (const auto& reaction : all_reactions) {
592 double flow = std::abs(m_baseEngine.calculateMolarReactionFlow(reaction, Y, T9, rho));
593 reaction_flows.push_back(flow);
594 if (flow > 1e-99) { // Avoid log(0)
595 double log_flow = std::log10(flow);
596 min_log_flow = std::min(min_log_flow, log_flow);
597 max_log_flow = std::max(max_log_flow, log_flow);
598 }
599 }
600 const double log_flow_range = (max_log_flow > min_log_flow) ? (max_log_flow - min_log_flow) : 1.0;
601
602 // --- 2. Write DOT file content ---
603
604 dotFile << "digraph PartitionedNetwork {\n";
605 dotFile << " graph [rankdir=TB, splines=true, overlap=false, bgcolor=\"#f8fafc\", label=\"Multiscale Partitioned Network View\", fontname=\"Helvetica\", fontsize=16, labeljust=l];\n";
606 dotFile << " node [shape=circle, style=filled, fontname=\"Helvetica\", width=0.8, fixedsize=true];\n";
607 dotFile << " edge [fontname=\"Helvetica\", fontsize=10];\n\n";
608
609 // --- Node Definitions ---
610 // Define all species nodes first, so they can be referenced by clusters and ranks later.
611 dotFile << " // --- Species Nodes Definitions ---\n";
612 std::map<int, std::vector<std::string>> species_by_mass;
613 for (size_t i = 0; i < all_species.size(); ++i) {
614 const auto& species = all_species[i];
615 std::string fillcolor = "#f1f5f9"; // Default: Other/Uninvolved
616
617 // Determine color based on category. A species can be a seed and also in the core dynamic group.
618 // The more specific category (algebraic, then seed) takes precedence.
619 if (algebraic_indices.contains(i)) {
620 fillcolor = "#e0f2fe"; // Light Blue: Algebraic (in QSE)
621 } else if (seed_indices.contains(i)) {
622 fillcolor = "#a7f3d0"; // Light Green: Seed (Dynamic, feeds a QSE group)
623 } else if (std::ranges::contains(m_dynamic_species_indices, i)) {
624 fillcolor = "#dcfce7"; // Pale Green: Core Dynamic
625 }
626 dotFile << " \"" << species.name() << "\" [label=\"" << species.name() << "\", fillcolor=\"" << fillcolor << "\"];\n";
627
628 // Group species by mass number for ranked layout.
629 // If species.a() returns incorrect values (e.g., 0 for many species), they will be grouped together here.
630 species_by_mass[species.a()].push_back(std::string(species.name()));
631 }
632 dotFile << "\n";
633
634 // --- Layout and Ranking ---
635 // Enforce a top-down layout based on mass number.
636 dotFile << " // --- Layout using Ranks ---\n";
637 for (const auto &species_list: species_by_mass | std::views::values) {
638 dotFile << " { rank=same; ";
639 for (const auto& name : species_list) {
640 dotFile << "\"" << name << "\"; ";
641 }
642 dotFile << "}\n";
643 }
644 dotFile << "\n";
645
646 // Chain by mass to get top down ordering
647 dotFile << " // --- Chain by Mass ---\n";
648 for (const auto& [mass, species_list] : species_by_mass) {
649 // Find the next largest mass in the species list
650 int minLargestMass = std::numeric_limits<int>::max();
651 for (const auto &next_mass: species_by_mass | std::views::keys) {
652 if (next_mass > mass && next_mass < minLargestMass) {
653 minLargestMass = next_mass;
654 }
655 }
656 if (minLargestMass != std::numeric_limits<int>::max()) {
657 // Connect the current mass to the next largest mass
658 dotFile << " \"" << species_list[0] << "\" -> \"" << species_by_mass[minLargestMass][0] << "\" [style=invis];\n";
659 }
660 }
661
662 // --- QSE Group Clusters ---
663 // Draw a prominent box around the algebraic species of each valid QSE group.
664 dotFile << " // --- QSE Group Clusters ---\n";
665 int group_counter = 0;
666 for (const auto& group : m_qse_groups) {
667 if (!group.is_in_equilibrium || group.algebraic_indices.empty()) {
668 continue;
669 }
670 dotFile << " subgraph cluster_qse_" << group_counter++ << " {\n";
671 dotFile << " label = \"QSE Group " << group_counter << "\";\n";
672 dotFile << " style = \"filled,rounded\";\n";
673 dotFile << " color = \"#38bdf8\";\n"; // A bright, visible blue for the border
674 dotFile << " penwidth = 2.0;\n"; // Thicker border
675 dotFile << " bgcolor = \"#f0f9ff80\";\n"; // Light blue fill with transparency
676 dotFile << " subgraph cluster_seed_" << group_counter << " {\n";
677 dotFile << " label = \"Seed Species\";\n";
678 dotFile << " style = \"filled,rounded\";\n";
679 dotFile << " color = \"#a7f3d0\";\n"; // Light green for seed species
680 dotFile << " penwidth = 1.5;\n"; // Thinner border for seed cluster
681 std::vector<std::string> seed_node_ids;
682 seed_node_ids.reserve(group.seed_indices.size());
683 for (const size_t species_idx : group.seed_indices) {
684 std::stringstream ss;
685 ss << "node_" << group_counter << "_seed_" << species_idx;
686 dotFile << " " << ss.str() << " [label=\"" << all_species[species_idx].name() << "\"];\n";
687 seed_node_ids.push_back(ss.str());
688 }
689 for (size_t i = 0; i < seed_node_ids.size() - 1; ++i) {
690 dotFile << " " << seed_node_ids[i] << " -> " << seed_node_ids[i + 1] << " [style=invis];\n";
691 }
692 dotFile << " }\n";
693 dotFile << " subgraph cluster_algebraic_" << group_counter << " {\n";
694 dotFile << " label = \"Algebraic Species\";\n";
695 dotFile << " style = \"filled,rounded\";\n";
696 dotFile << " color = \"#e0f2fe\";\n"; // Light blue for algebraic species
697 dotFile << " penwidth = 1.5;\n"; // Thinner border for algebraic cluster
698 std::vector<std::string> algebraic_node_ids;
699 algebraic_node_ids.reserve(group.algebraic_indices.size());
700 for (const size_t species_idx : group.algebraic_indices) {
701 std::stringstream ss;
702 ss << "node_" << group_counter << "_algebraic_" << species_idx;
703 dotFile << " " << ss.str() << " [label=\"" << all_species[species_idx].name() << "\"];\n";
704 algebraic_node_ids.push_back(ss.str());
705 }
706 // Make invisible edges between algebraic indices to keep them in top-down order
707 for (size_t i = 0; i < algebraic_node_ids.size() - 1; ++i) {
708 dotFile << " " << algebraic_node_ids[i] << " -> " << algebraic_node_ids[i + 1] << " [style=invis];\n";
709 }
710 dotFile << " }\n";
711 dotFile << " }\n";
712 }
713 dotFile << "\n";
714
715
716 // --- Legend ---
717 // Add a legend to explain colors and conventions.
718 dotFile << " // --- Legend ---\n";
719 dotFile << " subgraph cluster_legend {\n";
720 dotFile << " rank = sink"; // Try to push the legend to the bottom
721 dotFile << " label = \"Legend\";\n";
722 dotFile << " bgcolor = \"#ffffff\";\n";
723 dotFile << " color = \"#e2e8f0\";\n";
724 dotFile << " node [shape=box, style=filled, fontname=\"Helvetica\"];\n";
725 dotFile << " key_core [label=\"Core Dynamic\", fillcolor=\"#dcfce7\"];\n";
726 dotFile << " key_seed [label=\"Seed (Dynamic)\", fillcolor=\"#a7f3d0\"];\n";
727 dotFile << " key_qse [label=\"Algebraic (QSE)\", fillcolor=\"#e0f2fe\"];\n";
728 dotFile << " key_other [label=\"Other\", fillcolor=\"#f1f5f9\"];\n";
729 dotFile << " key_info [label=\"Edge Opacity ~ log(Reaction Flow)\", shape=plaintext];\n";
730 dotFile << " ";// Use invisible edges to stack legend items vertically
731 dotFile << " key_core -> key_seed -> key_qse -> key_other -> key_info [style=invis];\n";
732 dotFile << " }\n\n";
733
734 // --- Reaction Edges ---
735 // Draw edges with transparency scaled by the log of the molar reaction flow.
736 dotFile << " // --- Reaction Edges ---\n";
737 for (size_t i = 0; i < all_reactions.size(); ++i) {
738 const auto& reaction = all_reactions[i];
739 const double flow = reaction_flows[i];
740
741 if (flow < 1e-99) continue; // Don't draw edges for negligible flows
742
743 double log_flow_val = std::log10(flow);
744 double norm_alpha = (log_flow_val - min_log_flow) / log_flow_range;
745 int alpha_val = 0x30 + static_cast<int>(norm_alpha * (0xFF - 0x30)); // Scale from ~20% to 100% opacity
746 alpha_val = std::clamp(alpha_val, 0x00, 0xFF);
747
748 std::stringstream alpha_hex;
749 alpha_hex << std::setw(2) << std::setfill('0') << std::hex << alpha_val;
750 std::string edge_color = "#475569" + alpha_hex.str();
751
752 std::string reactionNodeId = "reaction_" + std::string(reaction.id());
753 dotFile << " \"" << reactionNodeId << "\" [shape=point, fillcolor=black, width=0.05, height=0.05];\n";
754 for (const auto& reactant : reaction.reactants()) {
755 dotFile << " \"" << reactant.name() << "\" -> \"" << reactionNodeId << "\" [color=\"" << edge_color << "\", arrowhead=none];\n";
756 }
757 for (const auto& product : reaction.products()) {
758 dotFile << " \"" << reactionNodeId << "\" -> \"" << product.name() << "\" [color=\"" << edge_color << "\"];\n";
759 }
760 dotFile << "\n";
761 }
762
763 dotFile << "}\n";
764 dotFile.close();
765 }
766
767
769 std::vector<double> Y(m_dynamic_species.size(), 0.0); // Initialize with zeros
770 for (const auto& [symbol, entry] : netIn.composition) {
771 Y[getSpeciesIndex(entry.isotope())] = netIn.composition.getMolarAbundance(symbol); // Map species to their molar abundance
772 }
773 return Y; // Return the vector of molar abundances
774 }
775
777 const auto& all_species = m_baseEngine.getNetworkSpecies();
778 std::vector<Species> fast_species;
779 fast_species.reserve(all_species.size() - m_dynamic_species.size());
780 for (const auto& species : all_species) {
781 auto it = std::ranges::find(m_dynamic_species, species);
782 if (it == m_dynamic_species.end()) {
783 fast_species.push_back(species);
784 }
785 }
786 return fast_species;
787 }
788
789 const std::vector<Species> & MultiscalePartitioningEngineView::getDynamicSpecies() const {
790 return m_dynamic_species;
791 }
792
794 return m_baseEngine.primeEngine(netIn);
795 }
796
798 const std::vector<double> &Y,
799 const double T9,
800 const double rho
801 ) {
802 partitionNetwork(Y, T9, rho);
803 const std::vector<double> Y_equilibrated = solveQSEAbundances(Y, T9, rho);
804 fourdst::composition::Composition composition;
805
806 std::vector<std::string> symbols;
807 symbols.reserve(m_baseEngine.getNetworkSpecies().size());
808 for (const auto& species : m_baseEngine.getNetworkSpecies()) {
809 symbols.emplace_back(species.name());
810 }
811 composition.registerSymbol(symbols);
812
813 std::vector<double> X;
814 X.reserve(Y_equilibrated.size());
815 for (size_t i = 0; i < Y_equilibrated.size(); ++i) {
816 const double molarMass = m_baseEngine.getNetworkSpecies()[i].mass();
817 X.push_back(Y_equilibrated[i] * molarMass); // Convert from molar abundance to mass fraction
818 }
819
820 for (size_t i = 0; i < Y_equilibrated.size(); ++i) {
821 const auto& species = m_baseEngine.getNetworkSpecies()[i];
822 if (X[i] < 0.0 && std::abs(X[i]) < 1e-20) {
823 composition.setMassFraction(std::string(species.name()), 0.0); // Avoid negative mass fractions
824 } else {
825 composition.setMassFraction(std::string(species.name()), X[i]);
826 }
827 }
828
829 composition.finalize(true);
830
831 return composition;
832 }
833
835 const NetIn &netIn
836 ) {
837 const PrimingReport primingReport = m_baseEngine.primeEngine(netIn);
838 const std::vector<double> Y = packCompositionToVector(primingReport.primedComposition, m_baseEngine);
839
840 const double T9 = netIn.temperature / 1e9; // Convert temperature from Kelvin to T9 (T9 = T / 1e9)
841 const double rho = netIn.density; // Density in g/cm^3
842
843 return equilibrateNetwork(Y, T9, rho);
844 }
845
846 int MultiscalePartitioningEngineView::getSpeciesIndex(const fourdst::atomic::Species &species) const {
847 return m_baseEngine.getSpeciesIndex(species);
848 }
849
851 const std::vector<double>& Y_full,
852 const double T9,
853 const double rho
854 ) const {
855 LOG_TRACE_L1(m_logger, "Partitioning by timescale...");
856 const auto result= m_baseEngine.getSpeciesDestructionTimescales(Y_full, T9, rho);
857 if (!result) {
858 LOG_ERROR(m_logger, "Failed to get species timescales due to stale engine state");
859 m_logger->flush_log();
860 throw exceptions::StaleEngineError("Failed to get species timescales due to stale engine state");
861 }
862 std::unordered_map<Species, double> all_timescales = result.value();
863 const auto& all_species = m_baseEngine.getNetworkSpecies();
864
865 std::vector<std::pair<double, size_t>> sorted_timescales;
866 for (size_t i = 0; i < all_species.size(); ++i) {
867 double timescale = all_timescales.at(all_species[i]);
868 if (std::isfinite(timescale) && timescale > 0) {
869 sorted_timescales.push_back({timescale, i});
870 }
871 }
872
873 std::ranges::sort(
874 sorted_timescales,
875 [](const auto& a, const auto& b)
876 {
877 return a.first > b.first;
878 }
879 );
880
881 std::vector<std::vector<size_t>> final_pools;
882 if (sorted_timescales.empty()) {
883 return final_pools;
884 }
885
886 constexpr double ABSOLUTE_QSE_TIMESCALE_THRESHOLD = 3.156e7; // Absolute threshold for QSE timescale (1 yr)
887 constexpr double MIN_GAP_THRESHOLD = 2.0; // Require a 2 order of magnitude gap
888
889 LOG_TRACE_L1(m_logger, "Found {} species with finite timescales.", sorted_timescales.size());
890 LOG_TRACE_L1(m_logger, "Absolute QSE timescale threshold: {} seconds ({} years).",
891 ABSOLUTE_QSE_TIMESCALE_THRESHOLD, ABSOLUTE_QSE_TIMESCALE_THRESHOLD / 3.156e7);
892 LOG_TRACE_L1(m_logger, "Minimum gap threshold: {} orders of magnitude.", MIN_GAP_THRESHOLD);
893
894 std::vector<size_t> dynamic_pool_indices;
895 std::vector<std::pair<double, size_t>> fast_candidates;
896
897 // 1. First Pass: Absolute Timescale Cutoff
898 for (const auto& ts_pair : sorted_timescales) {
899 if (ts_pair.first > ABSOLUTE_QSE_TIMESCALE_THRESHOLD) {
900 LOG_TRACE_L3(m_logger, "Species {} with timescale {} is considered dynamic (slower than qse timescale threshold).",
901 all_species[ts_pair.second].name(), ts_pair.first);
902 dynamic_pool_indices.push_back(ts_pair.second);
903 } else {
904 LOG_TRACE_L3(m_logger, "Species {} with timescale {} is a candidate fast species (faster than qse timescale threshold).",
905 all_species[ts_pair.second].name(), ts_pair.first);
906 fast_candidates.push_back(ts_pair);
907 }
908 }
909
910 if (!dynamic_pool_indices.empty()) {
911 LOG_TRACE_L1(m_logger, "Found {} dynamic species (slower than QSE timescale threshold).", dynamic_pool_indices.size());
912 final_pools.push_back(dynamic_pool_indices);
913 }
914
915 if (fast_candidates.empty()) {
916 LOG_TRACE_L1(m_logger, "No fast candidates found.");
917 return final_pools;
918 }
919
920 // 2. Second Pass: Gap Detection on the remaining "fast" candidates
921 std::vector<size_t> split_points;
922 for (size_t i = 0; i < fast_candidates.size() - 1; ++i) {
923 const double t1 = fast_candidates[i].first;
924 const double t2 = fast_candidates[i+1].first;
925 if (std::log10(t1) - std::log10(t2) > MIN_GAP_THRESHOLD) {
926 LOG_TRACE_L3(m_logger, "Detected gap between species {} (timescale {:0.2E}) and {} (timescale {:0.2E}).",
927 all_species[fast_candidates[i].second].name(), t1,
928 all_species[fast_candidates[i+1].second].name(), t2);
929 split_points.push_back(i + 1);
930 }
931 }
932
933 size_t last_split = 0;
934 for (const size_t split : split_points) {
935 std::vector<size_t> pool;
936 for (size_t i = last_split; i < split; ++i) {
937 pool.push_back(fast_candidates[i].second);
938 }
939 final_pools.push_back(pool);
940 last_split = split;
941 }
942
943 std::vector<size_t> final_fast_pool;
944 for (size_t i = last_split; i < fast_candidates.size(); ++i) {
945 final_fast_pool.push_back(fast_candidates[i].second);
946 }
947 final_pools.push_back(final_fast_pool);
948
949 LOG_TRACE_L1(m_logger, "Final partitioned pools: {}",
950 [&]() -> std::string {
951 std::stringstream ss;
952 int oc = 0;
953 for (const auto& pool : final_pools) {
954 ss << "[";
955 int ic = 0;
956 for (const auto& idx : pool) {
957 ss << all_species[idx].name();
958 if (ic < pool.size() - 1) {
959 ss << ", ";
960 }
961 ic++;
962 }
963 ss << "]";
964 if (oc < final_pools.size() - 1) {
965 ss << ", ";
966 }
967 oc++;
968 }
969 return ss.str();
970 }());
971 return final_pools;
972
973 }
974
975 // std::unordered_map<size_t, std::vector<size_t>> MultiscalePartitioningEngineView::buildConnectivityGraph(
976 // const std::unordered_set<size_t> &fast_reaction_indices
977 // ) const {
978 // const auto& all_reactions = m_baseEngine.getNetworkReactions();
979 // std::unordered_map<size_t, std::vector<size_t>> connectivity;
980 // for (const size_t reaction_idx : fast_reaction_indices) {
981 // const auto& reaction = all_reactions[reaction_idx];
982 // const auto& reactants = reaction.reactants();
983 // const auto& products = reaction.products();
984 //
985 // // For each fast reaction, create edges between all reactants and all products.
986 // // This represents that nucleons can flow quickly between these species.
987 // for (const auto& reactant : reactants) {
988 // const size_t reactant_idx = m_baseEngine.getSpeciesIndex(reactant);
989 // for (const auto& product : products) {
990 // const size_t product_idx = m_baseEngine.getSpeciesIndex(product);
991 //
992 // // Add a two-way edge to the adjacency list.
993 // connectivity[reactant_idx].push_back(product_idx);
994 // connectivity[product_idx].push_back(reactant_idx);
995 // }
996 // }
997 // }
998 // return connectivity;
999 // }
1000
1001 std::vector<MultiscalePartitioningEngineView::QSEGroup>
1003 const std::vector<QSEGroup> &candidate_groups,
1004 const std::vector<double> &Y,
1005 const double T9, const double rho
1006 ) const {
1007 constexpr double FLUX_RATIO_THRESHOLD = 100;
1008 std::vector<QSEGroup> validated_groups = candidate_groups;
1009 for (auto& group : validated_groups) {
1010 double internal_flux = 0.0;
1011 double external_flux = 0.0;
1012
1013 const std::unordered_set<size_t> group_members(
1014 group.species_indices.begin(),
1015 group.species_indices.end()
1016 );
1017
1018 for (const auto& reaction: m_baseEngine.getNetworkReactions()) {
1019 const double flow = std::abs(m_baseEngine.calculateMolarReactionFlow(reaction, Y, T9, rho));
1020 if (flow == 0.0) {
1021 continue; // Skip reactions with zero flow
1022 }
1023 bool has_internal_reactant = false;
1024 bool has_external_reactant = false;
1025
1026 for (const auto& reactant : reaction.reactants()) {
1027 if (group_members.contains(m_baseEngine.getSpeciesIndex(reactant))) {
1028 has_internal_reactant = true;
1029 } else {
1030 has_external_reactant = true;
1031 }
1032 }
1033
1034 bool has_internal_product = false;
1035 bool has_external_product = false;
1036
1037 for (const auto& product : reaction.products()) {
1038 if (group_members.contains(m_baseEngine.getSpeciesIndex(product))) {
1039 has_internal_product = true;
1040 } else {
1041 has_external_product = true;
1042 }
1043 }
1044
1045 // Classify the reaction based on its participants
1046 if ((has_internal_reactant && has_internal_product) && !(has_external_reactant || has_external_product)) {
1047 LOG_TRACE_L3(
1048 m_logger,
1049 "Reaction {} is internal to the group containing {} and contributes to internal flux by {}",
1050 reaction.id(),
1051 [&]() -> std::string {
1052 std::stringstream ss;
1053 int count = 0;
1054 for (const auto& idx : group.algebraic_indices) {
1055 ss << m_baseEngine.getNetworkSpecies()[idx].name();
1056 if (count < group.species_indices.size() - 1) {
1057 ss << ", ";
1058 }
1059 count++;
1060 }
1061 return ss.str();
1062 }(),
1063 flow
1064 );
1065 internal_flux += flow; // Internal flux within the group
1066 } else if ((has_internal_reactant || has_internal_product) && (has_external_reactant || has_external_product)) {
1067 LOG_TRACE_L3(
1068 m_logger,
1069 "Reaction {} is external to the group containing {} and contributes to external flux by {}",
1070 reaction.id(),
1071 [&]() -> std::string {
1072 std::stringstream ss;
1073 int count = 0;
1074 for (const auto& idx : group.algebraic_indices) {
1075 ss << m_baseEngine.getNetworkSpecies()[idx].name();
1076 if (count < group.species_indices.size() - 1) {
1077 ss << ", ";
1078 }
1079 count++;
1080 }
1081 return ss.str();
1082 }(),
1083 flow
1084 );
1085 external_flux += flow; // External flux to/from the group
1086 }
1087 // otherwise the reaction is fully decoupled from the QSE group and can be ignored.
1088 }
1089 if (internal_flux > FLUX_RATIO_THRESHOLD * external_flux) {
1090 LOG_TRACE_L1(
1091 m_logger,
1092 "Group containing {} is in equilibrium: internal flux = {}, external flux = {}, ratio = {}",
1093 [&]() -> std::string {
1094 std::stringstream ss;
1095 int count = 0;
1096 for (const auto& idx : group.algebraic_indices) {
1097 ss << m_baseEngine.getNetworkSpecies()[idx].name();
1098 if (count < group.species_indices.size() - 1) {
1099 ss << ", ";
1100 }
1101 count++;
1102 }
1103 return ss.str();
1104 }(),
1105 internal_flux,
1106 external_flux,
1107 internal_flux / external_flux
1108 );
1109 group.is_in_equilibrium = true; // This group is in equilibrium if internal flux is significantly larger than external flux.
1110 } else {
1111 LOG_TRACE_L1(
1112 m_logger,
1113 "Group containing {} is NOT in equilibrium: internal flux = {}, external flux = {}, ratio = {}",
1114 [&]() -> std::string {
1115 std::stringstream ss;
1116 int count = 0;
1117 for (const auto& idx : group.algebraic_indices) {
1118 ss << m_baseEngine.getNetworkSpecies()[idx].name();
1119 if (count < group.species_indices.size() - 1) {
1120 ss << ", ";
1121 }
1122 count++;
1123 }
1124 return ss.str();
1125 }(),
1126 internal_flux,
1127 external_flux,
1128 internal_flux / external_flux
1129 );
1130 group.is_in_equilibrium = false;
1131 }
1132 }
1133 return validated_groups;
1134 }
1135
1137 const std::vector<double> &Y_full,
1138 const double T9,
1139 const double rho
1140 ) {
1141 LOG_TRACE_L1(m_logger, "Solving for QSE abundances...");
1142 auto Y_out = Y_full;
1143
1144
1145 // Sort by timescale to solve fastest QSE groups first (these can seed slower groups)
1146 std::ranges::sort(m_qse_groups, [](const QSEGroup& a, const QSEGroup& b) {
1147 return a.mean_timescale < b.mean_timescale;
1148 });
1149
1150 for (const auto&[species_indices, is_in_equilibrium, algebraic_indices, seed_indices, mean_timescale] : m_qse_groups) {
1151 if (!is_in_equilibrium || species_indices.empty()) {
1152 LOG_TRACE_L1(
1153 m_logger,
1154 "Skipping QSE group with {} species ({} algebraic ({}), {} seeds ({})) as it is not in equilibrium.",
1155 species_indices.size(),
1156 algebraic_indices.size(),
1157 [&]() -> std::string {
1158 std::ostringstream os;
1159 int count = 0;
1160 for (const auto& idx : algebraic_indices) {
1161 os << m_baseEngine.getNetworkSpecies()[idx].name();
1162 if (count < algebraic_indices.size() - 1) {
1163 os << ", ";
1164 }
1165 count++;
1166 }
1167 return os.str();
1168 }(),
1169 seed_indices.size(),
1170 [&]() -> std::string {
1171 std::ostringstream os;
1172 int count = 0;
1173 for (const auto& idx : seed_indices) {
1174 os << m_baseEngine.getNetworkSpecies()[idx].name();
1175 if (count < seed_indices.size() - 1) {
1176 os << ", ";
1177 }
1178 count++;
1179 }
1180 return os.str();
1181 }()
1182 );
1183 continue;
1184 }
1185
1186 LOG_TRACE_L1(
1187 m_logger,
1188 "Solving for QSE abundances for group with {} species ([{}] algebraic, [{}] seeds).",
1189 species_indices.size(),
1190 [&]() -> std::string {
1191 std::stringstream ss;
1192 int count = 0;
1193 for (const auto& idx : algebraic_indices) {
1194 ss << m_baseEngine.getNetworkSpecies()[idx].name();
1195 if (count < algebraic_indices.size() - 1) {
1196 ss << ", ";
1197 }
1198 count++;
1199 }
1200 return ss.str();
1201 }(),
1202 [&]() -> std::string {
1203 std::stringstream ss;
1204 int count = 0;
1205 for (const auto& idx : seed_indices) {
1206 ss << m_baseEngine.getNetworkSpecies()[idx].name();
1207 if (count < seed_indices.size() - 1) {
1208 ss << ", ";
1209 }
1210 count++;
1211 }
1212 return ss.str();
1213 }()
1214 );
1215
1216 std::vector<size_t> qse_solve_indices;
1217 std::vector<size_t> seed_indices_vec;
1218
1219 seed_indices_vec.reserve(species_indices.size());
1220 qse_solve_indices.reserve(species_indices.size());
1221
1222 for (size_t seed_idx : seed_indices) {
1223 seed_indices_vec.emplace_back(seed_idx);
1224 }
1225
1226 for (size_t algebraic_idx : algebraic_indices) {
1227 qse_solve_indices.emplace_back(algebraic_idx);
1228 }
1229
1230 if (qse_solve_indices.empty()) continue;
1231
1232 Eigen::VectorXd Y_scale(qse_solve_indices.size());
1233 Eigen::VectorXd v_initial(qse_solve_indices.size());
1234 for (size_t i = 0; i < qse_solve_indices.size(); ++i) {
1235 constexpr double abundance_floor = 1.0e-15;
1236 const double initial_abundance = Y_full[qse_solve_indices[i]];
1237 Y_scale(i) = std::max(initial_abundance, abundance_floor);
1238 v_initial(i) = std::asinh(initial_abundance / Y_scale(i)); // Scale the initial abundances using asinh
1239 }
1240
1241 EigenFunctor functor(*this, qse_solve_indices, Y_full, T9, rho, Y_scale);
1242 Eigen::LevenbergMarquardt lm(functor);
1243 lm.parameters.ftol = 1.0e-10;
1244 lm.parameters.xtol = 1.0e-10;
1245
1246 LOG_TRACE_L1(m_logger, "Minimizing functor...");
1247 Eigen::LevenbergMarquardtSpace::Status status = lm.minimize(v_initial);
1248
1249 if (status <= 0 || status >= 4) {
1250 std::stringstream msg;
1251 msg << "QSE solver failed with status: " << status << " for group with seed ";
1252 if (seed_indices.size() == 1) {
1253 msg << "nucleus " << m_baseEngine.getNetworkSpecies()[seed_indices_vec[0]].name();
1254 } else {
1255 msg << "nuclei ";
1256 // TODO: Refactor nice list printing into utility function somewhere
1257 size_t counter = 0;
1258 for (const auto& seed_idx : seed_indices) {
1259 msg << m_baseEngine.getNetworkSpecies()[seed_idx].name();
1260 if (counter < seed_indices.size() - 2) {
1261 msg << ", ";
1262 } else if (counter == seed_indices.size() - 2) {
1263 if (seed_indices.size() < 2) {
1264 msg << " and ";
1265 } else {
1266 msg << ", and ";
1267 }
1268 }
1269 ++counter;
1270 }
1271 }
1272 throw std::runtime_error(msg.str());
1273 }
1274 LOG_TRACE_L1(m_logger, "Minimization succeeded!");
1275 Eigen::VectorXd Y_final_qse = Y_scale.array() * v_initial.array().sinh(); // Convert back to physical abundances using asinh scaling
1276 for (size_t i = 0; i < qse_solve_indices.size(); ++i) {
1277 LOG_TRACE_L1(
1278 m_logger,
1279 "Species {} (index {}) started with abundance {} and ended with {}.",
1280 m_baseEngine.getNetworkSpecies()[qse_solve_indices[i]].name(),
1281 qse_solve_indices[i],
1282 Y_full[qse_solve_indices[i]],
1283 Y_final_qse(i)
1284 );
1285 Y_out[qse_solve_indices[i]] = Y_final_qse(i);
1286 }
1287 }
1288 return Y_out;
1289 }
1290
1292 const std::vector<std::vector<size_t>> &pools,
1293 const std::vector<double> &Y,
1294 const double T9,
1295 const double rho
1296 ) const {
1297 const auto& result = m_baseEngine.getSpeciesDestructionTimescales(Y, T9, rho);
1298 if (!result) {
1299 LOG_ERROR(m_logger, "Failed to get species timescales due to stale engine state");
1300 m_logger->flush_log();
1301 throw exceptions::StaleEngineError("Failed to get species timescales due to stale engine state");
1302 }
1303 const std::unordered_map<Species, double> all_timescales = result.value();
1304
1305
1306 const auto& all_species = m_baseEngine.getNetworkSpecies();
1307
1308 size_t slowest_pool_index = 0; // Default to the first pool if no valid pool is found
1309 double slowest_mean_timescale = std::numeric_limits<double>::min();
1310 size_t count = 0;
1311 for (const auto& pool : pools) {
1312 double mean_timescale = 0.0;
1313 for (const auto& species_idx : pool) {
1314 const double timescale = all_timescales.at(all_species[species_idx]);
1315 mean_timescale += timescale;
1316 }
1317 mean_timescale = mean_timescale / pool.size();
1318 if (std::isinf(mean_timescale)) {
1319 LOG_CRITICAL(m_logger, "Encountered infinite mean timescale for pool {} with species: {}",
1320 count, [&]() -> std::string {
1321 std::stringstream ss;
1322 size_t iCount = 0;
1323 for (const auto& idx : pool) {
1324 ss << all_species[idx].name() << ": " << all_timescales.at(all_species[idx]);
1325 if (iCount < pool.size() - 1) {
1326 ss << ", ";
1327 }
1328 iCount++;
1329 }
1330 return ss.str();
1331 }()
1332 );
1333 m_logger->flush_log();
1334 throw std::logic_error("Encountered infinite mean destruction timescale for a pool while identifying the mean slowest pool set, indicating a potential issue with species timescales. Check log file for more details on specific pool composition...");
1335 }
1336 if (mean_timescale > slowest_mean_timescale) {
1337 slowest_mean_timescale = mean_timescale;
1338 slowest_pool_index = &pool - &pools[0]; // Get the index of the pool
1339 }
1340 }
1341 return slowest_pool_index;
1342 }
1343
1344 std::unordered_map<size_t, std::vector<size_t>> MultiscalePartitioningEngineView::buildConnectivityGraph(
1345 const std::vector<size_t> &species_pool
1346 ) const {
1347 std::unordered_map<size_t, std::vector<size_t>> connectivity_graph;
1348 const std::set<size_t> pool_set(species_pool.begin(), species_pool.end());
1349 const std::unordered_set<Species> pool_species = [&]() -> std::unordered_set<Species> {
1350 std::unordered_set<Species> result;
1351 for (const auto& species_idx : species_pool) {
1352 Species species = m_baseEngine.getNetworkSpecies()[species_idx];
1353 result.insert(species);
1354 }
1355 return result;
1356 }();
1357
1358 std::map<size_t, std::vector<reaction::LogicalReaction*>> speciesReactionMap;
1359 std::vector<const reaction::LogicalReaction*> candidate_reactions;
1360
1361 auto getSpeciesIdx = [&](const std::vector<Species> &species) -> std::vector<size_t> {
1362 std::vector<size_t> result;
1363 result.reserve(species.size());
1364 for (const auto& s : species) {
1365 size_t idx = m_baseEngine.getSpeciesIndex(s);
1366 result.push_back(idx);
1367 }
1368 return result;
1369 };
1370
1371 for (const auto& reaction : m_baseEngine.getNetworkReactions()) {
1372 const std::vector<Species> &reactants = reaction.reactants();
1373 const std::vector<Species> &products = reaction.products();
1374
1375 std::unordered_set<Species> reactant_set(reactants.begin(), reactants.end());
1376 std::unordered_set<Species> product_set(products.begin(), products.end());
1377
1378 // Only consider reactions where at least one distinct reactant and product are in the pool
1379 if (has_distinct_reactant_and_product_species(pool_species, reactant_set, product_set)) {
1380 std::vector<size_t> involvedIDs = getSpeciesIdx(reactants);
1381 std::vector<size_t> involvedProducts = getSpeciesIdx(products);
1382 involvedIDs.insert(involvedIDs.end(), involvedProducts.begin(), involvedProducts.end());
1383 std::set<size_t> involvedSet(involvedIDs.begin(), involvedIDs.end());
1384
1385 std::vector<size_t> intersection;
1386 intersection.reserve(involvedSet.size());
1387
1388 std::ranges::set_intersection(pool_set, involvedSet, std::back_inserter(intersection));
1389
1390 // Add clique
1391 for (const size_t& u : intersection) {
1392 for (const size_t& v : intersection) {
1393 if (u != v) { // Avoid self-loops
1394 connectivity_graph[u].push_back(v);
1395 }
1396 }
1397 }
1398 }
1399 }
1400
1401
1402 return connectivity_graph;
1403 }
1404
1405 std::vector<MultiscalePartitioningEngineView::QSEGroup> MultiscalePartitioningEngineView::constructCandidateGroups(
1406 const std::vector<std::vector<size_t>> &candidate_pools,
1407 const std::vector<double> &Y,
1408 const double T9, const double rho
1409 ) const {
1410 const auto& all_species = m_baseEngine.getNetworkSpecies();
1411 const auto& all_reactions = m_baseEngine.getNetworkReactions();
1412 const auto& result = m_baseEngine.getSpeciesDestructionTimescales(Y, T9, rho);
1413 if (!result) {
1414 LOG_ERROR(m_logger, "Failed to get species timescales due to stale engine state");
1415 m_logger->flush_log();
1416 throw exceptions::StaleEngineError("Failed to get species timescales due to stale engine state");
1417 }
1418 const std::unordered_map<Species, double> destruction_timescales = result.value();
1419
1420 std::vector<QSEGroup> candidate_groups;
1421 for (const auto& pool : candidate_pools) {
1422 if (pool.empty()) continue; // Skip empty pools
1423
1424 // For each pool first identify all topological bridge connections
1425 std::vector<std::pair<reaction::LogicalReaction, double>> bridge_reactions;
1426 for (const auto& species_idx : pool) {
1427 Species ash = all_species[species_idx];
1428 for (const auto& reaction : all_reactions) {
1429 if (reaction.contains(ash)) {
1430 // Check to make sure there is at least one reactant that is not in the pool
1431 // This lets seed nuclei bring mass into the QSE group.
1432 bool has_external_reactant = false;
1433 for (const auto& reactant : reaction.reactants()) {
1434 if (std::ranges::find(pool, m_baseEngine.getSpeciesIndex(reactant)) == pool.end()) {
1435 has_external_reactant = true;
1436 LOG_TRACE_L3(m_logger, "Found external reactant {} in reaction {} for species {}.", reactant.name(), reaction.id(), ash.name());
1437 break; // Found an external reactant, no need to check further
1438 }
1439 }
1440 if (has_external_reactant) {
1441 double flow = std::abs(m_baseEngine.calculateMolarReactionFlow(reaction, Y, T9, rho));
1442 LOG_TRACE_L3(m_logger, "Found bridge reaction {} with flow {} for species {}.", reaction.id(), flow, ash.name());
1443 bridge_reactions.push_back({reaction, flow});
1444 }
1445 }
1446 }
1447 }
1448 std::ranges::sort(
1449 bridge_reactions,
1450 [](const auto& a, const auto& b) {
1451 return a.second > b.second; // Sort by flow in descending order
1452 });
1453
1454 constexpr double MIN_GAP_THRESHOLD = 1; // Minimum logarithmic molar flow gap threshold for bridge reactions
1455 std::vector<size_t> split_points;
1456 for (size_t i = 0; i < bridge_reactions.size() - 1; ++i) {
1457 const double f1 = bridge_reactions[i].second;
1458 const double f2 = bridge_reactions[i + 1].second;
1459 if (std::log10(f1) - std::log10(f2) > MIN_GAP_THRESHOLD) {
1460 LOG_TRACE_L3(m_logger, "Detected gap between bridge reactions with flows {} and {}.", f1, f2);
1461 split_points.push_back(i + 1);
1462 }
1463 }
1464
1465 if (split_points.empty()) { // If no split points were found, we consider the whole set of bridge reactions as one group.
1466 split_points.push_back(bridge_reactions.size() - 1);
1467 }
1468
1469 std::vector<size_t> seed_indices;
1470 for (size_t i = 0; i < bridge_reactions.size(); ++i) {
1471 for (const auto& fuel : bridge_reactions[i].first.reactants()) {
1472 size_t fuel_idx = m_baseEngine.getSpeciesIndex(fuel);
1473 // Only add the fuel if it is not already in the pool
1474 if (std::ranges::find(pool, fuel_idx) == pool.end()) {
1475 seed_indices.push_back(fuel_idx);
1476 }
1477 }
1478 }
1479 std::set<size_t> all_indices(pool.begin(), pool.end());
1480 for (const auto& seed_idx : seed_indices) {
1481 all_indices.insert(seed_idx);
1482 }
1483 const std::set<size_t> poolSet(pool.begin(), pool.end());
1484 const std::set<size_t> seedSet(seed_indices.begin(), seed_indices.end());
1485
1486 double mean_timescale = 0.0;
1487 for (const auto& pool_idx : poolSet) {
1488 const auto& species = all_species[pool_idx];
1489 if (destruction_timescales.contains(species)) {
1490 mean_timescale += std::min(destruction_timescales.at(species), species.halfLife()); // Use the minimum of destruction timescale and half-life
1491 } else {
1492 mean_timescale += species.halfLife();
1493 }
1494 }
1495 mean_timescale /= poolSet.size();
1496 QSEGroup qse_group(all_indices, false, poolSet, seedSet, mean_timescale);
1497 candidate_groups.push_back(qse_group);
1498 }
1499 return candidate_groups;
1500 }
1501
1502
1505 std::vector<double> y_trial = m_Y_full_initial;
1506 Eigen::VectorXd y_qse = m_Y_scale.array() * v_qse.array().sinh(); // Convert to physical abundances using asinh scaling
1507
1508 for (size_t i = 0; i < m_qse_solve_indices.size(); ++i) {
1509 y_trial[m_qse_solve_indices[i]] = y_qse(i);
1510 }
1511
1512 const auto result = m_view->getBaseEngine().calculateRHSAndEnergy(y_trial, m_T9, m_rho);
1513 if (!result) {
1514 throw exceptions::StaleEngineError("Failed to calculate RHS and energy due to stale engine state");
1515 }
1516 const auto&[dydt, nuclearEnergyGenerationRate] = result.value();
1517 f_qse.resize(m_qse_solve_indices.size());
1518 for (size_t i = 0; i < m_qse_solve_indices.size(); ++i) {
1519 f_qse(i) = dydt[m_qse_solve_indices[i]];
1520 }
1521
1522 return 0; // Success
1523 }
1524
1526 std::vector<double> y_trial = m_Y_full_initial;
1527 Eigen::VectorXd y_qse = m_Y_scale.array() * v_qse.array().sinh(); // Convert to physical abundances using asinh scaling
1528
1529 for (size_t i = 0; i < m_qse_solve_indices.size(); ++i) {
1530 y_trial[m_qse_solve_indices[i]] = y_qse(i);
1531 }
1532
1533 // TODO: Think about if the jacobian matrix should be mutable so that generateJacobianMatrix can be const
1534 m_view->getBaseEngine().generateJacobianMatrix(y_trial, m_T9, m_rho);
1535
1536 // TODO: Think very carefully about the indices here.
1537 J_qse.resize(m_qse_solve_indices.size(), m_qse_solve_indices.size());
1538 for (size_t i = 0; i < m_qse_solve_indices.size(); ++i) {
1539 for (size_t j = 0; j < m_qse_solve_indices.size(); ++j) {
1540 J_qse(i, j) = m_view->getBaseEngine().getJacobianMatrixEntry(
1543 );
1544 }
1545 }
1546
1547 // Chain rule for asinh scaling:
1548 for (long j = 0; j < J_qse.cols(); ++j) {
1549 const double dY_dv = m_Y_scale(j) * std::cosh(v_qse(j));
1550 J_qse.col(j) *= dY_dv; // Scale the column by the derivative of the asinh scaling
1551 }
1552 return 0; // Success
1553 }
1554
1555
1557 const double T9,
1558 const double rho,
1559 const std::vector<double> &Y
1560 ) :
1561 m_T9(T9),
1562 m_rho(rho),
1563 m_Y(Y) {
1564 m_hash = hash();
1565 }
1566
1567 size_t QSECacheKey::hash() const {
1568 std::size_t seed = 0;
1569
1570 hash_combine(seed, m_Y.size());
1571
1572 hash_combine(seed, bin(m_T9, m_cacheConfig.T9_tol));
1573 hash_combine(seed, bin(m_rho, m_cacheConfig.rho_tol));
1574
1575 for (double Yi : m_Y) {
1576 if (Yi < 0.0 && std::abs(Yi) < 1e-20) {
1577 Yi = 0.0; // Avoid negative abundances
1578 } else if (Yi < 0.0 && std::abs(Yi) >= 1e-20) {
1579 throw std::invalid_argument("Yi should be positive for valid hashing (expected Yi > 0, received " + std::to_string(Yi) + ")");
1580 }
1581 hash_combine(seed, bin(Yi, m_cacheConfig.Yi_tol));
1582 }
1583
1584 return seed;
1585
1586 }
1587
1588 long QSECacheKey::bin(const double value, const double tol) {
1589 return static_cast<long>(std::floor(value / tol));
1590 }
1591
1592 bool QSECacheKey::operator==(const QSECacheKey &other) const {
1593 return m_hash == other.m_hash;
1594 }
1595
1599
1603
1607
1609 return !(*this == other);
1610 }
1611
1613 if (op == operators::All) {
1614 throw std::invalid_argument("Cannot use 'ALL' as an operator for a hit");
1615 }
1616
1617 m_hit ++;
1618 m_operatorHits[op]++;
1619 }
1621 if (op == operators::All) {
1622 throw std::invalid_argument("Cannot use 'ALL' as an operator for a miss");
1623 }
1624
1625 m_miss ++;
1626 m_operatorMisses[op]++;
1627 }
1628
1630 if (op == operators::All) {
1631 return m_hit;
1632 }
1633 return m_operatorHits.at(op);
1634 }
1635
1637 if (op == operators::All) {
1638 return m_miss;
1639 }
1640 return m_operatorMisses.at(op);
1641 }
1642
1643}
Abstract class for engines supporting Jacobian and stoichiometry operations.
A reaction network engine that uses a graph-based representation.
const std::vector< fourdst::atomic::Species > & getNetworkSpecies() const override
Gets the list of species in the network.
GraphEngine & m_baseEngine
The base engine to which this view delegates calculations.
PrimingReport primeEngine(const NetIn &netIn) override
Primes the engine with a specific species.
MultiscalePartitioningEngineView(GraphEngine &baseEngine)
Constructs a MultiscalePartitioningEngineView.
void setScreeningModel(screening::ScreeningType model) override
Sets the electron screening model.
std::vector< QSEGroup > m_qse_groups
The list of identified equilibrium groups.
const std::vector< fourdst::atomic::Species > & getDynamicSpecies() const
Gets the dynamic species in the network.
const DynamicEngine & getBaseEngine() const override
Gets the base engine.
std::vector< size_t > m_dynamic_species_indices
Indices mapping the dynamic species back to the base engine's full species list.
std::vector< double > solveQSEAbundances(const std::vector< double > &Y_full, double T9, double rho)
Solves for the QSE abundances of the algebraic species in a given state.
std::vector< fourdst::atomic::Species > getFastSpecies() const
Gets the fast species in the network.
std::vector< fourdst::atomic::Species > m_algebraic_species
Species that are treated as algebraic (in QSE) in the QSE groups.
fourdst::composition::Composition equilibrateNetwork(const std::vector< double > &Y, double T9, double rho)
Equilibrates the network by partitioning and solving for QSE abundances.
int getStoichiometryMatrixEntry(int speciesIndex, int reactionIndex) const override
Gets an entry from the stoichiometry matrix.
std::vector< size_t > m_algebraic_species_indices
Indices of algebraic species in the full network.
size_t identifyMeanSlowestPool(const std::vector< std::vector< size_t > > &pools, const std::vector< double > &Y, double T9, double rho) const
Identifies the pool with the slowest mean timescale.
std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > getSpeciesTimescales(const std::vector< double > &Y, double T9, double rho) const override
Computes timescales for all species in the network.
const std::vector< fourdst::atomic::Species > & getNetworkSpecies() const override
Gets the list of species in the network.
fourdst::composition::Composition update(const NetIn &netIn) override
Updates the internal state of the engine, performing partitioning and QSE equilibration.
std::unordered_map< QSECacheKey, std::vector< double > > m_qse_abundance_cache
Cache for QSE abundances based on T9, rho, and Y.
std::expected< StepDerivatives< double >, expectations::StaleEngineError > calculateRHSAndEnergy(const std::vector< double > &Y_full, double T9, double rho) const override
Calculates the right-hand side (dY/dt) and energy generation.
double calculateMolarReactionFlow(const reaction::Reaction &reaction, const std::vector< double > &Y_full, double T9, double rho) const override
Calculates the molar reaction flow for a given reaction.
screening::ScreeningType getScreeningModel() const override
Gets the current electron screening model.
void partitionNetwork(const std::vector< double > &Y, double T9, double rho)
Partitions the network into dynamic and algebraic (QSE) groups based on timescales.
quill::Logger * m_logger
Logger instance for logging messages.
int getSpeciesIndex(const fourdst::atomic::Species &species) const override
Gets the index of a species in the full network.
std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > getSpeciesDestructionTimescales(const std::vector< double > &Y, double T9, double rho) const override
Computes destruction timescales for all species in the network.
std::vector< QSEGroup > validateGroupsWithFluxAnalysis(const std::vector< QSEGroup > &candidate_groups, const std::vector< double > &Y, double T9, double rho) const
Validates candidate QSE groups using flux analysis.
CacheStats m_cacheStats
Statistics for the QSE abundance cache.
std::vector< double > mapNetInToMolarAbundanceVector(const NetIn &netIn) const override
Maps a NetIn struct to a molar abundance vector for the full network.
std::unordered_map< size_t, std::vector< size_t > > buildConnectivityGraph(const std::unordered_set< size_t > &fast_reaction_indices) const
Builds a connectivity graph from a set of fast reaction indices.
void generateStoichiometryMatrix() override
Generates the stoichiometry matrix for the network.
std::vector< QSEGroup > constructCandidateGroups(const std::vector< std::vector< size_t > > &candidate_pools, const std::vector< double > &Y, double T9, double rho) const
Constructs candidate QSE groups from connected timescale pools.
double getJacobianMatrixEntry(int i_full, int j_full) const override
Gets an entry from the previously generated Jacobian matrix.
void setNetworkReactions(const reaction::LogicalReactionSet &reactions) override
Sets the set of logical reactions in the network.
void generateJacobianMatrix(const std::vector< double > &Y_full, double T9, double rho) const override
Generates the Jacobian matrix for the current state.
void exportToDot(const std::string &filename, const std::vector< double > &Y, const double T9, const double rho) const
Exports the network to a DOT file for visualization.
std::vector< std::vector< size_t > > partitionByTimescale(const std::vector< double > &Y_full, double T9, double rho) const
Partitions the network by timescale.
const reaction::LogicalReactionSet & getNetworkReactions() const override
Gets the set of logical reactions in the network.
bool isStale(const NetIn &netIn) override
Checks if the engine's internal state is stale relative to the provided conditions.
std::vector< fourdst::atomic::Species > m_dynamic_species
The simplified set of species presented to the solver (the "slow" species).
std::vector< std::vector< size_t > > analyzeTimescalePoolConnectivity(const std::vector< std::vector< size_t > > &timescale_pools, const std::vector< double > &Y, double T9, double rho) const
Analyzes the connectivity of timescale pools.
Represents a single nuclear reaction from a specific data source.
Definition reaction.h:72
TemplatedReactionSet< LogicalReaction > LogicalReactionSet
A set of logical reactions.
Definition reaction.h:563
ScreeningType
Enumerates the available plasma screening models.
static int s_operator_parens_called
size_t misses(const operators op=operators::All) const
Gets the number of misses for a specific operator or all operators.
void hit(const operators op=operators::Other)
Increments the hit counter for a given operator.
size_t hits(const operators op=operators::All) const
Gets the number of hits for a specific operator or all operators.
std::map< operators,size_t > m_operatorHits
Map from operators to the number of cache hits for that operator.
void miss(const operators op=operators::Other)
Increments the miss counter for a given operator.
std::map< operators,size_t > m_operatorMisses
Map from operators to the number of cache misses for that operator.
const std::vector< double > & m_Y_full_initial
Initial abundances of all species in the full network.
Eigen::Matrix< double, Eigen::Dynamic, 1 > InputType
Eigen::Matrix< double, Eigen::Dynamic, 1 > OutputType
const std::vector< size_t > & m_qse_solve_indices
Indices of the species to solve for in the QSE group.
const double m_T9
Temperature in units of 10^9 K.
const Eigen::VectorXd & m_Y_scale
Scaling factors for the species abundances, used to improve solver stability.
int df(const InputType &v_qse, JacobianType &J_qse) const
Evaluates the Jacobian of the functor, J_qse = d(f_qse)/d(v_qse).
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > JacobianType
int operator()(const InputType &v_qse, OutputType &f_qse) const
Evaluates the functor's residual vector f_qse = dY_alg/dt.
MultiscalePartitioningEngineView * m_view
Pointer to the MultiscalePartitioningEngineView instance.
bool operator<(const QSEGroup &other) const
Less-than operator for QSEGroup, used for sorting.
bool operator>(const QSEGroup &other) const
Greater-than operator for QSEGroup.
bool operator==(const QSEGroup &other) const
Equality operator for QSEGroup.
bool operator!=(const QSEGroup &other) const
Inequality operator for QSEGroup.
double density
Density in g/cm^3.
Definition network.h:58
fourdst::composition::Composition composition
Composition of the network.
Definition network.h:54
double temperature
Temperature in Kelvin.
Definition network.h:57
Captures the result of a network priming operation.
Definition reporting.h:67
fourdst::composition::Composition primedComposition
Definition reporting.h:69
Key struct for the QSE abundance cache.
QSECacheKey(const double T9, const double rho, const std::vector< double > &Y)
Constructs a QSECacheKey.
QSECacheConfig m_cacheConfig
size_t hash() const
Computes the hash value for this key.
std::size_t m_hash
Precomputed hash value for this key.
static long bin(double value, double tol)
Converts a value to a discrete bin based on a tolerance.
bool operator==(const QSECacheKey &other) const
Equality operator for QSECacheKey.
std::vector< double > m_Y
Note that the ordering of Y must match the dynamic species indices in the view.