SERiF 0.0.1a
3+1D Stellar Structure and Evolution
|
Represents the coupled nonlinear operator for the polytropic system. More...
#include <polytropeOperator.h>
Public Member Functions | |
PolytropeOperator (std::unique_ptr< mfem::MixedBilinearForm > M, std::unique_ptr< mfem::MixedBilinearForm > Q, std::unique_ptr< mfem::BilinearForm > D, std::unique_ptr< mfem::BilinearForm > S, std::unique_ptr< mfem::NonlinearForm > f, const mfem::Array< int > &blockOffsets) | |
Constructs the PolytropeOperator. | |
~PolytropeOperator () override=default | |
Destructor. | |
void | Mult (const mfem::Vector &xFree, mfem::Vector &yFree) const override |
Applies the PolytropeOperator: y = F(x). This computes the residual of the nonlinear system. | |
mfem::Operator & | GetGradient (const mfem::Vector &xFree) const override |
Computes the Jacobian of the PolytropeOperator at a given state xFree. The Jacobian is of the form: J = [ G M ] [ -Q D ] where G is the gradient of f(θ). | |
void | finalize (const mfem::Vector &initTheta) |
Finalizes the operator setup. This must be called after setting essential true DOFs and before using Mult or GetGradient. It constructs reduced matrices, block offsets, and populates free DOFs. | |
bool | isFinalized () const |
Checks if the operator has been finalized. | |
void | set_essential_true_dofs (const serif::types::MFEMArrayPair &theta_ess_tdofs, const serif::types::MFEMArrayPair &phi_ess_tdofs) |
Sets the essential true degrees of freedom for both θ and φ variables. Marks the operator as not finalized. | |
void | set_essential_true_dofs (const serif::types::MFEMArrayPairSet &ess_tdof_pair_set) |
Sets the essential true degrees of freedom using a pair of pairs. Marks the operator as not finalized. | |
const mfem::Vector & | reconstruct_full_state_vector (const mfem::Vector &reducedState) const |
Reconstructs the full state vector (including essential DOFs) from a reduced state vector (free DOFs). | |
const mfem::BlockVector | reconstruct_full_block_state_vector (const mfem::Vector &reducedState) const |
void | populate_free_dof_array () |
Populates the internal array of free (non-essential) degree of freedom indices. This is called during finalize(). | |
const mfem::BlockOperator & | get_jacobian_operator () const |
— Getters for key internal state and operators — | |
mfem::BlockDiagonalPreconditioner & | get_preconditioner () const |
Gets the block diagonal preconditioner for the Schur complement system. Asserts that the operator is finalized and the preconditioner has been computed. | |
const mfem::Array< int > & | get_free_dofs () const |
— Getters for information on internal state — | |
int | get_reduced_system_size () const |
Gets the size of the reduced system (number of free DOFs). Asserts that the operator is finalized. | |
serif::types::MFEMArrayPairSet | get_essential_true_dofs () const |
Gets the currently set essential true degrees of freedom. | |
const mfem::Array< int > & | get_block_offsets () const |
Gets the block offsets for the full (unreduced) system. | |
const mfem::Array< int > & | get_reduced_block_offsets () const |
Gets the block offsets for the reduced system (after eliminating essential DOFs). | |
Private Member Functions | |
void | construct_matrix_representations () |
Constructs the sparse matrix representations of M, Q, and D, and their reduced forms. Called during finalize(). | |
void | construct_reduced_block_offsets () |
Constructs the block offsets for the reduced system. Called during finalize(). | |
void | construct_jacobian_constant_terms () |
Constructs the constant terms of the Jacobian operator (M, -Q, D). The (0,0) block (gradient of f) is set in GetGradient. Called during finalize(). | |
void | scatter_boundary_conditions () |
Scatters the values of essential boundary conditions into the full state vector. Called during finalize(). | |
void | update_inverse_nonlinear_jacobian (const mfem::Operator &grad) const |
Updates the solver for the inverse of the nonlinear Jacobian block (G_00). | |
void | update_inverse_schur_compliment () const |
Updates the inverse Schur complement operator and its components. This is typically called after the nonlinear Jacobian part has been updated. | |
void | update_preconditioner (const mfem::Operator &grad) const |
Updates the preconditioner components. This involves updating the inverse nonlinear Jacobian and then the inverse Schur complement. | |
Private Attributes | |
serif::probe::LogManager & | m_logManager = serif::probe::LogManager::getInstance() |
Reference to the global log manager. | |
quill::Logger * | m_logger = m_logManager.getLogger("log") |
Pointer to the specific logger instance. | |
std::unique_ptr< mfem::MixedBilinearForm > | m_M |
Bilinear form M, coupling θ and φ. | |
std::unique_ptr< mfem::MixedBilinearForm > | m_Q |
Bilinear form Q, coupling φ and θ. | |
std::unique_ptr< mfem::BilinearForm > | m_D |
Bilinear form D, acting on φ. | |
std::unique_ptr< mfem::BilinearForm > | m_S |
Bilinear form S, used for least squares stabilization. | |
std::unique_ptr< mfem::NonlinearForm > | m_f |
Nonlinear form f, acting on θ. | |
std::unique_ptr< mfem::SparseMatrix > | m_Mmat |
Sparse matrix representation of M. | |
std::unique_ptr< mfem::SparseMatrix > | m_Qmat |
Sparse matrix representation of Q. | |
std::unique_ptr< mfem::SparseMatrix > | m_Dmat |
Sparse matrix representation of D. | |
std::unique_ptr< mfem::SparseMatrix > | m_Smat |
std::unique_ptr< mfem::SparseMatrix > | m_MReduced |
Reduced M matrix (free DOFs only). | |
std::unique_ptr< mfem::SparseMatrix > | m_QReduced |
Reduced Q matrix (free DOFs only). | |
std::unique_ptr< mfem::SparseMatrix > | m_DReduced |
Reduced D matrix (free DOFs only). | |
std::unique_ptr< mfem::SparseMatrix > | m_SReduced |
Reduced S matrix (free DOFs only, used for least squares stabilization). | |
std::unique_ptr< mfem::SparseMatrix > | m_gradReduced |
Reduced gradient operator (G) for the nonlinear part f(θ). | |
std::unique_ptr< mfem::SparseMatrix > | m_MScaledReduced |
Scaled M matrix (free DOFs only, scaled by stabilization coefficient). | |
std::unique_ptr< mfem::SparseMatrix > | m_QScaledReduced |
Scaled Q matrix (free DOFs only, scaled by stabilization coefficient). | |
std::unique_ptr< mfem::SparseMatrix > | m_DScaledReduced |
Scaled D matrix (free DOFs only, scaled by stabilization coefficient). | |
std::unique_ptr< mfem::SparseMatrix > | m_ScaledSReduced |
Scaled S matrix (free DOFs only, scaled by stabilization coefficient). | |
mfem::Vector | m_state |
Full state vector [θ, φ]^T, including essential DOFs. | |
mfem::Array< int > | m_freeDofs |
Array of indices for free (non-essential) DOFs. | |
const mfem::Array< int > | m_blockOffsets |
Block offsets for the full system [0, size(θ), size(θ)+size(φ)]. | |
mfem::Array< int > | m_reducedBlockOffsets |
Block offsets for the reduced system (free DOFs). | |
serif::types::MFEMArrayPair | m_theta_ess_tdofs |
Essential true DOFs for θ (indices and values). | |
serif::types::MFEMArrayPair | m_phi_ess_tdofs |
Essential true DOFs for φ (indices and values). | |
std::unique_ptr< mfem::ScaledOperator > | m_negQ_mat |
Scaled operator for -Q_reduced. | |
std::unique_ptr< mfem::BlockOperator > | m_jacobian |
Jacobian operator J = [G M; -Q D]_reduced. | |
std::unique_ptr< SchurCompliment > | m_schurCompliment |
Schur complement S = D_reduced - Q_reduced * G_inv_reduced * M_reduced. | |
std::unique_ptr< GMRESInverter > | m_invSchurCompliment |
Approximate inverse of the Schur complement. | |
std::unique_ptr< mfem::Solver > | m_invNonlinearJacobian |
Solver for the inverse of the G block (gradient of f(θ)_reduced). | |
std::unique_ptr< mfem::BlockDiagonalPreconditioner > | m_schurPreconditioner |
Block diagonal preconditioner for the system. | |
bool | m_isFinalized = false |
Flag indicating if finalize() has been called. | |
Static Private Attributes | |
static constexpr double | m_stabilizationCoefficient = -2.0 |
Stabilization coefficient for the system, used to more tightly couple ∇θ and φ. | |
static constexpr double | m_IncrementedStabilizationCoefficient = 1.0 + m_stabilizationCoefficient |
1 + Stabilization coefficient for the system, used to more tightly couple ∇θ and φ. | |
Represents the coupled nonlinear operator for the polytropic system.
This operator defines the system F(X) = 0, where X = [θ, φ]^T, and F(X) = [ f(θ) + Mφ ] [ Dφ - Qθ ]. It handles essential boundary conditions and the construction of the Jacobian.
Definition at line 152 of file polytropeOperator.h.
serif::polytrope::PolytropeOperator::PolytropeOperator | ( | std::unique_ptr< mfem::MixedBilinearForm > | M, |
std::unique_ptr< mfem::MixedBilinearForm > | Q, | ||
std::unique_ptr< mfem::BilinearForm > | D, | ||
std::unique_ptr< mfem::BilinearForm > | S, | ||
std::unique_ptr< mfem::NonlinearForm > | f, | ||
const mfem::Array< int > & | blockOffsets ) |
Constructs the PolytropeOperator.
M | The M bilinear form (coupling θ and φ). |
Q | The Q bilinear form (coupling φ and θ). |
D | The D bilinear form (acting on φ). |
f | The nonlinear form f(θ). |
blockOffsets | Offsets defining the blocks for θ and φ variables. |
Definition at line 135 of file polytropeOperator.cpp.
|
overridedefault |
Destructor.
|
private |
Constructs the constant terms of the Jacobian operator (M, -Q, D). The (0,0) block (gradient of f) is set in GetGradient. Called during finalize().
Definition at line 429 of file polytropeOperator.cpp.
|
private |
Constructs the sparse matrix representations of M, Q, and D, and their reduced forms. Called during finalize().
Definition at line 367 of file polytropeOperator.cpp.
|
private |
Constructs the block offsets for the reduced system. Called during finalize().
Definition at line 422 of file polytropeOperator.cpp.
void serif::polytrope::PolytropeOperator::finalize | ( | const mfem::Vector & | initTheta | ) |
Finalizes the operator setup. This must be called after setting essential true DOFs and before using Mult or GetGradient. It constructs reduced matrices, block offsets, and populates free DOFs.
initTheta | Initial guess for θ, used to evaluate the nonlinear form if necessary during setup. |
Definition at line 262 of file polytropeOperator.cpp.
|
inline |
Gets the block offsets for the full (unreduced) system.
Definition at line 284 of file polytropeOperator.h.
serif::types::MFEMArrayPairSet serif::polytrope::PolytropeOperator::get_essential_true_dofs | ( | ) | const |
Gets the currently set essential true degrees of freedom.
Definition at line 303 of file polytropeOperator.cpp.
|
inline |
— Getters for information on internal state —
Gets the full state vector, including essential DOFs.
Definition at line 265 of file polytropeOperator.h.
const mfem::BlockOperator & serif::polytrope::PolytropeOperator::get_jacobian_operator | ( | ) | const |
— Getters for key internal state and operators —
Gets the Jacobian operator. Asserts that the operator is finalized and the Jacobian has been computed.
Definition at line 308 of file polytropeOperator.cpp.
mfem::BlockDiagonalPreconditioner & serif::polytrope::PolytropeOperator::get_preconditioner | ( | ) | const |
Gets the block diagonal preconditioner for the Schur complement system. Asserts that the operator is finalized and the preconditioner has been computed.
Definition at line 319 of file polytropeOperator.cpp.
|
inline |
Gets the block offsets for the reduced system (after eliminating essential DOFs).
Definition at line 290 of file polytropeOperator.h.
int serif::polytrope::PolytropeOperator::get_reduced_system_size | ( | ) | const |
Gets the size of the reduced system (number of free DOFs). Asserts that the operator is finalized.
Definition at line 329 of file polytropeOperator.cpp.
|
override |
Computes the Jacobian of the PolytropeOperator at a given state xFree. The Jacobian is of the form: J = [ G M ] [ -Q D ] where G is the gradient of f(θ).
xFree | The vector of free DOFs at which to evaluate the gradient. |
Definition at line 225 of file polytropeOperator.cpp.
|
inline |
Checks if the operator has been finalized.
Definition at line 206 of file polytropeOperator.h.
|
override |
Applies the PolytropeOperator: y = F(x). This computes the residual of the nonlinear system.
xFree | The vector of free (non-essential) degrees of freedom. |
yFree | The resulting residual vector corresponding to free DOFs. |
Definition at line 161 of file polytropeOperator.cpp.
void serif::polytrope::PolytropeOperator::populate_free_dof_array | ( | ) |
Populates the internal array of free (non-essential) degree of freedom indices. This is called during finalize().
Definition at line 349 of file polytropeOperator.cpp.
|
nodiscard |
@breif Reconstruct the full state vector (including essential DOFs) from a reduced state vector (free DOFs) as well as the block offsets.
reducedState | The vector containing only the free degrees of freedom. |
Definition at line 342 of file polytropeOperator.cpp.
|
nodiscard |
Reconstructs the full state vector (including essential DOFs) from a reduced state vector (free DOFs).
reducedState | The vector containing only the free degrees of freedom. |
Definition at line 337 of file polytropeOperator.cpp.
|
private |
Scatters the values of essential boundary conditions into the full state vector. Called during finalize().
Definition at line 437 of file polytropeOperator.cpp.
void serif::polytrope::PolytropeOperator::set_essential_true_dofs | ( | const serif::types::MFEMArrayPair & | theta_ess_tdofs, |
const serif::types::MFEMArrayPair & | phi_ess_tdofs ) |
Sets the essential true degrees of freedom for both θ and φ variables. Marks the operator as not finalized.
theta_ess_tdofs | Pair of arrays: (indices of essential DOFs for θ, values at these DOFs). |
phi_ess_tdofs | Pair of arrays: (indices of essential DOFs for φ, values at these DOFs). |
Definition at line 287 of file polytropeOperator.cpp.
void serif::polytrope::PolytropeOperator::set_essential_true_dofs | ( | const serif::types::MFEMArrayPairSet & | ess_tdof_pair_set | ) |
Sets the essential true degrees of freedom using a pair of pairs. Marks the operator as not finalized.
ess_tdof_pair_set | A pair containing the essential TDOF pairs for theta and phi. |
Definition at line 299 of file polytropeOperator.cpp.
|
private |
Updates the solver for the inverse of the nonlinear Jacobian block (G_00).
grad | The gradient operator (G_00) of the nonlinear part f(θ). |
Definition at line 460 of file polytropeOperator.cpp.
|
private |
Updates the inverse Schur complement operator and its components. This is typically called after the nonlinear Jacobian part has been updated.
Definition at line 464 of file polytropeOperator.cpp.
|
private |
Updates the preconditioner components. This involves updating the inverse nonlinear Jacobian and then the inverse Schur complement.
grad | The gradient operator (G_00) of the nonlinear part f(θ). |
Definition at line 487 of file polytropeOperator.cpp.
|
private |
Block offsets for the full system [0, size(θ), size(θ)+size(φ)].
Definition at line 333 of file polytropeOperator.h.
|
private |
Bilinear form D, acting on φ.
Definition at line 300 of file polytropeOperator.h.
|
private |
Sparse matrix representation of D.
Definition at line 307 of file polytropeOperator.h.
|
private |
Reduced D matrix (free DOFs only).
Definition at line 313 of file polytropeOperator.h.
|
private |
Scaled D matrix (free DOFs only, scaled by stabilization coefficient).
Definition at line 320 of file polytropeOperator.h.
|
private |
Nonlinear form f, acting on θ.
Definition at line 302 of file polytropeOperator.h.
|
private |
Array of indices for free (non-essential) DOFs.
Definition at line 330 of file polytropeOperator.h.
|
mutableprivate |
Reduced gradient operator (G) for the nonlinear part f(θ).
Definition at line 315 of file polytropeOperator.h.
|
staticconstexprprivate |
1 + Stabilization coefficient for the system, used to more tightly couple ∇θ and φ.
Definition at line 326 of file polytropeOperator.h.
|
mutableprivate |
Solver for the inverse of the G block (gradient of f(θ)_reduced).
Definition at line 345 of file polytropeOperator.h.
|
mutableprivate |
Approximate inverse of the Schur complement.
Definition at line 344 of file polytropeOperator.h.
|
private |
Flag indicating if finalize() has been called.
Definition at line 349 of file polytropeOperator.h.
|
mutableprivate |
Jacobian operator J = [G M; -Q D]_reduced.
Definition at line 342 of file polytropeOperator.h.
|
private |
Pointer to the specific logger instance.
Definition at line 295 of file polytropeOperator.h.
|
private |
Reference to the global log manager.
Definition at line 294 of file polytropeOperator.h.
|
private |
Bilinear form M, coupling θ and φ.
Definition at line 298 of file polytropeOperator.h.
|
private |
Sparse matrix representation of M.
Definition at line 305 of file polytropeOperator.h.
|
private |
Reduced M matrix (free DOFs only).
Definition at line 311 of file polytropeOperator.h.
|
private |
Scaled M matrix (free DOFs only, scaled by stabilization coefficient).
Definition at line 318 of file polytropeOperator.h.
|
private |
Scaled operator for -Q_reduced.
Definition at line 341 of file polytropeOperator.h.
|
private |
Essential true DOFs for φ (indices and values).
Definition at line 338 of file polytropeOperator.h.
|
private |
Bilinear form Q, coupling φ and θ.
Definition at line 299 of file polytropeOperator.h.
|
private |
Sparse matrix representation of Q.
Definition at line 306 of file polytropeOperator.h.
|
private |
Reduced Q matrix (free DOFs only).
Definition at line 312 of file polytropeOperator.h.
|
private |
Scaled Q matrix (free DOFs only, scaled by stabilization coefficient).
Definition at line 319 of file polytropeOperator.h.
|
private |
Block offsets for the reduced system (free DOFs).
Definition at line 334 of file polytropeOperator.h.
|
private |
Bilinear form S, used for least squares stabilization.
Definition at line 301 of file polytropeOperator.h.
|
private |
Scaled S matrix (free DOFs only, scaled by stabilization coefficient).
Definition at line 321 of file polytropeOperator.h.
|
mutableprivate |
Schur complement S = D_reduced - Q_reduced * G_inv_reduced * M_reduced.
Definition at line 343 of file polytropeOperator.h.
|
mutableprivate |
Block diagonal preconditioner for the system.
Definition at line 346 of file polytropeOperator.h.
|
private |
Definition at line 308 of file polytropeOperator.h.
|
private |
Reduced S matrix (free DOFs only, used for least squares stabilization).
Definition at line 314 of file polytropeOperator.h.
|
staticconstexprprivate |
Stabilization coefficient for the system, used to more tightly couple ∇θ and φ.
Definition at line 325 of file polytropeOperator.h.
|
mutableprivate |
Full state vector [θ, φ]^T, including essential DOFs.
Definition at line 329 of file polytropeOperator.h.
|
private |
Essential true DOFs for θ (indices and values).
Definition at line 337 of file polytropeOperator.h.