feat(reverse-reactions): major work towrds detailed balance calculations
This commit is contained in:
@@ -10,6 +10,8 @@
|
|||||||
#include "gridfire/engine/engine_abstract.h"
|
#include "gridfire/engine/engine_abstract.h"
|
||||||
#include "gridfire/screening/screening_abstract.h"
|
#include "gridfire/screening/screening_abstract.h"
|
||||||
#include "gridfire/screening/screening_types.h"
|
#include "gridfire/screening/screening_types.h"
|
||||||
|
#include "gridfire/partition/partition_abstract.h"
|
||||||
|
#include "gridfire/partition/partition_ground.h"
|
||||||
|
|
||||||
#include <string>
|
#include <string>
|
||||||
#include <unordered_map>
|
#include <unordered_map>
|
||||||
@@ -101,6 +103,10 @@ namespace gridfire {
|
|||||||
*/
|
*/
|
||||||
explicit GraphEngine(const fourdst::composition::Composition &composition);
|
explicit GraphEngine(const fourdst::composition::Composition &composition);
|
||||||
|
|
||||||
|
explicit GraphEngine(const fourdst::composition::Composition &composition,
|
||||||
|
const partition::PartitionFunction& partitionFunction
|
||||||
|
);
|
||||||
|
|
||||||
/**
|
/**
|
||||||
* @brief Constructs a GraphEngine from a set of reactions.
|
* @brief Constructs a GraphEngine from a set of reactions.
|
||||||
*
|
*
|
||||||
@@ -307,6 +313,14 @@ namespace gridfire {
|
|||||||
|
|
||||||
[[nodiscard]] bool isPrecomputationEnabled() const;
|
[[nodiscard]] bool isPrecomputationEnabled() const;
|
||||||
|
|
||||||
|
[[nodiscard]] const partition::PartitionFunction& getPartitionFunction() const;
|
||||||
|
|
||||||
|
[[nodiscard]] double calculateReverseRate(
|
||||||
|
const reaction::Reaction &reaction,
|
||||||
|
double T9,
|
||||||
|
double expFactor
|
||||||
|
) const;
|
||||||
|
|
||||||
private:
|
private:
|
||||||
struct PrecomputedReaction {
|
struct PrecomputedReaction {
|
||||||
size_t reaction_index;
|
size_t reaction_index;
|
||||||
@@ -321,6 +335,7 @@ namespace gridfire {
|
|||||||
const double u = Constants::getInstance().get("u").value; ///< Atomic mass unit in g.
|
const double u = Constants::getInstance().get("u").value; ///< Atomic mass unit in g.
|
||||||
const double Na = Constants::getInstance().get("N_a").value; ///< Avogadro's number.
|
const double Na = Constants::getInstance().get("N_a").value; ///< Avogadro's number.
|
||||||
const double c = Constants::getInstance().get("c").value; ///< Speed of light in cm/s.
|
const double c = Constants::getInstance().get("c").value; ///< Speed of light in cm/s.
|
||||||
|
const double kB = Constants::getInstance().get("kB").value; ///< Boltzmann constant in erg/K.
|
||||||
};
|
};
|
||||||
|
|
||||||
private:
|
private:
|
||||||
@@ -347,6 +362,7 @@ namespace gridfire {
|
|||||||
bool m_usePrecomputation = true; ///< Flag to enable or disable using precomputed reactions for efficiency. Mathematically, this should not change the results. Generally end users should not need to change this.
|
bool m_usePrecomputation = true; ///< Flag to enable or disable using precomputed reactions for efficiency. Mathematically, this should not change the results. Generally end users should not need to change this.
|
||||||
|
|
||||||
std::vector<PrecomputedReaction> m_precomputedReactions; ///< Precomputed reactions for efficiency.
|
std::vector<PrecomputedReaction> m_precomputedReactions; ///< Precomputed reactions for efficiency.
|
||||||
|
std::unique_ptr<partition::PartitionFunction> m_partitionFunction; ///< Partition function for the network.
|
||||||
|
|
||||||
private:
|
private:
|
||||||
/**
|
/**
|
||||||
@@ -432,6 +448,20 @@ namespace gridfire {
|
|||||||
double T9
|
double T9
|
||||||
);
|
);
|
||||||
|
|
||||||
|
|
||||||
|
double calculateReverseRateTwoBody(
|
||||||
|
const reaction::Reaction &reaction,
|
||||||
|
const double T9,
|
||||||
|
const double forwardRate,
|
||||||
|
const double expFactor
|
||||||
|
) const;
|
||||||
|
|
||||||
|
double GraphEngine::calculateReverseRateTwoBodyDerivative(
|
||||||
|
const reaction::Reaction &reaction,
|
||||||
|
const double T9,
|
||||||
|
const double reverseRate
|
||||||
|
) const;
|
||||||
|
|
||||||
[[nodiscard]] StepDerivatives<double> calculateAllDerivativesUsingPrecomputation(
|
[[nodiscard]] StepDerivatives<double> calculateAllDerivativesUsingPrecomputation(
|
||||||
const std::vector<double> &Y_in,
|
const std::vector<double> &Y_in,
|
||||||
const std::vector<double>& bare_rates,
|
const std::vector<double>& bare_rates,
|
||||||
@@ -647,4 +677,36 @@ namespace gridfire {
|
|||||||
// the entire network.
|
// the entire network.
|
||||||
return molar_concentration_product * k_reaction * threshold_flag;
|
return molar_concentration_product * k_reaction * threshold_flag;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
class AtomicReverseRate final : public CppAD::atomic_base<double> {
|
||||||
|
public:
|
||||||
|
AtomicReverseRate(
|
||||||
|
const reaction::Reaction& reaction,
|
||||||
|
const GraphEngine& engine
|
||||||
|
):
|
||||||
|
atomic_base<double>("AtomicReverseRate"),
|
||||||
|
m_reaction(reaction),
|
||||||
|
m_engine(engine) {}
|
||||||
|
|
||||||
|
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 reverse(
|
||||||
|
size_t id,
|
||||||
|
size_t an,
|
||||||
|
const CppAD::vector<double>& tx,
|
||||||
|
const CppAD::vector<double>& ty,
|
||||||
|
CppAD::vector<double>& px,
|
||||||
|
const CppAD::vector<double>& py
|
||||||
|
);
|
||||||
|
private:
|
||||||
|
const double m_kB = Constants::getInstance().get("k_b").value; ///< Boltzmann constant in erg/K.
|
||||||
|
const reaction::Reaction& m_reaction;
|
||||||
|
const GraphEngine& m_engine;
|
||||||
|
};
|
||||||
};
|
};
|
||||||
@@ -1,6 +1,7 @@
|
|||||||
#pragma once
|
#pragma once
|
||||||
|
|
||||||
#include "gridfire/partition/partition_abstract.h"
|
#include "gridfire/partition/partition_abstract.h"
|
||||||
|
#include "gridfire/partition/partition_types.h"
|
||||||
|
|
||||||
#include "fourdst/logging/logging.h"
|
#include "fourdst/logging/logging.h"
|
||||||
|
|
||||||
@@ -12,28 +13,17 @@
|
|||||||
|
|
||||||
|
|
||||||
namespace gridfire::partition {
|
namespace gridfire::partition {
|
||||||
enum BasePartitionType {
|
|
||||||
RauscherThielemann, ///< Rauscher-Thielemann partition function
|
|
||||||
GroundState, ///< Ground state partition function
|
|
||||||
};
|
|
||||||
|
|
||||||
inline std::unordered_map<BasePartitionType, std::string> basePartitionTypeToString = {
|
|
||||||
{RauscherThielemann, "RauscherThielemann"},
|
|
||||||
{GroundState, "GroundState"}
|
|
||||||
};
|
|
||||||
|
|
||||||
inline std::unordered_map<std::string, BasePartitionType> stringToBasePartitionType = {
|
|
||||||
{"RauscherThielemann", RauscherThielemann},
|
|
||||||
{"GroundState", GroundState}
|
|
||||||
};
|
|
||||||
|
|
||||||
class CompositePartitionFunction final : public PartitionFunction {
|
class CompositePartitionFunction final : public PartitionFunction {
|
||||||
public:
|
public:
|
||||||
explicit CompositePartitionFunction(const std::vector<BasePartitionType>& partitionFunctions);
|
explicit CompositePartitionFunction(const std::vector<BasePartitionType>& partitionFunctions);
|
||||||
double evaluate(int z, int a, double T9) const override;
|
CompositePartitionFunction(const CompositePartitionFunction& other);
|
||||||
double evaluateDerivative(int z, int a, double T9) const override;
|
[[nodiscard]] double evaluate(int z, int a, double T9) const override;
|
||||||
bool supports(int z, int a) const override;
|
[[nodiscard]] double evaluateDerivative(int z, int a, double T9) const override;
|
||||||
std::string type() const override;
|
[[nodiscard]] bool supports(int z, int a) const override;
|
||||||
|
[[nodiscard]] std::string type() const override;
|
||||||
|
[[nodiscard]] std::unique_ptr<PartitionFunction> clone() const override {
|
||||||
|
return std::make_unique<CompositePartitionFunction>(*this);
|
||||||
|
}
|
||||||
private:
|
private:
|
||||||
quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
|
quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
|
||||||
std::vector<std::unique_ptr<PartitionFunction>> m_partitionFunctions; ///< Set of partition functions to use in the composite partition function.
|
std::vector<std::unique_ptr<PartitionFunction>> m_partitionFunctions; ///< Set of partition functions to use in the composite partition function.
|
||||||
|
|||||||
@@ -1,14 +1,16 @@
|
|||||||
#pragma once
|
#pragma once
|
||||||
|
|
||||||
#include <string>
|
#include <string>
|
||||||
|
#include <memory>
|
||||||
|
|
||||||
namespace gridfire::partition {
|
namespace gridfire::partition {
|
||||||
class PartitionFunction {
|
class PartitionFunction {
|
||||||
public:
|
public:
|
||||||
virtual ~PartitionFunction() = default;
|
virtual ~PartitionFunction() = default;
|
||||||
virtual double evaluate(int z, int a, double T9) const = 0;
|
[[nodiscard]] virtual double evaluate(int z, int a, double T9) const = 0;
|
||||||
virtual double evaluateDerivative(int z, int a, double T9) const = 0;
|
[[nodiscard]] virtual double evaluateDerivative(int z, int a, double T9) const = 0;
|
||||||
virtual bool supports(int z, int a) const = 0;
|
[[nodiscard]] virtual bool supports(int z, int a) const = 0;
|
||||||
virtual std::string type() const = 0;
|
[[nodiscard]] virtual std::string type() const = 0;
|
||||||
|
[[nodiscard]] virtual std::unique_ptr<PartitionFunction> clone() const = 0;
|
||||||
};
|
};
|
||||||
}
|
}
|
||||||
@@ -5,6 +5,7 @@
|
|||||||
#include "fourdst/logging/logging.h"
|
#include "fourdst/logging/logging.h"
|
||||||
|
|
||||||
#include <unordered_map>
|
#include <unordered_map>
|
||||||
|
#include <memory>
|
||||||
|
|
||||||
#include "quill/Logger.h"
|
#include "quill/Logger.h"
|
||||||
|
|
||||||
@@ -27,6 +28,9 @@ namespace gridfire::partition {
|
|||||||
const int a
|
const int a
|
||||||
) const override;
|
) const override;
|
||||||
std::string type() const override { return "GroundState"; }
|
std::string type() const override { return "GroundState"; }
|
||||||
|
std::unique_ptr<PartitionFunction> clone() const override {
|
||||||
|
return std::make_unique<GroundStatePartitionFunction>(*this);
|
||||||
|
}
|
||||||
private:
|
private:
|
||||||
quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
|
quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
|
||||||
std::unordered_map<int, double> m_ground_state_spin;
|
std::unordered_map<int, double> m_ground_state_spin;
|
||||||
|
|||||||
@@ -9,6 +9,7 @@
|
|||||||
#include <unordered_map>
|
#include <unordered_map>
|
||||||
#include <array>
|
#include <array>
|
||||||
#include <utility>
|
#include <utility>
|
||||||
|
#include <memory>
|
||||||
|
|
||||||
namespace gridfire::partition {
|
namespace gridfire::partition {
|
||||||
class RauscherThielemannPartitionFunction final : public PartitionFunction {
|
class RauscherThielemannPartitionFunction final : public PartitionFunction {
|
||||||
@@ -41,6 +42,9 @@ namespace gridfire::partition {
|
|||||||
size_t upperIndex;
|
size_t upperIndex;
|
||||||
size_t lowerIndex;
|
size_t lowerIndex;
|
||||||
};
|
};
|
||||||
|
std::unique_ptr<PartitionFunction> clone() const override {
|
||||||
|
return std::make_unique<RauscherThielemannPartitionFunction>(*this);
|
||||||
|
}
|
||||||
private:
|
private:
|
||||||
quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
|
quill::Logger* m_logger = fourdst::logging::LogManager::getInstance().getLogger("log");
|
||||||
std::unordered_map<int, IsotopeData> m_partitionData;
|
std::unordered_map<int, IsotopeData> m_partitionData;
|
||||||
|
|||||||
21
src/network/include/gridfire/partition/partition_types.h
Normal file
21
src/network/include/gridfire/partition/partition_types.h
Normal file
@@ -0,0 +1,21 @@
|
|||||||
|
#pragma once
|
||||||
|
|
||||||
|
#include <unordered_map>
|
||||||
|
#include <string>
|
||||||
|
|
||||||
|
namespace gridfire::partition {
|
||||||
|
enum BasePartitionType {
|
||||||
|
RauscherThielemann, ///< Rauscher-Thielemann partition function
|
||||||
|
GroundState, ///< Ground state partition function
|
||||||
|
};
|
||||||
|
|
||||||
|
inline std::unordered_map<BasePartitionType, std::string> basePartitionTypeToString = {
|
||||||
|
{RauscherThielemann, "RauscherThielemann"},
|
||||||
|
{GroundState, "GroundState"}
|
||||||
|
};
|
||||||
|
|
||||||
|
inline std::unordered_map<std::string, BasePartitionType> stringToBasePartitionType = {
|
||||||
|
{"RauscherThielemann", RauscherThielemann},
|
||||||
|
{"GroundState", GroundState}
|
||||||
|
};
|
||||||
|
}
|
||||||
@@ -8,7 +8,6 @@ namespace gridfire::partition::record {
|
|||||||
uint32_t z; ///< Atomic number
|
uint32_t z; ///< Atomic number
|
||||||
uint32_t a; ///< Mass number
|
uint32_t a; ///< Mass number
|
||||||
double ground_state_spin; ///< Ground state spin
|
double ground_state_spin; ///< Ground state spin
|
||||||
double partition_function; ///< Partition function value
|
|
||||||
double normalized_g_values[24]; ///< Normalized g-values for the first 24 energy levels
|
double normalized_g_values[24]; ///< Normalized g-values for the first 24 energy levels
|
||||||
};
|
};
|
||||||
#pragma pack(pop)
|
#pragma pack(pop)
|
||||||
|
|||||||
@@ -113,6 +113,8 @@ namespace gridfire::reaction {
|
|||||||
*/
|
*/
|
||||||
[[nodiscard]] virtual CppAD::AD<double> calculate_rate(const CppAD::AD<double> T9) const;
|
[[nodiscard]] virtual CppAD::AD<double> calculate_rate(const CppAD::AD<double> T9) const;
|
||||||
|
|
||||||
|
[[nodiscard]] virtual double calculate_forward_rate_log_derivative(const double T9) const;
|
||||||
|
|
||||||
/**
|
/**
|
||||||
* @brief Gets the reaction name in (projectile, ejectile) notation.
|
* @brief Gets the reaction name in (projectile, ejectile) notation.
|
||||||
* @return The reaction name (e.g., "p(p,g)d").
|
* @return The reaction name (e.g., "p(p,g)d").
|
||||||
@@ -341,6 +343,8 @@ namespace gridfire::reaction {
|
|||||||
*/
|
*/
|
||||||
[[nodiscard]] double calculate_rate(const double T9) const override;
|
[[nodiscard]] double calculate_rate(const double T9) const override;
|
||||||
|
|
||||||
|
[[nodiscard]] virtual double calculate_forward_rate_log_derivative(const double T9) const override;
|
||||||
|
|
||||||
/**
|
/**
|
||||||
* @brief Calculates the total reaction rate using CppAD types.
|
* @brief Calculates the total reaction rate using CppAD types.
|
||||||
* @param T9 The temperature in units of 10^9 K, as a CppAD::AD<double>.
|
* @param T9 The temperature in units of 10^9 K, as a CppAD::AD<double>.
|
||||||
|
|||||||
@@ -18,6 +18,7 @@
|
|||||||
#include <utility>
|
#include <utility>
|
||||||
#include <vector>
|
#include <vector>
|
||||||
#include <fstream>
|
#include <fstream>
|
||||||
|
#include <ranges>
|
||||||
|
|
||||||
#include <boost/numeric/odeint.hpp>
|
#include <boost/numeric/odeint.hpp>
|
||||||
|
|
||||||
@@ -26,7 +27,18 @@ namespace gridfire {
|
|||||||
GraphEngine::GraphEngine(
|
GraphEngine::GraphEngine(
|
||||||
const fourdst::composition::Composition &composition
|
const fourdst::composition::Composition &composition
|
||||||
):
|
):
|
||||||
m_reactions(build_reaclib_nuclear_network(composition, false)) {
|
m_reactions(build_reaclib_nuclear_network(composition, false)),
|
||||||
|
m_partitionFunction(std::make_unique<partition::GroundStatePartitionFunction>()){
|
||||||
|
syncInternalMaps();
|
||||||
|
precomputeNetwork();
|
||||||
|
}
|
||||||
|
|
||||||
|
GraphEngine::GraphEngine(
|
||||||
|
const fourdst::composition::Composition &composition,
|
||||||
|
const partition::PartitionFunction& partitionFunction) :
|
||||||
|
m_reactions(build_reaclib_nuclear_network(composition, false)),
|
||||||
|
m_partitionFunction(partitionFunction.clone()) // Clone the partition function to ensure ownership
|
||||||
|
{
|
||||||
syncInternalMaps();
|
syncInternalMaps();
|
||||||
precomputeNetwork();
|
precomputeNetwork();
|
||||||
}
|
}
|
||||||
@@ -220,6 +232,129 @@ namespace gridfire {
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
double GraphEngine::calculateReverseRate(
|
||||||
|
const reaction::Reaction &reaction,
|
||||||
|
const double T9,
|
||||||
|
const double expFactor
|
||||||
|
) const {
|
||||||
|
double reverseRate = 0.0;
|
||||||
|
const double forwardRate = reaction.calculate_rate(T9);
|
||||||
|
|
||||||
|
if (reaction.reactants().size() == 2 && reaction.products().size() == 2) {
|
||||||
|
reverseRate = calculateReverseRateTwoBody(reaction, T9, forwardRate, expFactor);
|
||||||
|
} else {
|
||||||
|
LOG_WARNING(m_logger, "Reverse rate calculation for reactions with more than two reactants or products is not implemented.");
|
||||||
|
}
|
||||||
|
return reverseRate;
|
||||||
|
}
|
||||||
|
|
||||||
|
double GraphEngine::calculateReverseRateTwoBody(
|
||||||
|
const reaction::Reaction &reaction,
|
||||||
|
const double T9,
|
||||||
|
const double forwardRate,
|
||||||
|
const double expFactor
|
||||||
|
) const {
|
||||||
|
std::vector<double> reactantPartitionFunctions;
|
||||||
|
std::vector<double> productPartitionFunctions;
|
||||||
|
|
||||||
|
reactantPartitionFunctions.reserve(reaction.reactants().size());
|
||||||
|
productPartitionFunctions.reserve(reaction.products().size());
|
||||||
|
|
||||||
|
std::unordered_map<fourdst::atomic::Species, int> reactantMultiplicity;
|
||||||
|
std::unordered_map<fourdst::atomic::Species, int> productMultiplicity;
|
||||||
|
|
||||||
|
reactantMultiplicity.reserve(reaction.reactants().size());
|
||||||
|
productMultiplicity.reserve(reaction.products().size());
|
||||||
|
|
||||||
|
for (const auto& reactant : reaction.reactants()) {
|
||||||
|
reactantMultiplicity[reactant] += 1;
|
||||||
|
}
|
||||||
|
for (const auto& product : reaction.products()) {
|
||||||
|
productMultiplicity[product] += 1;
|
||||||
|
}
|
||||||
|
double reactantSymmetryFactor = 1.0;
|
||||||
|
double productSymmetryFactor = 1.0;
|
||||||
|
for (const auto& count : reactantMultiplicity | std::views::values) {
|
||||||
|
reactantSymmetryFactor *= std::tgamma(count + 1);
|
||||||
|
}
|
||||||
|
for (const auto& count : productMultiplicity | std::views::values) {
|
||||||
|
productSymmetryFactor *= std::tgamma(count + 1);
|
||||||
|
}
|
||||||
|
const double symmetryFactor = reactantSymmetryFactor / productSymmetryFactor;
|
||||||
|
|
||||||
|
// Accumulate mass terms
|
||||||
|
auto mass_op = [](double acc, const auto& species) { return acc * species.a(); };
|
||||||
|
const double massNumerator = std::accumulate(
|
||||||
|
reaction.reactants().begin(),
|
||||||
|
reaction.reactants().end(),
|
||||||
|
1.0,
|
||||||
|
mass_op
|
||||||
|
);
|
||||||
|
const double massDenominator = std::accumulate(
|
||||||
|
reaction.products().begin(),
|
||||||
|
reaction.products().end(),
|
||||||
|
1.0,
|
||||||
|
mass_op
|
||||||
|
);
|
||||||
|
|
||||||
|
// Accumulate partition functions
|
||||||
|
auto pf_op = [&](double acc, const auto& species) { return acc * m_partitionFunction->evaluate(species.z(), species.a(), T9); };
|
||||||
|
const double partitionFunctionNumerator = std::accumulate(
|
||||||
|
reaction.reactants().begin(),
|
||||||
|
reaction.reactants().end(),
|
||||||
|
1.0,
|
||||||
|
pf_op
|
||||||
|
);
|
||||||
|
const double partitionFunctionDenominator = std::accumulate(
|
||||||
|
reaction.products().begin(),
|
||||||
|
reaction.products().end(),
|
||||||
|
1.0,
|
||||||
|
pf_op
|
||||||
|
);
|
||||||
|
|
||||||
|
const double CT = std::pow(massNumerator/massDenominator, 1.5) * (partitionFunctionNumerator/partitionFunctionDenominator);
|
||||||
|
|
||||||
|
return forwardRate * symmetryFactor * CT * expFactor;
|
||||||
|
|
||||||
|
}
|
||||||
|
|
||||||
|
double GraphEngine::calculateReverseRateTwoBodyDerivative(
|
||||||
|
const reaction::Reaction &reaction,
|
||||||
|
const double T9,
|
||||||
|
const double reverseRate
|
||||||
|
) const {
|
||||||
|
const double d_log_kFwd = reaction.calculate_forward_rate_log_derivative(T9);
|
||||||
|
|
||||||
|
auto log_deriv_pf_op = [&](double acc, const auto& species) {
|
||||||
|
const double g = m_partitionFunction->evaluate(species.z(), species.a(), T9);
|
||||||
|
const double dg_dT = m_partitionFunction->evaluateDerivative(species.z(), species.a(), T9);
|
||||||
|
return (g == 0.0) ? acc : acc + (dg_dT / g);
|
||||||
|
};
|
||||||
|
|
||||||
|
const double reactant_log_derivative_sum = std::accumulate(
|
||||||
|
reaction.reactants().begin(),
|
||||||
|
reaction.reactants().end(),
|
||||||
|
0.0,
|
||||||
|
log_deriv_pf_op
|
||||||
|
);
|
||||||
|
|
||||||
|
const double product_log_derivative_sum = std::accumulate(
|
||||||
|
reaction.products().begin(),
|
||||||
|
reaction.products().end(),
|
||||||
|
0.0,
|
||||||
|
log_deriv_pf_op
|
||||||
|
);
|
||||||
|
|
||||||
|
const double d_log_C = reactant_log_derivative_sum - product_log_derivative_sum;
|
||||||
|
|
||||||
|
const double d_log_exp = reaction.qValue() / (m_constants.kB * T9 * T9);
|
||||||
|
|
||||||
|
const double log_total_derivative = d_log_kFwd + d_log_C + d_log_exp;
|
||||||
|
|
||||||
|
return reverseRate * log_total_derivative; // Return the derivative of the reverse rate with respect to T9
|
||||||
|
|
||||||
|
}
|
||||||
|
|
||||||
StepDerivatives<double> GraphEngine::calculateAllDerivativesUsingPrecomputation(
|
StepDerivatives<double> GraphEngine::calculateAllDerivativesUsingPrecomputation(
|
||||||
const std::vector<double> &Y_in,
|
const std::vector<double> &Y_in,
|
||||||
const std::vector<double> &bare_rates,
|
const std::vector<double> &bare_rates,
|
||||||
@@ -374,6 +509,10 @@ namespace gridfire {
|
|||||||
return m_usePrecomputation;
|
return m_usePrecomputation;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
const partition::PartitionFunction & GraphEngine::getPartitionFunction() const {
|
||||||
|
return *m_partitionFunction;
|
||||||
|
}
|
||||||
|
|
||||||
double GraphEngine::calculateMolarReactionFlow(
|
double GraphEngine::calculateMolarReactionFlow(
|
||||||
const reaction::Reaction &reaction,
|
const reaction::Reaction &reaction,
|
||||||
const std::vector<double> &Y,
|
const std::vector<double> &Y,
|
||||||
@@ -648,4 +787,21 @@ namespace gridfire {
|
|||||||
m_precomputedReactions.push_back(std::move(precomp));
|
m_precomputedReactions.push_back(std::move(precomp));
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
bool AtomicReverseRate::forward(
|
||||||
|
const size_t p,
|
||||||
|
const size_t q,
|
||||||
|
const CppAD::vector<bool> &vx,
|
||||||
|
CppAD::vector<bool> &vy,
|
||||||
|
const CppAD::vector<double> &tx,
|
||||||
|
CppAD::vector<double> &ty
|
||||||
|
) {
|
||||||
|
const double T9 = tx[0];
|
||||||
|
const double expFactor = std::exp(-m_reaction.qValue() / (m_kB * T9 * 1e9)); // Convert MeV to erg
|
||||||
|
if (p == 0) {
|
||||||
|
// --- Zeroth order forward sweep ---
|
||||||
|
const auto k_rev = m_engine.calculateReverseRate(m_reaction, T9, expFactor);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
}
|
}
|
||||||
|
|||||||
@@ -17,6 +17,13 @@ namespace gridfire::partition {
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
CompositePartitionFunction::CompositePartitionFunction(const CompositePartitionFunction &other) {
|
||||||
|
m_partitionFunctions.reserve(other.m_partitionFunctions.size());
|
||||||
|
for (const auto& pf : other.m_partitionFunctions) {
|
||||||
|
m_partitionFunctions.push_back(pf->clone());
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
double CompositePartitionFunction::evaluate(int z, int a, double T9) const {
|
double CompositePartitionFunction::evaluate(int z, int a, double T9) const {
|
||||||
LOG_TRACE_L1(m_logger, "Evaluating partition function for Z={} A={} T9={}", z, a, T9);
|
LOG_TRACE_L1(m_logger, "Evaluating partition function for Z={} A={} T9={}", z, a, T9);
|
||||||
for (const auto& partitionFunction : m_partitionFunctions) {
|
for (const auto& partitionFunction : m_partitionFunctions) {
|
||||||
|
|||||||
@@ -9,6 +9,7 @@
|
|||||||
#include <algorithm>
|
#include <algorithm>
|
||||||
#include <vector>
|
#include <vector>
|
||||||
#include <array>
|
#include <array>
|
||||||
|
#include <iostream>
|
||||||
|
|
||||||
namespace gridfire::partition {
|
namespace gridfire::partition {
|
||||||
static constexpr std::array<double, 24> RT_TEMPERATURE_GRID_T9 = {
|
static constexpr std::array<double, 24> RT_TEMPERATURE_GRID_T9 = {
|
||||||
@@ -25,8 +26,17 @@ namespace gridfire::partition {
|
|||||||
IsotopeData data;
|
IsotopeData data;
|
||||||
data.ground_state_spin = record.ground_state_spin;
|
data.ground_state_spin = record.ground_state_spin;
|
||||||
std::ranges::copy(record.normalized_g_values, data.normalized_g_values.begin());
|
std::ranges::copy(record.normalized_g_values, data.normalized_g_values.begin());
|
||||||
|
const int key = make_key(record.z, record.a);
|
||||||
|
LOG_TRACE_L3_LIMIT_EVERY_N(
|
||||||
|
100,
|
||||||
|
m_logger,
|
||||||
|
"(EVERY 100) Adding Rauscher-Thielemann partition data for Z={} A={} (key={})",
|
||||||
|
record.z,
|
||||||
|
record.a,
|
||||||
|
key
|
||||||
|
);
|
||||||
|
|
||||||
m_partitionData[make_key(record.z, record.a)] = std::move(data);
|
m_partitionData[key] = std::move(data);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|||||||
@@ -44,6 +44,27 @@ namespace gridfire::reaction {
|
|||||||
return calculate_rate<CppAD::AD<double>>(T9);
|
return calculate_rate<CppAD::AD<double>>(T9);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
double Reaction::calculate_forward_rate_log_derivative(const double T9) const {
|
||||||
|
constexpr double r_p13 = 1.0 / 3.0;
|
||||||
|
constexpr double r_p53 = 5.0 / 3.0;
|
||||||
|
constexpr double r_p23 = 2.0 / 3.0;
|
||||||
|
constexpr double r_p43 = 4.0 / 3.0;
|
||||||
|
|
||||||
|
const double T9_m1 = 1.0 / T9;
|
||||||
|
const double T9_m23 = std::pow(T9, -r_p23);
|
||||||
|
const double T9_m43 = std::pow(T9, -r_p43);
|
||||||
|
|
||||||
|
const double d_log_k_fwd_dT9 =
|
||||||
|
-m_rateCoefficients.a1 * T9_m1 * T9_m1
|
||||||
|
- r_p13 * m_rateCoefficients.a2 * T9_m43
|
||||||
|
+ r_p13 * m_rateCoefficients.a3 * T9_m23
|
||||||
|
+ m_rateCoefficients.a4
|
||||||
|
+ r_p53 * m_rateCoefficients.a5 * std::pow(T9, r_p23)
|
||||||
|
+ m_rateCoefficients.a6 * T9_m1;
|
||||||
|
|
||||||
|
return d_log_k_fwd_dT9; // Return the derivative of the log rate with respect to T9
|
||||||
|
}
|
||||||
|
|
||||||
bool Reaction::contains(const Species &species) const {
|
bool Reaction::contains(const Species &species) const {
|
||||||
return contains_reactant(species) || contains_product(species);
|
return contains_reactant(species) || contains_product(species);
|
||||||
}
|
}
|
||||||
@@ -194,6 +215,57 @@ namespace gridfire::reaction {
|
|||||||
return calculate_rate<double>(T9);
|
return calculate_rate<double>(T9);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
double LogicalReaction::calculate_forward_rate_log_derivative(const double T9) const {
|
||||||
|
constexpr double r_p13 = 1.0 / 3.0;
|
||||||
|
constexpr double r_p53 = 5.0 / 3.0;
|
||||||
|
constexpr double r_p23 = 2.0 / 3.0;
|
||||||
|
constexpr double r_p43 = 4.0 / 3.0;
|
||||||
|
|
||||||
|
double totalRate = 0.0;
|
||||||
|
double totalRateDerivative = 0.0;
|
||||||
|
|
||||||
|
|
||||||
|
const double T9_m1 = 1.0 / T9;
|
||||||
|
const double T913 = std::pow(T9, r_p13);
|
||||||
|
const double T953 = std::pow(T9, r_p53);
|
||||||
|
const double logT9 = std::log(T9);
|
||||||
|
|
||||||
|
const double T9_m1_sq = T9_m1 * T9_m1;
|
||||||
|
const double T9_m23 = std::pow(T9, -r_p23);
|
||||||
|
const double T9_m43 = std::pow(T9, -r_p43);
|
||||||
|
const double T9_p23 = std::pow(T9, r_p23);
|
||||||
|
|
||||||
|
|
||||||
|
for (const auto& coeffs : m_rates) {
|
||||||
|
const double exponent = coeffs.a0 +
|
||||||
|
coeffs.a1 * T9_m1 +
|
||||||
|
coeffs.a2 / T913 +
|
||||||
|
coeffs.a3 * T913 +
|
||||||
|
coeffs.a4 * T9 +
|
||||||
|
coeffs.a5 * T953 +
|
||||||
|
coeffs.a6 * logT9;
|
||||||
|
const double individualRate = std::exp(exponent);
|
||||||
|
|
||||||
|
const double d_exponent_T9 =
|
||||||
|
-coeffs.a1 * T9_m1_sq
|
||||||
|
- r_p13 * coeffs.a2 * T9_m43
|
||||||
|
+ r_p13 * coeffs.a3 * T9_m23
|
||||||
|
+ coeffs.a4
|
||||||
|
+ r_p53 * coeffs.a5 * T9_p23
|
||||||
|
+ coeffs.a6 * T9_m1;
|
||||||
|
|
||||||
|
const double individualRateDerivative = individualRate * d_exponent_T9;
|
||||||
|
totalRate += individualRate;
|
||||||
|
totalRateDerivative += individualRateDerivative;
|
||||||
|
}
|
||||||
|
|
||||||
|
if (totalRate == 0.0) {
|
||||||
|
return 0.0; // Avoid division by zero
|
||||||
|
}
|
||||||
|
|
||||||
|
return totalRateDerivative / totalRate;
|
||||||
|
}
|
||||||
|
|
||||||
CppAD::AD<double> LogicalReaction::calculate_rate(const CppAD::AD<double> T9) const {
|
CppAD::AD<double> LogicalReaction::calculate_rate(const CppAD::AD<double> T9) const {
|
||||||
return calculate_rate<CppAD::AD<double>>(T9);
|
return calculate_rate<CppAD::AD<double>>(T9);
|
||||||
}
|
}
|
||||||
|
|||||||
@@ -4,6 +4,7 @@
|
|||||||
#include "gridfire/engine/engine_graph.h"
|
#include "gridfire/engine/engine_graph.h"
|
||||||
#include "gridfire/engine/engine_approx8.h"
|
#include "gridfire/engine/engine_approx8.h"
|
||||||
#include "gridfire/engine/views/engine_adaptive.h"
|
#include "gridfire/engine/views/engine_adaptive.h"
|
||||||
|
#include "gridfire/partition/partition_types.h"
|
||||||
#include "gridfire/engine/views/engine_defined.h"
|
#include "gridfire/engine/views/engine_defined.h"
|
||||||
#include "gridfire/io/network_file.h"
|
#include "gridfire/io/network_file.h"
|
||||||
|
|
||||||
@@ -57,7 +58,7 @@ void quill_terminate_handler()
|
|||||||
int main() {
|
int main() {
|
||||||
g_previousHandler = std::set_terminate(quill_terminate_handler);
|
g_previousHandler = std::set_terminate(quill_terminate_handler);
|
||||||
quill::Logger* logger = fourdst::logging::LogManager::getInstance().getLogger("log");
|
quill::Logger* logger = fourdst::logging::LogManager::getInstance().getLogger("log");
|
||||||
logger->set_log_level(quill::LogLevel::TraceL3);
|
logger->set_log_level(quill::LogLevel::TraceL2);
|
||||||
LOG_DEBUG(logger, "Starting Adaptive Engine View Example...");
|
LOG_DEBUG(logger, "Starting Adaptive Engine View Example...");
|
||||||
|
|
||||||
using namespace gridfire;
|
using namespace gridfire;
|
||||||
@@ -93,13 +94,14 @@ int main() {
|
|||||||
BasePartitionType::RauscherThielemann,
|
BasePartitionType::RauscherThielemann,
|
||||||
BasePartitionType::GroundState
|
BasePartitionType::GroundState
|
||||||
});
|
});
|
||||||
std::cout << "Partition Function for Mg-24: " << partitionFunction.evaluate(12, 24, 1.5) << std::endl;
|
std::cout << "Partition Function for Mg-24: " << partitionFunction.evaluate(12, 24, 8) << std::endl;
|
||||||
std::cout << "Partition Function for F-23: " << partitionFunction.evaluate(9, 23, 1.5) << std::endl;
|
std::cout << "Partition Function for F-23: " << partitionFunction.evaluate(9, 23, 8) << std::endl;
|
||||||
std::cout << "Partition Function for O-13: " << partitionFunction.evaluate(8, 13, 1.5) << std::endl;
|
std::cout << "Partition Function for O-13: " << partitionFunction.evaluate(8, 13, 8) << std::endl;
|
||||||
|
|
||||||
netIn.dt0 = 1e-15;
|
netIn.dt0 = 1e-15;
|
||||||
|
|
||||||
// GraphEngine ReaclibEngine(composition);
|
GraphEngine ReaclibEngine(composition, partitionFunction);
|
||||||
|
std::cout << ReaclibEngine.getPartitionFunction().type() << std::endl;
|
||||||
// ReaclibEngine.setPrecomputation(true);
|
// ReaclibEngine.setPrecomputation(true);
|
||||||
// // AdaptiveEngineView adaptiveEngine(ReaclibEngine);
|
// // AdaptiveEngineView adaptiveEngine(ReaclibEngine);
|
||||||
// io::SimpleReactionListFileParser parser{};
|
// io::SimpleReactionListFileParser parser{};
|
||||||
|
|||||||
@@ -31,6 +31,7 @@ def generate_binary_file(input_filepath, output_filepath):
|
|||||||
spin = float(match.group('spin'))
|
spin = float(match.group('spin'))
|
||||||
coeffs_str = match.group('coeffs')
|
coeffs_str = match.group('coeffs')
|
||||||
g_values = [float(val) for val in coeffs_str.split()]
|
g_values = [float(val) for val in coeffs_str.split()]
|
||||||
|
print(f"Processing Z={z}, A={a}, Spin={spin}")
|
||||||
|
|
||||||
if len(g_values) != 24:
|
if len(g_values) != 24:
|
||||||
print(f"Warning: Found {len(g_values)} coefficients for Z={z}, A={a}. Expected 24. Skipping.")
|
print(f"Warning: Found {len(g_values)} coefficients for Z={z}, A={a}. Expected 24. Skipping.")
|
||||||
|
|||||||
Reference in New Issue
Block a user