2025-06-26 15:13:46 -04:00
|
|
|
|
#pragma once
|
|
|
|
|
|
|
2025-06-29 14:53:39 -04:00
|
|
|
|
#include "gridfire/reaction/reaction.h"
|
2025-11-24 09:07:49 -05:00
|
|
|
|
#include "gridfire/types/types.h"
|
2025-07-01 11:40:03 -04:00
|
|
|
|
#include "gridfire/screening/screening_abstract.h"
|
|
|
|
|
|
#include "gridfire/screening/screening_types.h"
|
2025-06-26 15:13:46 -04:00
|
|
|
|
|
2025-07-14 14:50:49 -04:00
|
|
|
|
#include "gridfire/engine/types/reporting.h"
|
2025-11-14 10:51:40 -05:00
|
|
|
|
#include "gridfire/engine/types/jacobian.h"
|
2025-07-14 14:50:49 -04:00
|
|
|
|
|
2026-04-13 07:17:14 -04:00
|
|
|
|
#include "gridfire/exceptions/error_engine.h"
|
|
|
|
|
|
|
2025-12-12 12:08:47 -05:00
|
|
|
|
#include "gridfire/engine/scratchpads/blob.h"
|
|
|
|
|
|
|
2025-11-10 10:40:03 -05:00
|
|
|
|
#include "fourdst/composition/composition_abstract.h"
|
|
|
|
|
|
|
2025-06-26 15:13:46 -04:00
|
|
|
|
#include <vector>
|
|
|
|
|
|
#include <unordered_map>
|
2025-07-14 14:50:49 -04:00
|
|
|
|
#include <utility>
|
2025-07-18 15:23:43 -04:00
|
|
|
|
#include <expected>
|
2025-06-26 15:13:46 -04:00
|
|
|
|
|
2025-06-29 14:53:39 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @file engine_abstract.h
|
|
|
|
|
|
* @brief Abstract interfaces for reaction network engines in GridFire.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This header defines the abstract base classes and concepts for implementing
|
|
|
|
|
|
* reaction network solvers in the GridFire framework. It provides the contract
|
|
|
|
|
|
* for calculating right-hand sides, energy generation, Jacobians, stoichiometry,
|
|
|
|
|
|
* and other core operations required for time integration of nuclear reaction networks.
|
|
|
|
|
|
*
|
|
|
|
|
|
* @author
|
|
|
|
|
|
* Emily M. Boudreaux
|
|
|
|
|
|
*/
|
|
|
|
|
|
|
2025-11-24 09:07:49 -05:00
|
|
|
|
namespace gridfire::engine {
|
2025-06-26 15:13:46 -04:00
|
|
|
|
|
|
|
|
|
|
|
2025-06-29 14:53:39 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Structure holding derivatives and energy generation for a network step.
|
|
|
|
|
|
*
|
|
|
|
|
|
* @tparam T Numeric type (double or CppAD::AD<double>).
|
|
|
|
|
|
*
|
|
|
|
|
|
* This struct is used to return both the time derivatives of all species abundances
|
|
|
|
|
|
* and the specific nuclear energy generation rate for a single network evaluation.
|
|
|
|
|
|
*
|
|
|
|
|
|
* Example usage:
|
|
|
|
|
|
* @code
|
|
|
|
|
|
* StepDerivatives<double> result = engine.calculateRHSAndEnergy(Y, T9, rho);
|
|
|
|
|
|
* for (double dydt_i : result.dydt) {
|
|
|
|
|
|
* // Use derivative
|
|
|
|
|
|
* }
|
|
|
|
|
|
* double energyRate = result.nuclearEnergyGenerationRate;
|
|
|
|
|
|
* @endcode
|
|
|
|
|
|
*/
|
2025-06-26 15:13:46 -04:00
|
|
|
|
template <IsArithmeticOrAD T>
|
|
|
|
|
|
struct StepDerivatives {
|
2025-11-10 10:40:03 -05:00
|
|
|
|
std::map<fourdst::atomic::Species, T> dydt{}; ///< Derivatives of abundances (dY/dt for each species).
|
2025-06-29 14:53:39 -04:00
|
|
|
|
T nuclearEnergyGenerationRate = T(0.0); ///< Specific energy generation rate (e.g., erg/g/s).
|
2025-12-02 13:09:19 -05:00
|
|
|
|
std::optional<std::map<fourdst::atomic::Species, std::unordered_map<std::string, T>>> reactionContributions = std::nullopt;
|
2025-11-27 14:34:20 -05:00
|
|
|
|
T neutrinoEnergyLossRate = T(0.0); // (erg/g/s)
|
|
|
|
|
|
T totalNeutrinoFlux = T(0.0); // (neutrinos/g/s)
|
|
|
|
|
|
|
2025-10-07 15:16:03 -04:00
|
|
|
|
|
|
|
|
|
|
StepDerivatives() : dydt(), nuclearEnergyGenerationRate(T(0.0)) {}
|
2025-06-26 15:13:46 -04:00
|
|
|
|
};
|
|
|
|
|
|
|
2025-11-24 09:07:49 -05:00
|
|
|
|
using SparsityPattern = std::vector<std::pair<size_t, size_t>>; ///< Type alias for sparsity pattern representation.
|
2025-07-14 14:50:49 -04:00
|
|
|
|
|
2025-11-24 09:07:49 -05:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Structure holding derivatives of energy generation rate with respect to T and rho.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This struct encapsulates the partial derivatives of the specific nuclear energy
|
|
|
|
|
|
* generation rate with respect to temperature and density.
|
|
|
|
|
|
*/
|
2025-09-19 15:14:46 -04:00
|
|
|
|
struct EnergyDerivatives {
|
|
|
|
|
|
double dEps_dT = 0.0; ///< Partial derivative of energy generation rate with respect to temperature.
|
|
|
|
|
|
double dEps_dRho = 0.0;///< Partial derivative of energy generation rate with respect to density.
|
|
|
|
|
|
|
|
|
|
|
|
friend std::ostream& operator<<(std::ostream& os, const EnergyDerivatives& ed) {
|
|
|
|
|
|
os << "<dε/dT: " << ed.dEps_dT << ", dε/dρ: " << ed.dEps_dRho << ">";
|
|
|
|
|
|
return os;
|
|
|
|
|
|
}
|
|
|
|
|
|
};
|
|
|
|
|
|
|
2025-11-24 09:07:49 -05:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Enumeration of possible engine statuses.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This enum defines the various states an engine can be in after performing
|
|
|
|
|
|
* calculations, such as being up-to-date (OKAY), needing an update (STALE),
|
|
|
|
|
|
* or encountering an error (ERROR).
|
|
|
|
|
|
*/
|
|
|
|
|
|
enum class EngineStatus {
|
|
|
|
|
|
OKAY,
|
|
|
|
|
|
STALE,
|
|
|
|
|
|
ERROR,
|
|
|
|
|
|
COUNT
|
|
|
|
|
|
};
|
|
|
|
|
|
|
|
|
|
|
|
/**
|
|
|
|
|
|
* @brief Convert EngineStatus enum to string representation.
|
|
|
|
|
|
*
|
|
|
|
|
|
* @param status The EngineStatus value to convert.
|
|
|
|
|
|
* @return A string_view representing the name of the EngineStatus.
|
|
|
|
|
|
*/
|
|
|
|
|
|
constexpr std::string_view EngineStatus_to_string(const EngineStatus status) {
|
|
|
|
|
|
constexpr std::array<std::string_view, static_cast<size_t>(EngineStatus::COUNT)> names = {
|
|
|
|
|
|
"OKAY",
|
|
|
|
|
|
"STALE",
|
|
|
|
|
|
"ERROR"
|
|
|
|
|
|
};
|
|
|
|
|
|
return names[static_cast<size_t>(status)];
|
|
|
|
|
|
}
|
|
|
|
|
|
|
2025-06-29 14:53:39 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Abstract base class for a reaction network engine.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This class defines the minimal interface for a reaction network engine,
|
|
|
|
|
|
* which is responsible for evaluating the right-hand side (dY/dt) and
|
|
|
|
|
|
* energy generation for a given set of abundances, temperature, and density.
|
|
|
|
|
|
*
|
|
|
|
|
|
* Intended usage: Derive from this class to implement a concrete engine
|
|
|
|
|
|
* for a specific network or integration method.
|
|
|
|
|
|
*
|
|
|
|
|
|
* Example:
|
|
|
|
|
|
* @code
|
|
|
|
|
|
* class MyEngine : public gridfire::Engine {
|
|
|
|
|
|
* // Implement required methods...
|
|
|
|
|
|
* };
|
|
|
|
|
|
* @endcode
|
|
|
|
|
|
*/
|
2025-06-26 15:13:46 -04:00
|
|
|
|
class Engine {
|
|
|
|
|
|
public:
|
2025-06-29 14:53:39 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Virtual destructor.
|
|
|
|
|
|
*/
|
2025-06-26 15:13:46 -04:00
|
|
|
|
virtual ~Engine() = default;
|
2025-06-29 14:53:39 -04:00
|
|
|
|
|
|
|
|
|
|
/**
|
|
|
|
|
|
* @brief Get the list of species in the network.
|
|
|
|
|
|
* @return Vector of Species objects representing all network species.
|
|
|
|
|
|
*/
|
2025-12-12 12:08:47 -05:00
|
|
|
|
[[nodiscard]] virtual const std::vector<fourdst::atomic::Species>& getNetworkSpecies(
|
|
|
|
|
|
scratch::StateBlob& ctx
|
|
|
|
|
|
) const = 0;
|
2025-06-29 14:53:39 -04:00
|
|
|
|
|
|
|
|
|
|
/**
|
|
|
|
|
|
* @brief Calculate the right-hand side (dY/dt) and energy generation.
|
|
|
|
|
|
*
|
2025-10-07 15:16:03 -04:00
|
|
|
|
* @param comp Composition object containing current abundances.
|
2025-06-29 14:53:39 -04:00
|
|
|
|
* @param T9 Temperature in units of 10^9 K.
|
|
|
|
|
|
* @param rho Density in g/cm^3.
|
2025-12-07 12:34:12 -05:00
|
|
|
|
* @param trust If true, indicates that the engine should trust the passed composition has already been collected.
|
2025-10-30 15:05:08 -04:00
|
|
|
|
* @return expected<StepDerivatives<double>> containing either dY/dt and energy generation rate or a stale engine
|
|
|
|
|
|
* error indicating that the engine must be updated
|
2025-06-29 14:53:39 -04:00
|
|
|
|
*
|
|
|
|
|
|
* This function must be implemented by derived classes to compute the
|
|
|
|
|
|
* time derivatives of all species and the specific nuclear energy generation
|
|
|
|
|
|
* rate for the current state.
|
|
|
|
|
|
*/
|
2025-11-24 09:07:49 -05:00
|
|
|
|
[[nodiscard]] virtual std::expected<StepDerivatives<double>, EngineStatus> calculateRHSAndEnergy(
|
2025-12-12 12:08:47 -05:00
|
|
|
|
scratch::StateBlob&,
|
2025-11-10 10:40:03 -05:00
|
|
|
|
const fourdst::composition::CompositionAbstract &comp,
|
2025-06-26 15:13:46 -04:00
|
|
|
|
double T9,
|
2025-12-07 12:34:12 -05:00
|
|
|
|
double rho,
|
|
|
|
|
|
bool trust
|
2025-06-26 15:13:46 -04:00
|
|
|
|
) const = 0;
|
|
|
|
|
|
};
|
|
|
|
|
|
|
2025-06-29 14:53:39 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Abstract class for engines supporting Jacobian and stoichiometry operations.
|
|
|
|
|
|
*
|
|
|
|
|
|
* Extends Engine with additional methods for:
|
|
|
|
|
|
* - Generating and accessing the Jacobian matrix (for implicit solvers).
|
|
|
|
|
|
* - Generating and accessing the stoichiometry matrix.
|
|
|
|
|
|
* - Calculating molar reaction flows for individual reactions.
|
|
|
|
|
|
* - Accessing the set of logical reactions in the network.
|
|
|
|
|
|
* - Computing timescales for each species.
|
|
|
|
|
|
*
|
|
|
|
|
|
* Intended usage: Derive from this class to implement engines that support
|
|
|
|
|
|
* advanced solver features such as implicit integration, sensitivity analysis,
|
2025-11-24 09:07:49 -05:00
|
|
|
|
* QSE (Quasi-Steady-State Equilibrium) handling, and more. Generally this will be the main engine type
|
2025-06-29 14:53:39 -04:00
|
|
|
|
*/
|
2025-06-26 15:13:46 -04:00
|
|
|
|
class DynamicEngine : public Engine {
|
|
|
|
|
|
public:
|
2025-06-29 14:53:39 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Generate the Jacobian matrix for the current state.
|
|
|
|
|
|
*
|
2026-04-13 07:17:14 -04:00
|
|
|
|
* @param ctx The scratchpad context for the current state.
|
2025-10-07 15:16:03 -04:00
|
|
|
|
* @param comp Composition object containing current abundances.
|
2025-06-29 14:53:39 -04:00
|
|
|
|
* @param T9 Temperature in units of 10^9 K.
|
|
|
|
|
|
* @param rho Density in g/cm^3.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method must compute and store the Jacobian matrix (∂(dY/dt)_i/∂Y_j)
|
|
|
|
|
|
* for the current state. The matrix can then be accessed via getJacobianMatrixEntry().
|
|
|
|
|
|
*/
|
2025-11-14 10:51:40 -05:00
|
|
|
|
[[nodiscard]] virtual NetworkJacobian generateJacobianMatrix(
|
2025-12-12 12:08:47 -05:00
|
|
|
|
scratch::StateBlob& ctx,
|
2025-11-10 10:40:03 -05:00
|
|
|
|
const fourdst::composition::CompositionAbstract &comp,
|
2025-07-14 14:50:49 -04:00
|
|
|
|
double T9,
|
|
|
|
|
|
double rho
|
2025-07-18 15:23:43 -04:00
|
|
|
|
) const = 0;
|
2025-06-29 14:53:39 -04:00
|
|
|
|
|
2025-11-24 09:07:49 -05:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Generate the Jacobian matrix for the current state using a subset of active species.
|
|
|
|
|
|
*
|
2026-04-13 07:17:14 -04:00
|
|
|
|
* @param ctx The scratchpad context for the current state.
|
2025-11-24 09:07:49 -05:00
|
|
|
|
* @param comp Composition object containing current abundances.
|
|
|
|
|
|
* @param T9 Temperature in units of 10^9 K.
|
|
|
|
|
|
* @param rho Density in g/cm^3.
|
|
|
|
|
|
* @param activeSpecies The set of species to include in the Jacobian calculation.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method must compute and store the Jacobian matrix (∂(dY/dt)_i/∂Y_j)
|
|
|
|
|
|
* for the current state, considering only the specified subset of active species.
|
|
|
|
|
|
* The matrix can then be accessed via getJacobianMatrixEntry().
|
|
|
|
|
|
*/
|
2025-11-14 10:51:40 -05:00
|
|
|
|
[[nodiscard]] virtual NetworkJacobian generateJacobianMatrix(
|
2025-12-12 12:08:47 -05:00
|
|
|
|
scratch::StateBlob& ctx,
|
2025-11-10 10:40:03 -05:00
|
|
|
|
const fourdst::composition::CompositionAbstract &comp,
|
2025-10-24 11:17:22 -04:00
|
|
|
|
double T9,
|
|
|
|
|
|
double rho,
|
|
|
|
|
|
const std::vector<fourdst::atomic::Species>& activeSpecies
|
|
|
|
|
|
) const = 0;
|
|
|
|
|
|
|
2025-11-24 09:07:49 -05:00
|
|
|
|
|
|
|
|
|
|
/**
|
|
|
|
|
|
* @brief Generate the Jacobian matrix for the current state with a specified sparsity pattern.
|
|
|
|
|
|
*
|
2026-04-13 07:17:14 -04:00
|
|
|
|
* @param ctx Get the scratchpad context for the current state.
|
2025-11-24 09:07:49 -05:00
|
|
|
|
* @param comp Composition object containing current abundances.
|
|
|
|
|
|
* @param T9 Temperature in units of 10^9 K.
|
|
|
|
|
|
* @param rho Density in g/cm^3.
|
|
|
|
|
|
* @param sparsityPattern The sparsity pattern to use for the Jacobian matrix.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method must compute and store the Jacobian matrix (∂(dY/dt)_i/∂Y_j)
|
|
|
|
|
|
* for the current state using automatic differentiation, taking into
|
|
|
|
|
|
* account the provided sparsity pattern. The matrix can then be accessed
|
|
|
|
|
|
* via `getJacobianMatrixEntry()`.
|
|
|
|
|
|
*
|
|
|
|
|
|
* @see getJacobianMatrixEntry()
|
|
|
|
|
|
*/
|
2025-11-14 10:51:40 -05:00
|
|
|
|
[[nodiscard]] virtual NetworkJacobian generateJacobianMatrix(
|
2025-12-12 12:08:47 -05:00
|
|
|
|
scratch::StateBlob& ctx,
|
2025-11-10 10:40:03 -05:00
|
|
|
|
const fourdst::composition::CompositionAbstract &comp,
|
2025-07-14 14:50:49 -04:00
|
|
|
|
double T9,
|
|
|
|
|
|
double rho,
|
|
|
|
|
|
const SparsityPattern& sparsityPattern
|
2025-10-24 11:17:22 -04:00
|
|
|
|
) const = 0;
|
2025-07-14 14:50:49 -04:00
|
|
|
|
|
|
|
|
|
|
|
2025-06-29 14:53:39 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Calculate the molar reaction flow for a given reaction.
|
|
|
|
|
|
*
|
2026-04-13 07:17:14 -04:00
|
|
|
|
* @param ctx The scratchpad context for the current state.
|
2025-06-29 14:53:39 -04:00
|
|
|
|
* @param reaction The reaction for which to calculate the flow.
|
2025-10-07 15:16:03 -04:00
|
|
|
|
* @param comp Composition object containing current abundances.
|
2025-06-29 14:53:39 -04:00
|
|
|
|
* @param T9 Temperature in units of 10^9 K.
|
|
|
|
|
|
* @param rho Density in g/cm^3.
|
|
|
|
|
|
* @return Molar flow rate for the reaction (e.g., mol/g/s).
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method computes the net rate at which the given reaction proceeds
|
|
|
|
|
|
* under the current state.
|
|
|
|
|
|
*/
|
2025-07-01 11:40:03 -04:00
|
|
|
|
[[nodiscard]] virtual double calculateMolarReactionFlow(
|
2025-12-12 12:08:47 -05:00
|
|
|
|
scratch::StateBlob& ctx,
|
2025-06-26 15:13:46 -04:00
|
|
|
|
const reaction::Reaction& reaction,
|
2025-11-10 10:40:03 -05:00
|
|
|
|
const fourdst::composition::CompositionAbstract &comp,
|
2025-06-26 15:13:46 -04:00
|
|
|
|
double T9,
|
|
|
|
|
|
double rho
|
|
|
|
|
|
) const = 0;
|
|
|
|
|
|
|
2025-09-19 15:14:46 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Calculate the derivatives of the energy generation rate with respect to T and rho.
|
|
|
|
|
|
*
|
2025-10-07 15:16:03 -04:00
|
|
|
|
* @param comp Composition object containing current abundances.
|
2025-09-19 15:14:46 -04:00
|
|
|
|
* @param T9 Temperature in units of 10^9 K.
|
|
|
|
|
|
* @param rho Density in g/cm^3.
|
|
|
|
|
|
* @return EnergyDerivatives containing dEps/dT and dEps/dRho.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method computes the partial derivatives of the specific nuclear energy
|
|
|
|
|
|
* generation rate with respect to temperature and density for the current state.
|
|
|
|
|
|
*/
|
|
|
|
|
|
[[nodiscard]] virtual EnergyDerivatives calculateEpsDerivatives(
|
2025-12-12 12:08:47 -05:00
|
|
|
|
scratch::StateBlob& ctx,
|
2025-11-10 10:40:03 -05:00
|
|
|
|
const fourdst::composition::CompositionAbstract &comp,
|
2025-10-07 15:16:03 -04:00
|
|
|
|
double T9,
|
|
|
|
|
|
double rho
|
2025-09-19 15:14:46 -04:00
|
|
|
|
) const = 0;
|
|
|
|
|
|
|
2025-06-29 14:53:39 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Get the set of logical reactions in the network.
|
|
|
|
|
|
*
|
|
|
|
|
|
* @return Reference to the LogicalReactionSet containing all reactions.
|
|
|
|
|
|
*/
|
2025-12-12 12:08:47 -05:00
|
|
|
|
[[nodiscard]] virtual const reaction::ReactionSet& getNetworkReactions(
|
|
|
|
|
|
scratch::StateBlob& ctx
|
|
|
|
|
|
) const = 0;
|
2025-06-29 14:53:39 -04:00
|
|
|
|
|
2026-04-13 07:17:14 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Get the set of inactive reactions in the network.
|
|
|
|
|
|
*
|
|
|
|
|
|
* @return ReactionSet containing all inactive reactions.
|
|
|
|
|
|
*
|
|
|
|
|
|
* By default, this method returns an empty set. Derived classes can override
|
|
|
|
|
|
* this method to provide the actual set of inactive reactions based on their
|
|
|
|
|
|
* internal logic (e.g., reaction flow culling, QSE partitioning).
|
|
|
|
|
|
*/
|
|
|
|
|
|
[[nodiscard]] virtual reaction::ReactionSet getInactiveNetworkReactions(
|
|
|
|
|
|
scratch::StateBlob &ctx
|
|
|
|
|
|
) const {
|
|
|
|
|
|
return reaction::ReactionSet{};
|
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
[[nodiscard]] virtual double getInactiveReactionMolarReactionFlow(
|
|
|
|
|
|
scratch::StateBlob& ctx,
|
|
|
|
|
|
const reaction::Reaction &reaction,
|
|
|
|
|
|
const fourdst::composition::CompositionAbstract &comp,
|
|
|
|
|
|
const double T9,
|
|
|
|
|
|
const double rho
|
|
|
|
|
|
) const {
|
|
|
|
|
|
std::string warning_msg = std::format(
|
|
|
|
|
|
"[GridFire Warning ({}, {}, {})]: Engine of type '{}' does not implement getInactiveReactionMolarReactionFlow. Returning 0.0 flow for reaction '{}'.",
|
|
|
|
|
|
__FILE__,
|
|
|
|
|
|
__LINE__,
|
|
|
|
|
|
__FUNCTION__,
|
|
|
|
|
|
typeid(*this).name(),
|
|
|
|
|
|
reaction.id()
|
|
|
|
|
|
);
|
|
|
|
|
|
return 0.0;
|
|
|
|
|
|
}
|
|
|
|
|
|
|
2025-07-22 12:48:24 -04:00
|
|
|
|
|
2025-06-29 14:53:39 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Compute timescales for all species in the network.
|
|
|
|
|
|
*
|
2025-10-07 15:16:03 -04:00
|
|
|
|
* @param comp Composition object containing current abundances.
|
2025-06-29 14:53:39 -04:00
|
|
|
|
* @param T9 Temperature in units of 10^9 K.
|
|
|
|
|
|
* @param rho Density in g/cm^3.
|
|
|
|
|
|
* @return Map from Species to their characteristic timescales (s).
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method estimates the timescale for abundance change of each species,
|
|
|
|
|
|
* which can be used for timestep control, diagnostics, and reaction network culling.
|
|
|
|
|
|
*/
|
2025-11-24 09:07:49 -05:00
|
|
|
|
[[nodiscard]] virtual std::expected<std::unordered_map<fourdst::atomic::Species, double>, EngineStatus> getSpeciesTimescales(
|
2025-12-12 12:08:47 -05:00
|
|
|
|
scratch::StateBlob& ctx,
|
2025-11-10 10:40:03 -05:00
|
|
|
|
const fourdst::composition::CompositionAbstract &comp,
|
2025-07-22 12:48:24 -04:00
|
|
|
|
double T9,
|
|
|
|
|
|
double rho
|
|
|
|
|
|
) const = 0;
|
|
|
|
|
|
|
2025-11-24 09:07:49 -05:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Compute destruction timescales for all species in the network.
|
|
|
|
|
|
*
|
2026-04-13 07:17:14 -04:00
|
|
|
|
* @param ctx The scratchpad context for the current state.
|
2025-11-24 09:07:49 -05:00
|
|
|
|
* @param comp Composition object containing current abundances.
|
|
|
|
|
|
* @param T9 Temperature in units of 10^9 K.
|
|
|
|
|
|
* @param rho Density in g/cm^3.
|
|
|
|
|
|
* @return Map from Species to their destruction timescales (s).
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method estimates the destruction timescale for each species,
|
|
|
|
|
|
* which can be useful for understanding reaction flows and equilibrium states.
|
|
|
|
|
|
*/
|
|
|
|
|
|
[[nodiscard]] virtual std::expected<std::unordered_map<fourdst::atomic::Species, double>, EngineStatus> getSpeciesDestructionTimescales(
|
2025-12-12 12:08:47 -05:00
|
|
|
|
scratch::StateBlob& ctx,
|
2025-11-10 10:40:03 -05:00
|
|
|
|
const fourdst::composition::CompositionAbstract &comp,
|
2025-06-26 15:13:46 -04:00
|
|
|
|
double T9,
|
|
|
|
|
|
double rho
|
|
|
|
|
|
) const = 0;
|
2025-07-01 11:40:03 -04:00
|
|
|
|
|
2025-07-01 15:06:22 -04:00
|
|
|
|
/**
|
2025-12-12 12:08:47 -05:00
|
|
|
|
* @brief Update the thread local scratch pad state of a network.
|
2025-07-01 15:06:22 -04:00
|
|
|
|
*
|
2026-04-13 07:17:14 -04:00
|
|
|
|
* @param ctx The scratchpad context for the current state.
|
2025-07-01 15:06:22 -04:00
|
|
|
|
* @param netIn A struct containing the current network input, such as
|
|
|
|
|
|
* temperature, density, and composition.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method is intended to be implemented by derived classes to update
|
|
|
|
|
|
* their internal state based on the provided network conditions. For example,
|
|
|
|
|
|
* an adaptive engine might use this to re-evaluate which reactions and species
|
|
|
|
|
|
* are active. For other engines that do not support manually updating, this
|
|
|
|
|
|
* method might do nothing.
|
|
|
|
|
|
*
|
|
|
|
|
|
* @par Usage Example:
|
|
|
|
|
|
* @code
|
|
|
|
|
|
* NetIn input = { ... };
|
|
|
|
|
|
* myEngine.update(input);
|
|
|
|
|
|
* @endcode
|
|
|
|
|
|
*
|
|
|
|
|
|
* @post The internal state of the engine is updated to reflect the new conditions.
|
|
|
|
|
|
*/
|
2025-12-12 12:08:47 -05:00
|
|
|
|
[[nodiscard]] virtual fourdst::composition::Composition project(
|
|
|
|
|
|
scratch::StateBlob& ctx,
|
|
|
|
|
|
const NetIn &netIn
|
|
|
|
|
|
) const = 0;
|
2025-07-01 11:40:03 -04:00
|
|
|
|
|
2025-07-01 15:06:22 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Get the current electron screening model.
|
|
|
|
|
|
*
|
2026-04-13 07:17:14 -04:00
|
|
|
|
* @param ctx The scratchpad context for the current state.
|
|
|
|
|
|
*
|
2025-07-01 15:06:22 -04:00
|
|
|
|
* @return The currently active screening model type.
|
|
|
|
|
|
*
|
|
|
|
|
|
* @par Usage Example:
|
|
|
|
|
|
* @code
|
|
|
|
|
|
* screening::ScreeningType currentModel = myEngine.getScreeningModel();
|
|
|
|
|
|
* @endcode
|
|
|
|
|
|
*/
|
2025-12-12 12:08:47 -05:00
|
|
|
|
[[nodiscard]] virtual screening::ScreeningType getScreeningModel(
|
|
|
|
|
|
scratch::StateBlob& ctx
|
|
|
|
|
|
) const = 0;
|
2025-07-10 09:36:05 -04:00
|
|
|
|
|
2025-07-31 15:04:57 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Get the index of a species in the network.
|
|
|
|
|
|
*
|
2026-04-13 07:17:14 -04:00
|
|
|
|
* @param ctx The scratchpad context for the current state.
|
2025-07-31 15:04:57 -04:00
|
|
|
|
* @param species The species to look up.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method allows querying the index of a specific species in the
|
|
|
|
|
|
* engine's internal representation. It is useful for accessing species
|
|
|
|
|
|
* data efficiently.
|
|
|
|
|
|
*/
|
2025-12-12 12:08:47 -05:00
|
|
|
|
[[nodiscard]] virtual size_t getSpeciesIndex(
|
|
|
|
|
|
scratch::StateBlob& ctx,
|
|
|
|
|
|
const fourdst::atomic::Species &species
|
|
|
|
|
|
) const = 0;
|
2025-07-14 14:50:49 -04:00
|
|
|
|
|
2025-07-31 15:04:57 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Prime the engine with initial conditions.
|
|
|
|
|
|
*
|
2026-04-13 07:17:14 -04:00
|
|
|
|
* @param ctx The scratchpad context for the current state.
|
2025-07-31 15:04:57 -04:00
|
|
|
|
* @param netIn The input conditions for the network.
|
|
|
|
|
|
* @return PrimingReport containing information about the priming process.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method is used to prepare the engine for calculations by setting up
|
|
|
|
|
|
* initial conditions, reactions, and species. It may involve compiling reaction
|
|
|
|
|
|
* rates, initializing internal data structures, and performing any necessary
|
|
|
|
|
|
* pre-computation.
|
|
|
|
|
|
*/
|
2025-12-12 12:08:47 -05:00
|
|
|
|
[[nodiscard]] virtual PrimingReport primeEngine(
|
|
|
|
|
|
scratch::StateBlob& ctx,
|
|
|
|
|
|
const NetIn &netIn
|
|
|
|
|
|
) const = 0;
|
2025-07-18 15:23:43 -04:00
|
|
|
|
|
2025-10-30 15:05:08 -04:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Recursively collect composition from current engine and any sub engines if they exist.
|
|
|
|
|
|
* @details If species i is defined in comp and in any sub engine or self composition then the molar abundance of
|
|
|
|
|
|
* species i in the returned composition will be that defined in comp. If there are species defined in sub engine
|
|
|
|
|
|
* compositions which are not defined in comp then their molar abundances will be based on the reported values
|
|
|
|
|
|
* from each sub engine.
|
|
|
|
|
|
* @note It is up to each engine to decide how to handle filling in the return composition.
|
|
|
|
|
|
* @note These methods return an unfinalized composition which must then be finalized by the caller
|
2026-04-13 07:17:14 -04:00
|
|
|
|
* @param ctx The scratchpad context for the current state.
|
2025-10-30 15:05:08 -04:00
|
|
|
|
* @param comp Input composition to "normalize".
|
2025-11-12 16:54:12 -05:00
|
|
|
|
* @param T9
|
|
|
|
|
|
* @param rho
|
2025-10-30 15:05:08 -04:00
|
|
|
|
* @return An updated composition which is a superset of comp. This may contain species which were culled, for
|
|
|
|
|
|
* example, by either QSE partitioning or reaction flow rate culling
|
|
|
|
|
|
*/
|
2025-12-12 12:08:47 -05:00
|
|
|
|
[[nodiscard]] virtual fourdst::composition::Composition collectComposition(
|
|
|
|
|
|
scratch::StateBlob& ctx,
|
2025-11-12 16:54:12 -05:00
|
|
|
|
const fourdst::composition::CompositionAbstract &comp,
|
|
|
|
|
|
double T9,
|
|
|
|
|
|
double rho
|
2025-10-30 15:05:08 -04:00
|
|
|
|
) const = 0;
|
|
|
|
|
|
|
2025-11-24 09:07:49 -05:00
|
|
|
|
/**
|
|
|
|
|
|
* @brief Get the status of a species in the network.
|
|
|
|
|
|
*
|
|
|
|
|
|
* @param species The species to check.
|
|
|
|
|
|
* @return SpeciesStatus indicating whether the species is active, inactive, or culled.
|
|
|
|
|
|
*
|
|
|
|
|
|
* This method allows querying the current status of a specific species
|
|
|
|
|
|
* within the engine's network.
|
|
|
|
|
|
*/
|
2025-12-12 12:08:47 -05:00
|
|
|
|
[[nodiscard]] virtual SpeciesStatus getSpeciesStatus(
|
|
|
|
|
|
scratch::StateBlob& ctx,
|
|
|
|
|
|
const fourdst::atomic::Species& species
|
|
|
|
|
|
) const = 0;
|
2025-11-14 10:51:40 -05:00
|
|
|
|
|
2025-12-12 12:08:47 -05:00
|
|
|
|
[[nodiscard]] virtual std::optional<StepDerivatives<double>> getMostRecentRHSCalculation(
|
|
|
|
|
|
scratch::StateBlob& ctx
|
|
|
|
|
|
) const = 0;
|
2025-12-10 12:50:35 -05:00
|
|
|
|
|
2026-04-13 07:17:14 -04:00
|
|
|
|
[[nodiscard]] virtual std::unique_ptr<scratch::StateBlob> constructStateBlob(const scratch::StateBlob *blob) const = 0;
|
|
|
|
|
|
|
2025-06-26 15:13:46 -04:00
|
|
|
|
};
|
|
|
|
|
|
}
|