|
GridFire v0.7.0-alpha
General Purpose Nuclear Network
|
Functor for solving QSE abundances using Eigen's nonlinear optimization. More...
Public Types | |
| enum | { InputsAtCompileTime = Eigen::Dynamic , ValuesAtCompileTime = Eigen::Dynamic } |
| using | InputType = Eigen::Matrix<double, Eigen::Dynamic, 1> |
| using | OutputType = Eigen::Matrix<double, Eigen::Dynamic, 1> |
| using | JacobianType = Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> |
Public Member Functions | |
| EigenFunctor (MultiscalePartitioningEngineView &view, const std::set< fourdst::atomic::Species > &qse_solve_species, const fourdst::composition::Composition &initial_comp, const double T9, const double rho, const Eigen::VectorXd &Y_scale, const std::unordered_map< fourdst::atomic::Species, size_t > &qse_solve_species_index_map) | |
| Constructs an EigenFunctor. | |
| size_t | values () const |
| Gets the number of output values from the functor (size of the residual vector). | |
| size_t | inputs () const |
| Gets the number of input values to the functor (size of the variable vector). | |
| int | operator() (const InputType &v_qse, OutputType &f_qse) const |
Evaluates the functor's residual vector f_qse = dY_alg/dt. | |
| int | df (const InputType &v_qse, JacobianType &J_qse) const |
Evaluates the Jacobian of the functor, J_qse = d(f_qse)/d(v_qse). | |
Public Attributes | |
| MultiscalePartitioningEngineView * | m_view |
| Pointer to the MultiscalePartitioningEngineView instance. | |
| const std::set< fourdst::atomic::Species > & | m_qse_solve_species |
| The set of species to solve for in the QSE group. | |
| const fourdst::composition::Composition & | m_initial_comp |
| Initial abundances of all species in the full network. | |
| const double | m_T9 |
| Temperature in units of 10^9 K. | |
| const double | m_rho |
| Density in g/cm^3. | |
| const Eigen::VectorXd & | m_Y_scale |
| Scaling factors for the species abundances, used to improve solver stability. | |
| const std::unordered_map< fourdst::atomic::Species, size_t > | m_qse_solve_species_index_map |
| Mapping from species to their indices in the QSE solve vector. | |
Functor for solving QSE abundances using Eigen's nonlinear optimization.
operator()) and its Jacobian (df) to Eigen's Levenberg-Marquardt solver. The goal is to find the abundances of algebraic species that make their time derivatives (dY/dt) equal to zero.@how
operator(): Takes a vector v_qse (scaled abundances of algebraic species) as input. It constructs a full trial abundance vector y_trial, calls the base engine's calculateRHSAndEnergy, and returns the dY/dt values for the algebraic species. The solver attempts to drive this return vector to zero.df: Computes the Jacobian of the objective function. It calls the base engine's generateJacobianMatrix and extracts the sub-matrix corresponding to the algebraic species. It applies the chain rule to account for the asinh scaling used on the abundances.The abundances are scaled using asinh to handle the large dynamic range and ensure positivity.
| using gridfire::MultiscalePartitioningEngineView::EigenFunctor::InputType = Eigen::Matrix<double, Eigen::Dynamic, 1> |
| using gridfire::MultiscalePartitioningEngineView::EigenFunctor::JacobianType = Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> |
| using gridfire::MultiscalePartitioningEngineView::EigenFunctor::OutputType = Eigen::Matrix<double, Eigen::Dynamic, 1> |
|
inline |
Constructs an EigenFunctor.
| view | The MultiscalePartitioningEngineView instance. |
| qse_solve_species | Species to solve for in the QSE group. |
| initial_comp | Initial abundances of all species in the full network. |
| T9 | Temperature in units of 10^9 K. |
| rho | Density in g/cm^3. |
| Y_scale | Scaling factors for the species abundances. |
| qse_solve_species_index_map | Mapping from species to their indices in the QSE solve vector. |
| int gridfire::MultiscalePartitioningEngineView::EigenFunctor::df | ( | const InputType & | v_qse, |
| JacobianType & | J_qse ) const |
Evaluates the Jacobian of the functor, J_qse = d(f_qse)/d(v_qse).
| v_qse | The input vector of scaled algebraic abundances. |
| J_qse | The output Jacobian matrix. |
|
inlinenodiscard |
Gets the number of input values to the functor (size of the variable vector).
| int gridfire::MultiscalePartitioningEngineView::EigenFunctor::operator() | ( | const InputType & | v_qse, |
| OutputType & | f_qse ) const |
Evaluates the functor's residual vector f_qse = dY_alg/dt.
| v_qse | The input vector of scaled algebraic abundances. |
| f_qse | The output residual vector. |
|
inlinenodiscard |
Gets the number of output values from the functor (size of the residual vector).
| const fourdst::composition::Composition& gridfire::MultiscalePartitioningEngineView::EigenFunctor::m_initial_comp |
Initial abundances of all species in the full network.
| const std::set<fourdst::atomic::Species>& gridfire::MultiscalePartitioningEngineView::EigenFunctor::m_qse_solve_species |
The set of species to solve for in the QSE group.
| const std::unordered_map<fourdst::atomic::Species, size_t> gridfire::MultiscalePartitioningEngineView::EigenFunctor::m_qse_solve_species_index_map |
Mapping from species to their indices in the QSE solve vector.
| const double gridfire::MultiscalePartitioningEngineView::EigenFunctor::m_rho |
Density in g/cm^3.
| const double gridfire::MultiscalePartitioningEngineView::EigenFunctor::m_T9 |
Temperature in units of 10^9 K.
| MultiscalePartitioningEngineView* gridfire::MultiscalePartitioningEngineView::EigenFunctor::m_view |
Pointer to the MultiscalePartitioningEngineView instance.
| const Eigen::VectorXd& gridfire::MultiscalePartitioningEngineView::EigenFunctor::m_Y_scale |
Scaling factors for the species abundances, used to improve solver stability.