/* *********************************************************************** // // Copyright (C) 2025 -- The 4D-STAR Collaboration // File Author: Emily Boudreaux // Last Modified: April 21, 2025 // // 4DSSE is free software; you can use it and/or modify // it under the terms and restrictions the GNU General Library Public // License version 3 (GPLv3) as published by the Free Software Foundation. // // 4DSSE is distributed in the hope that it will be useful, // but WITHOUT ANY WARRANTY; without even the implied warranty of // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. // See the GNU Library General Public License for more details. // // You should have received a copy of the GNU Library General Public License // along with this software; if not, write to the Free Software // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA // // *********************************************************************** */ #pragma once #include "mfem.hpp" #include "4DSTARTypes.h" #include #include "probe.h" class SchurCompliment final : public mfem::Operator { public: SchurCompliment(const mfem::SparseMatrix &QOp, const mfem::SparseMatrix &DOp, const mfem::SparseMatrix &MOp, const mfem::Solver &GradInvOp); SchurCompliment(const mfem::SparseMatrix &QOp, const mfem::SparseMatrix &DOp, const mfem::SparseMatrix &MOp); ~SchurCompliment() override = default; void Mult(const mfem::Vector &x, mfem::Vector &y) const override; void SetOperator(const mfem::SparseMatrix &QOp, const mfem::SparseMatrix &DOp, const mfem::SparseMatrix &MOp, const mfem::Solver & GradInvOp); void updateInverseNonlinearJacobian(const mfem::Solver &gradInv); private: void updateConstantTerms(const mfem::SparseMatrix &QOp, const mfem::SparseMatrix &DOp, const mfem::SparseMatrix &MOp); private: // Note that these are not owned by this class const mfem::SparseMatrix* m_QOp = nullptr; const mfem::SparseMatrix* m_DOp = nullptr; const mfem::SparseMatrix* m_MOp = nullptr; const mfem::Solver* m_GradInvOp = nullptr; int m_nPhi = 0; int m_nTheta = 0; mutable std::unique_ptr m_matrixForm; }; class GMRESInverter final : public mfem::Operator { public: explicit GMRESInverter(const SchurCompliment& op); ~GMRESInverter() override = default; void Mult(const mfem::Vector &x, mfem::Vector &y) const override; private: const SchurCompliment& m_op; mfem::GMRESSolver m_solver; }; class PolytropeOperator final : public mfem::Operator { public: PolytropeOperator( std::unique_ptr M, std::unique_ptr Q, std::unique_ptr D, std::unique_ptr f, const mfem::Array &blockOffsets); void populate_free_dof_array(); ~PolytropeOperator() override = default; void Mult(const mfem::Vector &xFree, mfem::Vector &yFree) const override; mfem::Operator& GetGradient(const mfem::Vector &xFree) const override; void set_essential_true_dofs(const SSE::MFEMArrayPair& theta_ess_tdofs, const SSE::MFEMArrayPair& phi_ess_tdofs); void set_essential_true_dofs(const SSE::MFEMArrayPairSet& ess_tdof_pair_set); SSE::MFEMArrayPairSet get_essential_true_dofs() const; bool isFinalized() const { return m_isFinalized; } void finalize(const mfem::Vector &initTheta); const mfem::Array& get_block_offsets() const { return m_blockOffsets; } const mfem::Array& get_reduced_block_offsets() const {return m_reducedBlockOffsets; } const mfem::BlockOperator &get_jacobian_operator() const; mfem::BlockDiagonalPreconditioner &get_preconditioner() const; int get_reduced_system_size() const; private: Probe::LogManager& m_logManager = Probe::LogManager::getInstance(); quill::Logger* m_logger = m_logManager.getLogger("log"); std::unique_ptr m_M; std::unique_ptr m_Q; std::unique_ptr m_D; std::unique_ptr m_f; std::unique_ptr m_Mmat; std::unique_ptr m_Qmat; std::unique_ptr m_Dmat; mutable mfem::Vector m_state; mfem::Array m_freeDofs; // These are used to store the reduced system matrices after essential dofs have been eliminated std::unique_ptr m_MReduced; std::unique_ptr m_QReduced; std::unique_ptr m_DReduced; const mfem::Array m_blockOffsets; mfem::Array m_reducedBlockOffsets; SSE::MFEMArrayPair m_theta_ess_tdofs; SSE::MFEMArrayPair m_phi_ess_tdofs; std::unique_ptr m_negQ_mat; mutable std::unique_ptr m_jacobian; mutable std::unique_ptr m_schurCompliment; mutable std::unique_ptr m_invSchurCompliment; mutable std::unique_ptr m_invNonlinearJacobian; mutable std::unique_ptr m_schurPreconditioner; bool m_isFinalized = false; private: void update_inverse_nonlinear_jacobian(const mfem::Operator &grad) const; void update_inverse_schur_compliment() const; void update_preconditioner(const mfem::Operator &grad) const; void scatter_boundary_conditions(); void construct_matrix_representations(); void construct_reduced_block_offsets(); void construct_jacobian_constant_terms(); };