10 #ifdef MOOSE_MFEM_ENABLED 20 params.
addClassDescription(
"Hypre auxiliary-space Maxwell solver and preconditioner for the " 21 "iterative solution of MFEM equation systems.");
22 params.
addParam<UserObjectName>(
"fespace",
"H(curl) FESpace to use in HypreAMS setup.");
25 "Declare that the system is singular; use when solving curl-curl problem " 26 "if mass term is zero");
27 params.
addParam<
int>(
"print_level", 2,
"Set the solver verbosity.");
42 if (getParam<bool>(
"singular"))
43 solver->SetSingularProblem();
45 solver->SetPrintLevel(getParam<int>(
"print_level"));
56 mooseError(
"Low-Order-Refined solver requires the FESpace closed_basis to be GaussLobatto " 57 "and the open-basis to be IntegratedGLL for ND and RT elements.");
60 mfem::Geometry::Type::CUBE)
61 mooseError(
"LOR HypreAMS Solver only supports hex meshes.");
63 auto lor_solver =
new mfem::LORSolver<mfem::HypreAMS>(a, tdofs);
64 lor_solver->GetSolver().SetPrintLevel(getParam<int>(
"print_level"));
65 if (getParam<bool>(
"singular"))
66 lor_solver->GetSolver().SetSingularProblem();
Wrapper for mfem::HypreAMS solver.
registerMooseObject("MooseApp", MFEMHypreAMS)
const MFEMFESpace & _mfem_fespace
const InputParameters & parameters() const
Get the parameters of the object.
static InputParameters validParams()
static InputParameters validParams()
Constructs and stores an mfem::ParFiniteElementSpace object.
virtual bool checkSpectralEquivalence(mfem::ParBilinearForm &blf) const
Checks for the correct configuration of quadrature bases for LOR spectral equivalence.
void updateSolver(mfem::ParBilinearForm &a, mfem::Array< int > &tdofs) override
Updates the solver with the bilinear form in case LOR solve is required.
MFEMHypreAMS(const InputParameters &)
bool _lor
Variable defining whether to use LOR solver.
Base class for wrapping mfem::Solver-derived classes.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
void constructSolver(const InputParameters ¶meters) override
Override in derived classes to construct and set the solver options.
std::shared_ptr< mfem::ParFiniteElementSpace > getFESpace() const
Returns a shared pointer to the constructed fespace.
std::unique_ptr< mfem::Solver > _solver
Solver to be used for the problem.