GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
engine_graph.cpp
Go to the documentation of this file.
3#include "gridfire/network.h"
8
9#include "fourdst/composition/species.h"
10#include "fourdst/composition/atomicSpecies.h"
11
12#include "quill/LogMacros.h"
13
14#include <cstdint>
15#include <iostream>
16#include <set>
17#include <stdexcept>
18#include <string>
19#include <string_view>
20#include <unordered_map>
21#include <utility>
22#include <vector>
23#include <fstream>
24#include <ranges>
25
26#include <boost/numeric/odeint.hpp>
27
28#include "cppad/cppad.hpp"
29#include "cppad/utility/sparse_rc.hpp"
30#include "cppad/utility/sparse_rcv.hpp"
31
32
33namespace gridfire {
35 const fourdst::composition::Composition &composition,
36 const BuildDepthType buildDepth
37 ): GraphEngine(composition, partition::GroundStatePartitionFunction(), buildDepth) {}
38
40 const fourdst::composition::Composition &composition,
41 const partition::PartitionFunction& partitionFunction,
42 const BuildDepthType buildDepth) :
43 m_reactions(build_reaclib_nuclear_network(composition, buildDepth, false)),
44 m_depth(buildDepth),
45 m_partitionFunction(partitionFunction.clone())
46 {
48 }
49
51 const reaction::LogicalReactionSet &reactions
52 ) :
53 m_reactions(reactions) {
55 }
56
58 const std::vector<double> &Y,
59 const double T9,
60 const double rho
61 ) const {
63 std::vector<double> bare_rates;
64 std::vector<double> bare_reverse_rates;
65 bare_rates.reserve(m_reactions.size());
66 bare_reverse_rates.reserve(m_reactions.size());
67
68 // TODO: Add cache to this
69 for (const auto& reaction: m_reactions) {
70 bare_rates.push_back(reaction.calculate_rate(T9));
71 bare_reverse_rates.push_back(calculateReverseRate(reaction, T9));
72 }
73
74 // --- The public facing interface can always use the precomputed version since taping is done internally ---
75 return calculateAllDerivativesUsingPrecomputation(Y, bare_rates, bare_reverse_rates, T9, rho);
76 } else {
77 return calculateAllDerivatives<double>(Y, T9, rho);
78 }
79 }
80
82 LOG_INFO(m_logger, "Synchronizing internal maps for REACLIB graph network (serif::network::GraphNetwork)...");
90
91 const size_t n = m_rhsADFun.Domain();
92 const size_t m = m_rhsADFun.Range();
93
94 const std::vector<bool> select_domain(n, true);
95 const std::vector<bool> select_range(m, true);
96
97 m_rhsADFun.subgraph_sparsity(select_domain, select_range, false, m_full_jacobian_sparsity_pattern);
98 m_jac_work.clear();
99
101 LOG_INFO(m_logger, "Internal maps synchronized. Network contains {} species and {} reactions.",
102 m_networkSpecies.size(), m_reactions.size());
103 }
104
105 // --- Network Graph Construction Methods ---
107 m_networkSpecies.clear();
108 m_networkSpeciesMap.clear();
109
110 std::set<std::string_view> uniqueSpeciesNames;
111
112 for (const auto& reaction: m_reactions) {
113 for (const auto& reactant: reaction.reactants()) {
114 uniqueSpeciesNames.insert(reactant.name());
115 }
116 for (const auto& product: reaction.products()) {
117 uniqueSpeciesNames.insert(product.name());
118 }
119 }
120
121 for (const auto& name: uniqueSpeciesNames) {
122 auto it = fourdst::atomic::species.find(std::string(name));
123 if (it != fourdst::atomic::species.end()) {
124 m_networkSpecies.push_back(it->second);
125 m_networkSpeciesMap.insert({name, it->second});
126 } else {
127 LOG_ERROR(m_logger, "Species '{}' not found in global atomic species database.", name);
128 m_logger->flush_log();
129 throw std::runtime_error("Species not found in global atomic species database: " + std::string(name));
130 }
131 }
132
133 }
134
136 LOG_TRACE_L1(m_logger, "Populating reaction ID map for REACLIB graph network (serif::network::GraphNetwork)...");
137 m_reactionIDMap.clear();
138 for (auto& reaction: m_reactions) {
139 m_reactionIDMap.emplace(reaction.id(), &reaction);
140 }
141 LOG_TRACE_L1(m_logger, "Populated {} reactions in the reaction ID map.", m_reactionIDMap.size());
142 }
143
145 m_speciesToIndexMap.clear();
146 for (size_t i = 0; i < m_networkSpecies.size(); ++i) {
148 }
149 }
150
152 // The implementation of this function (and others) constrains this nuclear network to a constant temperature and density during
153 // each evaluation.
154 const size_t numSpecies = m_networkSpecies.size();
155 m_jacobianMatrix.clear();
156 m_jacobianMatrix.resize(numSpecies, numSpecies, false); // Sparse matrix, no initial values
157 LOG_TRACE_L2(m_logger, "Jacobian matrix resized to {} rows and {} columns.",
158 m_jacobianMatrix.size1(), m_jacobianMatrix.size2());
159 }
160
161 // --- Basic Accessors and Queries ---
162 const std::vector<fourdst::atomic::Species>& GraphEngine::getNetworkSpecies() const {
163 // Returns a constant reference to the vector of unique species in the network.
164 LOG_TRACE_L3(m_logger, "Providing access to network species vector. Size: {}.", m_networkSpecies.size());
165 return m_networkSpecies;
166 }
167
169 // Returns a constant reference to the set of reactions in the network.
170 LOG_TRACE_L3(m_logger, "Providing access to network reactions set. Size: {}.", m_reactions.size());
171 return m_reactions;
172 }
173
178
179 bool GraphEngine::involvesSpecies(const fourdst::atomic::Species& species) const {
180 // Checks if a given species is present in the network's species map for efficient lookup.
181 const bool found = m_networkSpeciesMap.contains(species.name());
182 LOG_DEBUG(m_logger, "Checking if species '{}' is involved in the network: {}.", species.name(), found ? "Yes" : "No");
183 return found;
184 }
185
186 // --- Validation Methods ---
188 LOG_TRACE_L1(m_logger, "Validating mass (A) and charge (Z) conservation across all reactions in the network.");
189
190 for (const auto& reaction : m_reactions) {
191 uint64_t totalReactantA = 0;
192 uint64_t totalReactantZ = 0;
193 uint64_t totalProductA = 0;
194 uint64_t totalProductZ = 0;
195
196 // Calculate total A and Z for reactants
197 for (const auto& reactant : reaction.reactants()) {
198 auto it = m_networkSpeciesMap.find(reactant.name());
199 if (it != m_networkSpeciesMap.end()) {
200 totalReactantA += it->second.a();
201 totalReactantZ += it->second.z();
202 } else {
203 // This scenario indicates a severe data integrity issue:
204 // a reactant is part of a reaction but not in the network's species map.
205 LOG_ERROR(m_logger, "CRITICAL ERROR: Reactant species '{}' in reaction '{}' not found in network species map during conservation validation.",
206 reactant.name(), reaction.id());
207 return false;
208 }
209 }
210
211 // Calculate total A and Z for products
212 for (const auto& product : reaction.products()) {
213 auto it = m_networkSpeciesMap.find(product.name());
214 if (it != m_networkSpeciesMap.end()) {
215 totalProductA += it->second.a();
216 totalProductZ += it->second.z();
217 } else {
218 // Similar critical error for product species
219 LOG_ERROR(m_logger, "CRITICAL ERROR: Product species '{}' in reaction '{}' not found in network species map during conservation validation.",
220 product.name(), reaction.id());
221 return false;
222 }
223 }
224
225 // Compare totals for conservation
226 if (totalReactantA != totalProductA) {
227 LOG_ERROR(m_logger, "Mass number (A) not conserved for reaction '{}': Reactants A={} vs Products A={}.",
228 reaction.id(), totalReactantA, totalProductA);
229 return false;
230 }
231 if (totalReactantZ != totalProductZ) {
232 LOG_ERROR(m_logger, "Atomic number (Z) not conserved for reaction '{}': Reactants Z={} vs Products Z={}.",
233 reaction.id(), totalReactantZ, totalProductZ);
234 return false;
235 }
236 }
237
238 LOG_TRACE_L1(m_logger, "Mass (A) and charge (Z) conservation validated successfully for all reactions.");
239 return true; // All reactions passed the conservation check
240 }
241
244 const double T9
245 ) const {
247 LOG_TRACE_L3_LIMIT_EVERY_N(std::numeric_limits<int>::max(), m_logger, "Reverse reactions are disabled. Returning 0.0 for reverse rate of reaction '{}'.", reaction.id());
248 return 0.0; // If reverse reactions are not used, return 0.0
249 }
250 const double temp = T9 * 1e9; // Convert T9 to Kelvin
251
252 // In debug builds we check the units on kB to ensure it is in erg/K. This is removed in release builds to avoid overhead. (Note assert is a no-op in release builds)
253 assert(Constants::getInstance().get("kB").unit == "erg / K");
254
255 const double kBMeV = m_constants.kB * 624151; // Convert kB to MeV/K NOTE: This relies on the fact that m_constants.kB is in erg/K!
256 const double expFactor = std::exp(-reaction.qValue() / (kBMeV * temp));
257 double reverseRate = 0.0;
258 const double forwardRate = reaction.calculate_rate(T9);
259
260 if (reaction.reactants().size() == 2 && reaction.products().size() == 2) {
261 reverseRate = calculateReverseRateTwoBody(reaction, T9, forwardRate, expFactor);
262 } else {
263 LOG_WARNING_LIMIT_EVERY_N(1000000, m_logger, "Reverse rate calculation for reactions with more than two reactants or products is not implemented (reaction id {}).", reaction.peName());
264 }
265 LOG_TRACE_L2_LIMIT_EVERY_N(1000, m_logger, "Calculated reverse rate for reaction '{}': {:.3E} at T9={:.3E}.", reaction.id(), reverseRate, T9);
266 return reverseRate;
267 }
268
271 const double T9,
272 const double forwardRate,
273 const double expFactor
274 ) const {
275 std::vector<double> reactantPartitionFunctions;
276 std::vector<double> productPartitionFunctions;
277
278 reactantPartitionFunctions.reserve(reaction.reactants().size());
279 productPartitionFunctions.reserve(reaction.products().size());
280
281 std::unordered_map<fourdst::atomic::Species, int> reactantMultiplicity;
282 std::unordered_map<fourdst::atomic::Species, int> productMultiplicity;
283
284 reactantMultiplicity.reserve(reaction.reactants().size());
285 productMultiplicity.reserve(reaction.products().size());
286
287 for (const auto& reactant : reaction.reactants()) {
288 reactantMultiplicity[reactant] += 1;
289 }
290 for (const auto& product : reaction.products()) {
291 productMultiplicity[product] += 1;
292 }
293 double reactantSymmetryFactor = 1.0;
294 double productSymmetryFactor = 1.0;
295 for (const auto& count : reactantMultiplicity | std::views::values) {
296 reactantSymmetryFactor *= std::tgamma(count + 1);
297 }
298 for (const auto& count : productMultiplicity | std::views::values) {
299 productSymmetryFactor *= std::tgamma(count + 1);
300 }
301 const double symmetryFactor = reactantSymmetryFactor / productSymmetryFactor;
302
303 // Accumulate mass terms
304 auto mass_op = [](double acc, const auto& species) { return acc * species.a(); };
305 const double massNumerator = std::accumulate(
306 reaction.reactants().begin(),
307 reaction.reactants().end(),
308 1.0,
309 mass_op
310 );
311 const double massDenominator = std::accumulate(
312 reaction.products().begin(),
313 reaction.products().end(),
314 1.0,
315 mass_op
316 );
317
318 // Accumulate partition functions
319 auto pf_op = [&](double acc, const auto& species) {
320 return acc * m_partitionFunction->evaluate(species.z(), species.a(), T9);
321 };
322 const double partitionFunctionNumerator = std::accumulate(
323 reaction.reactants().begin(),
324 reaction.reactants().end(),
325 1.0,
326 pf_op
327 );
328 const double partitionFunctionDenominator = std::accumulate(
329 reaction.products().begin(),
330 reaction.products().end(),
331 1.0,
332 pf_op
333 );
334
335 const double CT = std::pow(massNumerator/massDenominator, 1.5) *
336 (partitionFunctionNumerator/partitionFunctionDenominator);
337
338 const double reverseRate = forwardRate * symmetryFactor * CT * expFactor;
339 if (!std::isfinite(reverseRate)) {
340 return 0.0; // If the reverse rate is not finite, return 0.0
341 }
342 return reverseRate; // Return the calculated reverse rate
343
344 }
345
348 const double T9,
349 const double reverseRate
350 ) const {
352 LOG_TRACE_L3_LIMIT_EVERY_N(std::numeric_limits<int>::max(), m_logger, "Reverse reactions are disabled. Returning 0.0 for reverse rate of reaction '{}'.", reaction.id());
353 return 0.0; // If reverse reactions are not used, return 0.0
354 }
355 const double d_log_kFwd = reaction.calculate_forward_rate_log_derivative(T9);
356
357 auto log_deriv_pf_op = [&](double acc, const auto& species) {
358 const double g = m_partitionFunction->evaluate(species.z(), species.a(), T9);
359 const double dg_dT = m_partitionFunction->evaluateDerivative(species.z(), species.a(), T9);
360 return (g == 0.0) ? acc : acc + (dg_dT / g);
361 };
362
363 const double reactant_log_derivative_sum = std::accumulate(
364 reaction.reactants().begin(),
365 reaction.reactants().end(),
366 0.0,
367 log_deriv_pf_op
368 );
369
370 const double product_log_derivative_sum = std::accumulate(
371 reaction.products().begin(),
372 reaction.products().end(),
373 0.0,
374 log_deriv_pf_op
375 );
376
377 const double d_log_C = reactant_log_derivative_sum - product_log_derivative_sum;
378
379 const double d_log_exp = reaction.qValue() / (m_constants.kB * T9 * T9);
380
381 const double log_total_derivative = d_log_kFwd + d_log_C + d_log_exp;
382
383 return reverseRate * log_total_derivative; // Return the derivative of the reverse rate with respect to T9
384
385 }
386
390
391 void GraphEngine::setUseReverseReactions(const bool useReverse) {
392 m_useReverseReactions = useReverse;
393 }
394
395 int GraphEngine::getSpeciesIndex(const fourdst::atomic::Species &species) const {
396 return m_speciesToIndexMap.at(species); // Returns the index of the species in the stoichiometry matrix
397 }
398
399 std::vector<double> GraphEngine::mapNetInToMolarAbundanceVector(const NetIn &netIn) const {
400 std::vector<double> Y(m_networkSpecies.size(), 0.0); // Initialize with zeros
401 for (const auto& [symbol, entry] : netIn.composition) {
402 Y[getSpeciesIndex(entry.isotope())] = netIn.composition.getMolarAbundance(symbol); // Map species to their molar abundance
403 }
404 return Y; // Return the vector of molar abundances
405 }
406
408 NetIn fullNetIn;
409 fourdst::composition::Composition composition;
410
411 std::vector<std::string> symbols;
412 symbols.reserve(m_networkSpecies.size());
413 for (const auto &symbol: m_networkSpecies) {
414 symbols.emplace_back(symbol.name());
415 }
416 composition.registerSymbol(symbols);
417 for (const auto& [symbol, entry] : netIn.composition) {
418 if (m_networkSpeciesMap.contains(symbol)) {
419 composition.setMassFraction(symbol, entry.mass_fraction());
420 } else {
421 composition.setMassFraction(symbol, 0.0);
422 }
423 }
424 composition.finalize(true);
425 fullNetIn.composition = composition;
426 fullNetIn.temperature = netIn.temperature;
427 fullNetIn.density = netIn.density;
428
429 auto primingReport = primeNetwork(fullNetIn, *this);
430
431 return primingReport;
432 }
433
435 return m_depth;
436 }
437
438 void GraphEngine::rebuild(const fourdst::composition::Composition& comp, const BuildDepthType depth) {
439 if (depth != m_depth) {
440 m_depth = depth;
442 syncInternalMaps(); // Resync internal maps after changing the depth
443 } else {
444 LOG_DEBUG(m_logger, "Rebuild requested with the same depth. No changes made to the network.");
445 }
446 }
447
449 const std::vector<double> &Y_in,
450 const std::vector<double> &bare_rates,
451 const std::vector<double> &bare_reverse_rates,
452 const double T9,
453 const double rho
454 ) const {
455 // --- Calculate screening factors ---
456 const std::vector<double> screeningFactors = m_screeningModel->calculateScreeningFactors(
459 Y_in,
460 T9,
461 rho
462 );
463
464 // --- Optimized loop ---
465 std::vector<double> molarReactionFlows;
466 molarReactionFlows.reserve(m_precomputedReactions.size());
467
468 for (const auto& precomp : m_precomputedReactions) {
469 double forwardAbundanceProduct = 1.0;
470 // bool below_threshold = false;
471 for (size_t i = 0; i < precomp.unique_reactant_indices.size(); ++i) {
472 const size_t reactantIndex = precomp.unique_reactant_indices[i];
473 const int power = precomp.reactant_powers[i];
474 // const double abundance = Y_in[reactantIndex];
475 // if (abundance < MIN_ABUNDANCE_THRESHOLD) {
476 // below_threshold = true;
477 // break;
478 // }
479
480 forwardAbundanceProduct *= std::pow(Y_in[reactantIndex], power);
481 }
482 // if (below_threshold) {
483 // molarReactionFlows.push_back(0.0);
484 // continue; // Skip this reaction if any reactant is below the abundance threshold
485 // }
486
487 const double bare_rate = bare_rates[precomp.reaction_index];
488 const double screeningFactor = screeningFactors[precomp.reaction_index];
489 const size_t numReactants = m_reactions[precomp.reaction_index].reactants().size();
490 const size_t numProducts = m_reactions[precomp.reaction_index].products().size();
491
492 const double forwardMolarReactionFlow =
493 screeningFactor *
494 bare_rate *
495 precomp.symmetry_factor *
496 forwardAbundanceProduct *
497 std::pow(rho, numReactants > 1 ? numReactants - 1 : 0.0);
498
499 double reverseMolarReactionFlow = 0.0;
500 if (precomp.reverse_symmetry_factor != 0.0 and m_useReverseReactions) {
501 const double bare_reverse_rate = bare_reverse_rates[precomp.reaction_index];
502 double reverseAbundanceProduct = 1.0;
503 for (size_t i = 0; i < precomp.unique_product_indices.size(); ++i) {
504 reverseAbundanceProduct *= std::pow(Y_in[precomp.unique_product_indices[i]], precomp.product_powers[i]);
505 }
506 reverseMolarReactionFlow = screeningFactor *
507 bare_reverse_rate *
508 precomp.reverse_symmetry_factor *
509 reverseAbundanceProduct *
510 std::pow(rho, numProducts > 1 ? numProducts - 1 : 0.0);
511 }
512
513 molarReactionFlows.push_back(forwardMolarReactionFlow - reverseMolarReactionFlow);
514
515 }
516
517 // --- Assemble molar abundance derivatives ---
519 result.dydt.assign(m_networkSpecies.size(), 0.0); // Initialize derivatives to zero
520 for (size_t j = 0; j < m_precomputedReactions.size(); ++j) {
521 const auto& precomp = m_precomputedReactions[j];
522 const double R_j = molarReactionFlows[j];
523
524 for (size_t i = 0; i < precomp.affected_species_indices.size(); ++i) {
525 const size_t speciesIndex = precomp.affected_species_indices[i];
526 const int stoichiometricCoefficient = precomp.stoichiometric_coefficients[i];
527
528 // Update the derivative for this species
529 result.dydt[speciesIndex] += static_cast<double>(stoichiometricCoefficient) * R_j;
530 }
531 }
532
533 // --- Calculate the nuclear energy generation rate ---
534 double massProductionRate = 0.0; // [mol][s^-1]
535 for (size_t i = 0; i < m_networkSpecies.size(); ++i) {
536 const auto& species = m_networkSpecies[i];
537 massProductionRate += result.dydt[i] * species.mass() * m_constants.u;
538 }
539 result.nuclearEnergyGenerationRate = -massProductionRate * m_constants.Na * m_constants.c * m_constants.c; // [erg][s^-1][g^-1]
540 return result;
541
542 }
543
544 // --- Generate Stoichiometry Matrix ---
546 LOG_TRACE_L1(m_logger, "Generating stoichiometry matrix...");
547
548 // Task 1: Set dimensions and initialize the matrix
549 size_t numSpecies = m_networkSpecies.size();
550 size_t numReactions = m_reactions.size();
551 m_stoichiometryMatrix.resize(numSpecies, numReactions, false);
552
553 LOG_TRACE_L1(m_logger, "Stoichiometry matrix initialized with dimensions: {} rows (species) x {} columns (reactions).",
554 numSpecies, numReactions);
555
556 // Task 2: Populate the stoichiometry matrix
557 // Iterate through all reactions, assign them a column index, and fill in their stoichiometric coefficients.
558 size_t reactionColumnIndex = 0;
559 for (const auto& reaction : m_reactions) {
560 // Get the net stoichiometry for the current reaction
561 std::unordered_map<fourdst::atomic::Species, int> netStoichiometry = reaction.stoichiometry();
562
563 // Iterate through the species and their coefficients in the stoichiometry map
564 for (const auto& [species, coefficient] : netStoichiometry) {
565 // Find the row index for this species
566 auto it = m_speciesToIndexMap.find(species);
567 if (it != m_speciesToIndexMap.end()) {
568 const size_t speciesRowIndex = it->second;
569 // Set the matrix element. Boost.uBLAS handles sparse insertion.
570 m_stoichiometryMatrix(speciesRowIndex, reactionColumnIndex) = coefficient;
571 } else {
572 // This scenario should ideally not happen if m_networkSpeciesMap and m_speciesToIndexMap are correctly synced
573 LOG_ERROR(m_logger, "CRITICAL ERROR: Species '{}' from reaction '{}' stoichiometry not found in species to index map.",
574 species.name(), reaction.id());
575 m_logger -> flush_log();
576 throw std::runtime_error("Species not found in species to index map: " + std::string(species.name()));
577 }
578 }
579 reactionColumnIndex++; // Move to the next column for the next reaction
580 }
581
582 LOG_TRACE_L1(m_logger, "Stoichiometry matrix population complete. Number of non-zero elements: {}.",
583 m_stoichiometryMatrix.nnz()); // Assuming nnz() exists for compressed_matrix
584 }
585
587 const std::vector<double> &Y_in,
588 const double T9,
589 const double rho
590 ) const {
591 return calculateAllDerivatives<double>(Y_in, T9, rho);
592 }
593
595 const std::vector<ADDouble> &Y_in,
596 const ADDouble &T9,
597 const ADDouble &rho
598 ) const {
599 return calculateAllDerivatives<ADDouble>(Y_in, T9, rho);
600 }
601
606
610
611 void GraphEngine::setPrecomputation(const bool precompute) {
612 m_usePrecomputation = precompute;
613 }
614
618
622
625 const std::vector<double> &Y,
626 const double T9,
627 const double rho
628 ) const {
630 }
631
633 const std::vector<double> &Y_dynamic,
634 const double T9,
635 const double rho
636 ) const {
637
638 LOG_TRACE_L1_LIMIT_EVERY_N(1000, m_logger, "Generating jacobian matrix for T9={}, rho={}..", T9, rho);
639 const size_t numSpecies = m_networkSpecies.size();
640
641 // 1. Pack the input variables into a vector for CppAD
642 std::vector<double> adInput(numSpecies + 2, 0.0); // +2 for T9 and rho
643 for (size_t i = 0; i < numSpecies; ++i) {
644 adInput[i] = std::max(Y_dynamic[i], 1e-99); // regularize the jacobian...
645 }
646 adInput[numSpecies] = T9; // T9
647 adInput[numSpecies + 1] = rho; // rho
648
649 // 2. Calculate the full jacobian
650 const std::vector<double> dotY = m_rhsADFun.Jacobian(adInput);
651
652 // 3. Pack jacobian vector into sparse matrix
653 m_jacobianMatrix.clear();
654 for (size_t i = 0; i < numSpecies; ++i) {
655 for (size_t j = 0; j < numSpecies; ++j) {
656 const double value = dotY[i * (numSpecies + 2) + j];
657 if (std::abs(value) > MIN_JACOBIAN_THRESHOLD || i == j) { // Always keep diagonal elements to avoid pathological stiffness
658 m_jacobianMatrix(i, j) = value;
659 }
660 }
661 }
662 LOG_TRACE_L1_LIMIT_EVERY_N(1000, m_logger, "Jacobian matrix generated with dimensions: {} rows x {} columns.", m_jacobianMatrix.size1(), m_jacobianMatrix.size2());
663 }
664
666 const std::vector<double> &Y_dynamic,
667 const double T9,
668 const double rho,
669 const SparsityPattern &sparsityPattern
670 ) const {
671 // --- Pack the input variables into a vector for CppAD ---
672 const size_t numSpecies = m_networkSpecies.size();
673 std::vector<double> x(numSpecies + 2, 0.0);
674 for (size_t i = 0; i < numSpecies; ++i) {
675 x[i] = Y_dynamic[i];
676 }
677 x[numSpecies] = T9;
678 x[numSpecies + 1] = rho;
679
680 // --- Convert into CppAD Sparsity pattern ---
681 const size_t nnz = sparsityPattern.size(); // Number of non-zero entries in the sparsity pattern
682 std::vector<size_t> row_indices(nnz);
683 std::vector<size_t> col_indices(nnz);
684
685 for (size_t k = 0; k < nnz; ++k) {
686 row_indices[k] = sparsityPattern[k].first;
687 col_indices[k] = sparsityPattern[k].second;
688 }
689
690 std::vector<double> values(nnz);
691 const size_t num_rows_jac = numSpecies;
692 const size_t num_cols_jac = numSpecies + 2; // +2 for T9 and rho
693
694 CppAD::sparse_rc<std::vector<size_t>> CppAD_sparsity_pattern(num_rows_jac, num_cols_jac, nnz);
695 for (size_t k = 0; k < nnz; ++k) {
696 CppAD_sparsity_pattern.set(k, sparsityPattern[k].first, sparsityPattern[k].second);
697 }
698
699 CppAD::sparse_rcv<std::vector<size_t>, std::vector<double>> jac_subset(CppAD_sparsity_pattern);
700
701 m_rhsADFun.sparse_jac_rev(
702 x,
703 jac_subset, // Sparse Jacobian output
705 "cppad",
706 m_jac_work // Work vector for CppAD
707 );
708
709 // --- Convert the sparse Jacobian back to the Boost uBLAS format ---
710 m_jacobianMatrix.clear();
711 for (size_t k = 0; k < nnz; ++k) {
712 const size_t row = jac_subset.row()[k];
713 const size_t col = jac_subset.col()[k];
714 const double value = jac_subset.val()[k];
715
716 if (std::abs(value) > MIN_JACOBIAN_THRESHOLD) {
717 m_jacobianMatrix(row, col) = value; // Insert into the sparse matrix
718 }
719 }
720 }
721
722 double GraphEngine::getJacobianMatrixEntry(const int i, const int j) const {
723 // LOG_TRACE_L3(m_logger, "Getting jacobian matrix entry for {},{} = {}", i, j, m_jacobianMatrix(i, j));
724 return m_jacobianMatrix(i, j);
725 }
726
727 std::unordered_map<fourdst::atomic::Species, int> GraphEngine::getNetReactionStoichiometry(
729 ) {
730 return reaction.stoichiometry();
731 }
732
734 const int speciesIndex,
735 const int reactionIndex
736 ) const {
737 return m_stoichiometryMatrix(speciesIndex, reactionIndex);
738 }
739
740 void GraphEngine::exportToDot(const std::string &filename) const {
741 LOG_TRACE_L1(m_logger, "Exporting network graph to DOT file: {}", filename);
742
743 std::ofstream dotFile(filename);
744 if (!dotFile.is_open()) {
745 LOG_ERROR(m_logger, "Failed to open file for writing: {}", filename);
746 m_logger->flush_log();
747 throw std::runtime_error("Failed to open file for writing: " + filename);
748 }
749
750 dotFile << "digraph NuclearReactionNetwork {\n";
751 dotFile << " graph [rankdir=LR, splines=true, overlap=false, bgcolor=\"#f0f0f0\"];\n";
752 dotFile << " node [shape=circle, style=filled, fillcolor=\"#a7c7e7\", fontname=\"Helvetica\"];\n";
753 dotFile << " edge [fontname=\"Helvetica\", fontsize=10];\n\n";
754
755 // 1. Define all species as nodes
756 dotFile << " // --- Species Nodes ---\n";
757 for (const auto& species : m_networkSpecies) {
758 dotFile << " \"" << species.name() << "\" [label=\"" << species.name() << "\"];\n";
759 }
760 dotFile << "\n";
761
762 // 2. Define all reactions as intermediate nodes and connect them
763 dotFile << " // --- Reaction Edges ---\n";
764 for (const auto& reaction : m_reactions) {
765 // Create a unique ID for the reaction node
766 std::string reactionNodeId = "reaction_" + std::string(reaction.id());
767
768 // Define the reaction node (small, black dot)
769 dotFile << " \"" << reactionNodeId << "\" [shape=point, fillcolor=black, width=0.1, height=0.1, label=\"\"];\n";
770
771 // Draw edges from reactants to the reaction node
772 for (const auto& reactant : reaction.reactants()) {
773 dotFile << " \"" << reactant.name() << "\" -> \"" << reactionNodeId << "\";\n";
774 }
775
776 // Draw edges from the reaction node to products
777 for (const auto& product : reaction.products()) {
778 dotFile << " \"" << reactionNodeId << "\" -> \"" << product.name() << "\" [label=\"" << reaction.qValue() << " MeV\"];\n";
779 }
780 dotFile << "\n";
781 }
782
783 dotFile << "}\n";
784 dotFile.close();
785 LOG_TRACE_L1(m_logger, "Successfully exported network to {}", filename);
786 }
787
788 void GraphEngine::exportToCSV(const std::string &filename) const {
789 LOG_TRACE_L1(m_logger, "Exporting network graph to CSV file: {}", filename);
790
791 std::ofstream csvFile(filename, std::ios::out | std::ios::trunc);
792 if (!csvFile.is_open()) {
793 LOG_ERROR(m_logger, "Failed to open file for writing: {}", filename);
794 m_logger->flush_log();
795 throw std::runtime_error("Failed to open file for writing: " + filename);
796 }
797 csvFile << "Reaction;Reactants;Products;Q-value;sources;rates\n";
798 for (const auto& reaction : m_reactions) {
799 // Dynamic cast to REACLIBReaction to access specific properties
800 csvFile << reaction.id() << ";";
801 // Reactants
802 size_t count = 0;
803 for (const auto& reactant : reaction.reactants()) {
804 csvFile << reactant.name();
805 if (++count < reaction.reactants().size()) {
806 csvFile << ",";
807 }
808 }
809 csvFile << ";";
810 count = 0;
811 for (const auto& product : reaction.products()) {
812 csvFile << product.name();
813 if (++count < reaction.products().size()) {
814 csvFile << ",";
815 }
816 }
817 csvFile << ";" << reaction.qValue() << ";";
818 // Reaction coefficients
819 auto sources = reaction.sources();
820 count = 0;
821 for (const auto& source : sources) {
822 csvFile << source;
823 if (++count < sources.size()) {
824 csvFile << ",";
825 }
826 }
827 csvFile << ";";
828 // Reaction coefficients
829 count = 0;
830 for (const auto& rates : reaction) {
831 csvFile << rates;
832 if (++count < reaction.size()) {
833 csvFile << ",";
834 }
835 }
836 csvFile << "\n";
837 }
838 csvFile.close();
839 LOG_TRACE_L1(m_logger, "Successfully exported network graph to {}", filename);
840 }
841
842 std::expected<std::unordered_map<fourdst::atomic::Species, double>, expectations::StaleEngineError> GraphEngine::getSpeciesTimescales(
843 const std::vector<double> &Y,
844 const double T9,
845 const double rho
846 ) const {
847 auto [dydt, _] = calculateAllDerivatives<double>(Y, T9, rho);
848 std::unordered_map<fourdst::atomic::Species, double> speciesTimescales;
849 speciesTimescales.reserve(m_networkSpecies.size());
850 for (size_t i = 0; i < m_networkSpecies.size(); ++i) {
851 double timescale = std::numeric_limits<double>::infinity();
852 const auto species = m_networkSpecies[i];
853 if (std::abs(dydt[i]) > 0.0) {
854 timescale = std::abs(Y[i] / dydt[i]);
855 }
856 speciesTimescales.emplace(species, timescale);
857 }
858 return speciesTimescales;
859 }
860
861 std::expected<std::unordered_map<fourdst::atomic::Species, double>, expectations::StaleEngineError> GraphEngine::getSpeciesDestructionTimescales(
862 const std::vector<double> &Y,
863 const double T9,
864 const double rho
865 ) const {
866 auto [dydt, _] = calculateAllDerivatives<double>(Y, T9, rho);
867 std::unordered_map<fourdst::atomic::Species, double> speciesDestructionTimescales;
868 speciesDestructionTimescales.reserve(m_networkSpecies.size());
869 for (const auto& species : m_networkSpecies) {
870 double netDestructionFlow = 0.0;
871 for (const auto& reaction : m_reactions) {
872 if (reaction.stoichiometry(species) < 0) {
873 const double flow = calculateMolarReactionFlow<double>(reaction, Y, T9, rho);
874 netDestructionFlow += flow;
875 }
876 }
877 double timescale = std::numeric_limits<double>::infinity();
878 if (netDestructionFlow != 0.0) {
879 timescale = Y[getSpeciesIndex(species)] / netDestructionFlow;
880 }
881 speciesDestructionTimescales.emplace(species, timescale);
882 }
883 return speciesDestructionTimescales;
884 }
885
886 fourdst::composition::Composition GraphEngine::update(const NetIn &netIn) {
887 fourdst::composition::Composition baseUpdatedComposition = netIn.composition;
888 for (const auto& species : m_networkSpecies) {
889 if (!netIn.composition.contains(species)) {
890 baseUpdatedComposition.registerSpecies(species);
891 baseUpdatedComposition.setMassFraction(species, 0.0);
892 }
893 }
894 baseUpdatedComposition.finalize(false);
895 return baseUpdatedComposition;
896 }
897
898 bool GraphEngine::isStale(const NetIn &netIn) {
899 return false;
900 }
901
903 LOG_TRACE_L1(m_logger, "Recording AD tape for the RHS calculation...");
904
905 // Task 1: Set dimensions and initialize the matrix
906 const size_t numSpecies = m_networkSpecies.size();
907 if (numSpecies == 0) {
908 LOG_ERROR(m_logger, "Cannot record AD tape: No species in the network.");
909 m_logger->flush_log();
910 throw std::runtime_error("Cannot record AD tape: No species in the network.");
911 }
912 const size_t numADInputs = numSpecies + 2; // Note here that by not letting T9 and rho be independent variables, we are constraining the network to a constant temperature and density during each evaluation.
913
914 // --- CppAD Tape Recording ---
915 // 1. Declare independent variable (adY)
916 // We also initialize the dummy variable for tape recording (these tell CppAD what the derivative chain looks like).
917 // Their numeric values are irrelevant except for in so far as they avoid numerical instabilities.
918
919 // Distribute total mass fraction uniformly between species in the dummy variable space
920 const auto uniformMassFraction = static_cast<CppAD::AD<double>>(1.0 / static_cast<double>(numSpecies));
921 std::vector<CppAD::AD<double>> adInput(numADInputs, uniformMassFraction);
922 adInput[numSpecies] = 1.0; // Dummy T9
923 adInput[numSpecies + 1] = 1.0; // Dummy rho
924
925 // 3. Declare independent variables (what CppAD will differentiate wrt.)
926 // This also beings the tape recording process.
927 CppAD::Independent(adInput);
928
929 std::vector<CppAD::AD<double>> adY(numSpecies);
930 for(size_t i = 0; i < numSpecies; ++i) {
931 adY[i] = adInput[i];
932 }
933 const CppAD::AD<double> adT9 = adInput[numSpecies];
934 const CppAD::AD<double> adRho = adInput[numSpecies + 1];
935
936
937 // 5. Call the actual templated function
938 // We let T9 and rho be constant, so we pass them as fixed values.
939 auto [dydt, nuclearEnergyGenerationRate] = calculateAllDerivatives<CppAD::AD<double>>(adY, adT9, adRho);
940
941 m_rhsADFun.Dependent(adInput, dydt);
942
943 LOG_TRACE_L1(m_logger, "AD tape recorded successfully for the RHS calculation. Number of independent variables: {}.",
944 adInput.size());
945 }
946
948 m_atomicReverseRates.clear();
949 m_atomicReverseRates.reserve(m_reactions.size());
950
951 for (const auto& reaction: m_reactions) {
952 if (reaction.qValue() != 0.0) {
953 m_atomicReverseRates.push_back(std::make_unique<AtomicReverseRate>(reaction, *this));
954 } else {
955 m_atomicReverseRates.push_back(nullptr);
956 }
957 }
958 }
959
961 LOG_TRACE_L1(m_logger, "Pre-computing constant components of GraphNetwork state...");
962
963 // --- Reverse map for fast species lookups ---
964 std::unordered_map<fourdst::atomic::Species, size_t> speciesIndexMap;
965 for (size_t i = 0; i < m_networkSpecies.size(); ++i) {
966 speciesIndexMap[m_networkSpecies[i]] = i;
967 }
968
970 m_precomputedReactions.reserve(m_reactions.size());
971
972 for (size_t i = 0; i < m_reactions.size(); ++i) {
973 const auto& reaction = m_reactions[i];
974 PrecomputedReaction precomp;
975 precomp.reaction_index = i;
976
977 // --- Precompute forward reaction information ---
978 // Count occurrences for each reactant to determine powers and symmetry
979 std::unordered_map<size_t, int> reactantCounts;
980 for (const auto& reactant: reaction.reactants()) {
981 size_t reactantIndex = speciesIndexMap.at(reactant);
982 reactantCounts[reactantIndex]++;
983 }
984
985 double symmetryDenominator = 1.0;
986 for (const auto& [index, count] : reactantCounts) {
987 precomp.unique_reactant_indices.push_back(index);
988 precomp.reactant_powers.push_back(count);
989
990 symmetryDenominator *= std::tgamma(count + 1);
991 }
992
993 precomp.symmetry_factor = 1.0/symmetryDenominator;
994
995 // --- Precompute reverse reaction information ---
996 if (reaction.qValue() != 0.0) {
997 std::unordered_map<size_t, int> productCounts;
998 for (const auto& product : reaction.products()) {
999 productCounts[speciesIndexMap.at(product)]++;
1000 }
1001 double reverseSymmetryDenominator = 1.0;
1002 for (const auto& [index, count] : productCounts) {
1003 precomp.unique_product_indices.push_back(index);
1004 precomp.product_powers.push_back(count);
1005 reverseSymmetryDenominator *= std::tgamma(count + 1);
1006 }
1007
1008 precomp.reverse_symmetry_factor = 1.0/reverseSymmetryDenominator;
1009 } else {
1010 precomp.unique_product_indices.clear();
1011 precomp.product_powers.clear();
1012 precomp.reverse_symmetry_factor = 0.0; // No reverse reaction for Q = 0 reactions
1013 }
1014
1015 // --- Precompute stoichiometry information ---
1016 const auto stoichiometryMap = reaction.stoichiometry();
1017 precomp.affected_species_indices.reserve(stoichiometryMap.size());
1018 precomp.stoichiometric_coefficients.reserve(stoichiometryMap.size());
1019
1020 for (const auto& [species, coeff] : stoichiometryMap) {
1021 precomp.affected_species_indices.push_back(speciesIndexMap.at(species));
1022 precomp.stoichiometric_coefficients.push_back(coeff);
1023 }
1024
1025 m_precomputedReactions.push_back(std::move(precomp));
1026 }
1027 }
1028
1030 const size_t p,
1031 const size_t q,
1032 const CppAD::vector<bool> &vx,
1033 CppAD::vector<bool> &vy,
1034 const CppAD::vector<double> &tx,
1035 CppAD::vector<double> &ty
1036 ) {
1037
1038 if ( p != 0) { return false; }
1039 const double T9 = tx[0];
1040
1041 const double reverseRate = m_engine.calculateReverseRate(m_reaction, T9);
1042 // std::cout << m_reaction.peName() << " reverseRate: " << reverseRate << " at T9: " << T9 << "\n";
1043 ty[0] = reverseRate; // Store the reverse rate in the output vector
1044
1045 if (vx.size() > 0) {
1046 vy[0] = vx[0];
1047 }
1048 return true;
1049 }
1050
1052 size_t q,
1053 const CppAD::vector<double> &tx,
1054 const CppAD::vector<double> &ty,
1055 CppAD::vector<double> &px,
1056 const CppAD::vector<double> &py
1057 ) {
1058 const double T9 = tx[0];
1059 const double reverseRate = ty[0];
1060
1061 const double derivative = m_engine.calculateReverseRateTwoBodyDerivative(m_reaction, T9, reverseRate);
1062 // std::cout << m_reaction.peName() << " reverseRate Derivative: " << derivative << "\n";
1063
1064 px[0] = py[0] * derivative; // Return the derivative of the reverse rate with respect to T9
1065
1066 return true;
1067 }
1068
1070 size_t q,
1071 const CppAD::vector<std::set<size_t>> &r,
1072 CppAD::vector<std::set<size_t>> &s
1073 ) {
1074 s[0] = r[0];
1075 return true;
1076 }
1077
1079 size_t q,
1080 const CppAD::vector<std::set<size_t>> &rt,
1081 CppAD::vector<std::set<size_t>> &st
1082 ) {
1083 st[0] = rt[0];
1084 return true;
1085 }
1086}
bool reverse(size_t q, const CppAD::vector< double > &tx, const CppAD::vector< double > &ty, CppAD::vector< double > &px, const CppAD::vector< double > &py) override
bool rev_sparse_jac(size_t q, const CppAD::vector< std::set< size_t > > &rt, CppAD::vector< std::set< size_t > > &st) override
const reaction::Reaction & m_reaction
bool forward(size_t p, size_t q, const CppAD::vector< bool > &vx, CppAD::vector< bool > &vy, const CppAD::vector< double > &tx, CppAD::vector< double > &ty) override
bool for_sparse_jac(size_t q, const CppAD::vector< std::set< size_t > > &r, CppAD::vector< std::set< size_t > > &s) override
bool isPrecomputationEnabled() const
double calculateReverseRateTwoBody(const reaction::Reaction &reaction, const double T9, const double forwardRate, const double expFactor) const
double calculateReverseRate(const reaction::Reaction &reaction, double T9) const
const std::vector< fourdst::atomic::Species > & getNetworkSpecies() const override
Gets the list of species in the network.
BuildDepthType getDepth() const override
bool m_usePrecomputation
Flag to enable or disable using precomputed reactions for efficiency. Mathematically,...
CppAD::sparse_rc< std::vector< size_t > > m_full_jacobian_sparsity_pattern
Full sparsity pattern for the Jacobian matrix.
CppAD::sparse_jac_work m_jac_work
Work object for sparse Jacobian calculations.
void populateReactionIDMap()
Populates the reaction ID map.
std::vector< double > mapNetInToMolarAbundanceVector(const NetIn &netIn) const override
void collectAtomicReverseRateAtomicBases()
CppAD::ADFun< double > m_rhsADFun
CppAD function for the right-hand side of the ODE.
boost::numeric::ublas::compressed_matrix< double > m_jacobianMatrix
Jacobian matrix (species x species).
double getJacobianMatrixEntry(const int i, const int j) const override
Gets an entry from the previously generated Jacobian matrix.
std::unordered_map< std::string_view, fourdst::atomic::Species > m_networkSpeciesMap
Map from species name to Species object.
bool m_useReverseReactions
Flag to enable or disable reverse reactions. If false, only forward reactions are considered.
std::unique_ptr< partition::PartitionFunction > m_partitionFunction
Partition function for the network.
void setNetworkReactions(const reaction::LogicalReactionSet &reactions) override
void setUseReverseReactions(bool useReverse)
void populateSpeciesToIndexMap()
Populates the species-to-index map.
quill::Logger * m_logger
screening::ScreeningType m_screeningType
Screening type for the reaction network. Default to no screening.
fourdst::composition::Composition update(const NetIn &netIn) override
Update the internal state of the engine.
std::vector< PrecomputedReaction > m_precomputedReactions
Precomputed reactions for efficiency.
std::unordered_map< std::string_view, reaction::Reaction * > m_reactionIDMap
Map from reaction ID to REACLIBReaction. //PERF: This makes copies of REACLIBReaction and could be a ...
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.
screening::ScreeningType getScreeningModel() const override
Get the current electron screening model.
int getStoichiometryMatrixEntry(const int speciesIndex, const int reactionIndex) const override
Gets an entry from the stoichiometry matrix.
void setPrecomputation(bool precompute)
BuildDepthType m_depth
void setScreeningModel(screening::ScreeningType) override
Set the electron screening model.
std::vector< std::unique_ptr< AtomicReverseRate > > m_atomicReverseRates
void exportToCSV(const std::string &filename) const
Exports the network to a CSV file for analysis.
static std::unordered_map< fourdst::atomic::Species, int > getNetReactionStoichiometry(const reaction::Reaction &reaction)
Gets the net stoichiometry for a given reaction.
void reserveJacobianMatrix() const
Reserves space for the Jacobian matrix.
int getSpeciesIndex(const fourdst::atomic::Species &species) const override
double calculateMolarReactionFlow(const reaction::Reaction &reaction, const std::vector< double > &Y, const double T9, const double rho) const override
Calculates the molar reaction flow for a given reaction.
std::vector< fourdst::atomic::Species > m_networkSpecies
Vector of unique species in the network.
void recordADTape()
Records the AD tape for the right-hand side of the ODE.
StepDerivatives< double > calculateAllDerivativesUsingPrecomputation(const std::vector< double > &Y_in, const std::vector< double > &bare_rates, const std::vector< double > &bare_reverse_rates, double T9, double rho) const
bool involvesSpecies(const fourdst::atomic::Species &species) const
Checks if a given species is involved in the network.
std::expected< StepDerivatives< double >, expectations::StaleEngineError > calculateRHSAndEnergy(const std::vector< double > &Y, const double T9, const double rho) const override
Calculates the right-hand side (dY/dt) and energy generation rate.
reaction::LogicalReactionSet m_reactions
Set of REACLIB reactions in the network.
void syncInternalMaps()
Synchronizes the internal maps.
bool validateConservation() const
Validates mass and charge conservation across all reactions.
void generateJacobianMatrix(const std::vector< double > &Y_dynamic, const double T9, const double rho) const override
Generates the Jacobian matrix for the current state.
boost::numeric::ublas::compressed_matrix< int > m_stoichiometryMatrix
Stoichiometry matrix (species x reactions).
const reaction::LogicalReactionSet & getNetworkReactions() const override
Gets the set of logical reactions in the network.
std::unordered_map< fourdst::atomic::Species, size_t > m_speciesToIndexMap
Map from species to their index in the stoichiometry matrix.
void rebuild(const fourdst::composition::Composition &comp, const BuildDepthType depth) override
void exportToDot(const std::string &filename) const
Exports the network to a DOT file for visualization.
const partition::PartitionFunction & getPartitionFunction() const
bool isUsingReverseReactions() const
PrimingReport primeEngine(const NetIn &netIn) override
void generateStoichiometryMatrix() override
Generates the stoichiometry matrix for the network.
void collectNetworkSpecies()
Collects the unique species in the network.
std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > getSpeciesDestructionTimescales(const std::vector< double > &Y, double T9, double rho) const override
bool isStale(const NetIn &netIn) override
std::unique_ptr< screening::ScreeningModel > m_screeningModel
double calculateReverseRateTwoBodyDerivative(const reaction::Reaction &reaction, const double T9, const double reverseRate) const
StepDerivatives< T > calculateAllDerivatives(const std::vector< T > &Y_in, T T9, T rho) const
Calculates all derivatives (dY/dt) and the energy generation rate.
GraphEngine(const fourdst::composition::Composition &composition, const BuildDepthType=NetworkBuildDepth::Full)
Constructs a GraphEngine from a composition.
Abstract interface for evaluating nuclear partition functions.
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
std::unique_ptr< ScreeningModel > selectScreeningModel(ScreeningType type)
A factory function to select and create a screening model.
ScreeningType
Enumerates the available plasma screening models.
CppAD::AD< double > ADDouble
Alias for CppAD AD type for double precision.
std::variant< NetworkBuildDepth, int > BuildDepthType
Variant specifying either a predefined NetworkBuildDepth or a custom integer depth.
Definition building.h:37
std::vector< std::pair< size_t, size_t > > SparsityPattern
PrimingReport primeNetwork(const NetIn &, DynamicEngine &engine)
Primes absent species in the network to their equilibrium abundances.
Definition priming.cpp:39
static constexpr double MIN_JACOBIAN_THRESHOLD
Minimum value for Jacobian matrix entries.
reaction::LogicalReactionSet build_reaclib_nuclear_network(const fourdst::composition::Composition &composition, BuildDepthType maxLayers=NetworkBuildDepth::Full, bool reverse=false)
Builds a nuclear reaction network from the Reaclib library based on an initial composition.
Defines classes for representing and managing nuclear reactions.
std::vector< int > product_powers
Powers of each unique product in the reverse reaction.
double reverse_symmetry_factor
Symmetry factor for reverse reactions.
std::vector< size_t > unique_product_indices
Unique product indices for reverse reactions.
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
Structure holding derivatives and energy generation for a network step.
T nuclearEnergyGenerationRate
Specific energy generation rate (e.g., erg/g/s).
std::vector< T > dydt
Derivatives of abundances (dY/dt for each species).