Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 3 : //* 4 : //* All rights reserved, see COPYRIGHT for full restrictions 5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT 6 : //* 7 : //* Licensed under LGPL 2.1, please see LICENSE for details 8 : //* https://www.gnu.org/licenses/lgpl-2.1.html 9 : 10 : #ifdef MFEM_ENABLED 11 : 12 : #include "MFEMSolverBase.h" 13 : 14 : InputParameters 15 86806 : MFEMSolverBase::validParams() 16 : { 17 86806 : InputParameters params = MFEMGeneralUserObject::validParams(); 18 86806 : params.addClassDescription("Base class for defining mfem::Solver derived classes for Moose."); 19 86806 : params.registerBase("MFEMSolverBase"); 20 86806 : params.addParam<bool>("low_order_refined", false, "Set usage of Low-Order Refined solver."); 21 : 22 86806 : return params; 23 0 : } 24 : 25 253 : MFEMSolverBase::MFEMSolverBase(const InputParameters & parameters) 26 : : MFEMGeneralUserObject(parameters), 27 253 : _lor{getParam<bool>("low_order_refined")}, 28 253 : _solver{nullptr}, 29 253 : _preconditioner{nullptr} 30 : { 31 253 : } 32 : 33 : template <typename T> 34 : void 35 254 : MFEMSolverBase::setPreconditioner(T & solver) 36 : { 37 254 : if (isParamSetByUser("preconditioner")) 38 : { 39 240 : if (!_preconditioner) 40 98 : _preconditioner = 41 98 : &const_cast<MFEMSolverBase &>(getUserObject<MFEMSolverBase>("preconditioner")); 42 : 43 240 : auto & mfem_pre = _preconditioner->getSolver(); 44 : if constexpr (std::is_base_of_v<mfem::HypreSolver, T>) 45 160 : if (auto * const hypre_pre = dynamic_cast<mfem::HypreSolver *>(&mfem_pre)) 46 160 : solver.SetPreconditioner(*hypre_pre); 47 : else 48 0 : mooseError("hypre solver preconditioners must themselves be hypre solvers"); 49 : else 50 80 : solver.SetPreconditioner(mfem_pre); 51 : } 52 254 : } 53 : 54 : template void MFEMSolverBase::setPreconditioner(mfem::CGSolver &); 55 : template void MFEMSolverBase::setPreconditioner(mfem::GMRESSolver &); 56 : template void MFEMSolverBase::setPreconditioner(mfem::HypreFGMRES &); 57 : template void MFEMSolverBase::setPreconditioner(mfem::HypreGMRES &); 58 : template void MFEMSolverBase::setPreconditioner(mfem::HyprePCG &); 59 : 60 : bool 61 34 : MFEMSolverBase::checkSpectralEquivalence(mfem::ParBilinearForm & blf) const 62 : { 63 34 : bool equiv = true; 64 : 65 34 : if (auto fec = dynamic_cast<const mfem::ND_FECollection *>(blf.FESpace()->FEColl())) 66 : { 67 8 : if (fec->GetClosedBasisType() != mfem::BasisType::GaussLobatto || 68 4 : fec->GetOpenBasisType() != mfem::BasisType::IntegratedGLL) 69 0 : equiv = false; 70 : } 71 30 : else if (auto fec = dynamic_cast<const mfem::RT_FECollection *>(blf.FESpace()->FEColl())) 72 : { 73 8 : if (fec->GetClosedBasisType() != mfem::BasisType::GaussLobatto || 74 4 : fec->GetOpenBasisType() != mfem::BasisType::IntegratedGLL) 75 0 : equiv = false; 76 : } 77 : 78 34 : return equiv; 79 : } 80 : 81 : #endif