feat(Jacobian): Jacobian generation is now stateless.

Previously Jacobians were stored by engines and accessed through engine
accessors (e.g getJacobianMatrixEntry); however, this resulted in
desynced jacobian states. We have changed to a pattern of Engine creates
a jacobian and returns it to the caller. The caller can then do what
they will with it. Because of this the getJacobianMatrixEntry method has
been removed.

BREAKING CHANGE:
    - There is no longer any getJacobianMatrixEntry method on
DynamicEngine classes
    - the generateJacobian method signature has changed to return a
NetworkJacobian object. Internally this uses an Eigen Sparse Matrix to
store its data.
This commit is contained in:
2025-11-14 10:51:40 -05:00
parent 1500f863b6
commit 9417b79a32
14 changed files with 352 additions and 440 deletions

View File

@@ -73,6 +73,7 @@ namespace gridfire {
/**
* @class GraphEngine
* @brief A reaction network engine that uses a graph-based representation.
@@ -236,7 +237,7 @@ namespace gridfire {
*
* @see getJacobianMatrixEntry()
*/
void generateJacobianMatrix(
[[nodiscard]] NetworkJacobian generateJacobianMatrix(
const fourdst::composition::CompositionAbstract &comp,
double T9,
double rho
@@ -254,7 +255,7 @@ namespace gridfire {
* @see getJacobianMatrixEntry()
* @see generateJacobianMatrix()
*/
void generateJacobianMatrix(
[[nodiscard]] NetworkJacobian generateJacobianMatrix(
const fourdst::composition::CompositionAbstract &comp,
double T9,
double rho,
@@ -276,7 +277,7 @@ namespace gridfire {
*
* @see getJacobianMatrixEntry()
*/
void generateJacobianMatrix(
[[nodiscard]] NetworkJacobian generateJacobianMatrix(
const fourdst::composition::CompositionAbstract &comp,
double T9,
double rho,
@@ -334,22 +335,6 @@ namespace gridfire {
*/
void setNetworkReactions(const reaction::ReactionSet& reactions) override;
/**
* @brief Gets an entry from the previously generated Jacobian matrix.
*
* @param rowSpecies Species corresponding to the row index (i).
* @param colSpecies Species corresponding to the column index (j).
* @return Value of the Jacobian matrix at (i, j).
*
* The Jacobian must have been generated by `generateJacobianMatrix()` before calling this.
*
* @see generateJacobianMatrix()
*/
[[nodiscard]] double getJacobianMatrixEntry(
const fourdst::atomic::Species& rowSpecies,
const fourdst::atomic::Species& colSpecies
) const override;
/**
* @brief Gets the net stoichiometry for a given reaction.
*
@@ -754,6 +739,8 @@ namespace gridfire {
*/
fourdst::composition::Composition collectComposition(const fourdst::composition::CompositionAbstract &comp, double T9, double rho) const override;
SpeciesStatus getSpeciesStatus(const fourdst::atomic::Species &species) const override;
private:
struct PrecomputedReaction {
@@ -864,9 +851,6 @@ namespace gridfire {
boost::numeric::ublas::compressed_matrix<int> m_stoichiometryMatrix; ///< Stoichiometry matrix (species x reactions).
mutable boost::numeric::ublas::compressed_matrix<double> m_jacobianMatrix; ///< Jacobian matrix (species x species).
mutable JacobianMatrixState m_jacobianMatrixState = JacobianMatrixState::UNINITIALIZED;
mutable CppAD::ADFun<double> m_rhsADFun; ///< CppAD function for the right-hand side of the ODE.
mutable CppAD::ADFun<double> m_epsADFun; ///< CppAD function for the energy generation rate.
mutable CppAD::sparse_jac_work m_jac_work; ///< Work object for sparse Jacobian calculations.
@@ -921,14 +905,6 @@ namespace gridfire {
*/
void populateSpeciesToIndexMap();
/**
* @brief Reserves space for the Jacobian matrix.
*
* This method reserves space for the Jacobian matrix, which is used
* to store the partial derivatives of the right-hand side of the ODE
* with respect to the species abundances.
*/
void reserveJacobianMatrix() const;
/**
* @brief Records the AD tape for the right-hand side of the ODE.
@@ -1254,7 +1230,7 @@ namespace gridfire {
std::unordered_map<fourdst::atomic::Species, int> reactant_counts;
reactant_counts.reserve(reaction.reactants().size());
for (const auto& reactant : reaction.reactants()) {
reactant_counts[reactant]++;
reactant_counts[reactant] = reaction.countReactantOccurrences(reactant);
}
const int totalReactants = static_cast<int>(reaction.reactants().size());
@@ -1285,7 +1261,6 @@ namespace gridfire {
// the tape more expensive to record, but it will also mean that we only need to record it once for
// the entire network.
const T densityTerm = CppAD::pow(rho, totalReactants > 1 ? static_cast<T>(totalReactants - 1) : zero); // Density raised to the power of (N-1) for N reactants
return molar_concentration_product * k_reaction * densityTerm;
}