https://mooseframework.inl.gov
MFEMSolverBase.C
Go to the documentation of this file.
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 MOOSE_MFEM_ENABLED
11 
12 #include "MFEMSolverBase.h"
13 
16 {
18  params.addClassDescription("Base class for defining mfem::Solver derived classes for Moose.");
19  params.registerBase("MFEMSolverBase");
20  params.addParam<bool>("low_order_refined", false, "Set usage of Low-Order Refined solver.");
21 
22  return params;
23 }
24 
26  : MFEMGeneralUserObject(parameters),
27  _lor{getParam<bool>("low_order_refined")},
28  _solver{nullptr},
29  _preconditioner{nullptr}
30 {
31 }
32 
33 template <typename T>
34 void
36 {
37  if (isParamSetByUser("preconditioner"))
38  {
39  if (!_preconditioner)
41  &const_cast<MFEMSolverBase &>(getUserObject<MFEMSolverBase>("preconditioner"));
42 
43  auto & mfem_pre = _preconditioner->getSolver();
44  if constexpr (std::is_base_of_v<mfem::HypreSolver, T>)
45  if (auto * const hypre_pre = dynamic_cast<mfem::HypreSolver *>(&mfem_pre))
46  solver.SetPreconditioner(*hypre_pre);
47  else
48  mooseError("hypre solver preconditioners must themselves be hypre solvers");
49  else
50  solver.SetPreconditioner(mfem_pre);
51  }
52 }
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 
62 
63 void
64 MFEMSolverBase::checkSpectralEquivalence(mfem::ParBilinearForm & blf) const
65 {
66  if (auto fec = dynamic_cast<const mfem::H1_FECollection *>(blf.FESpace()->FEColl()))
67  {
68  if (fec->GetBasisType() != mfem::BasisType::GaussLobatto)
69  mooseError("Low-Order-Refined solver requires the FESpace basis to be GaussLobatto "
70  "for H1 elements.");
71  }
72  else if (auto fec = dynamic_cast<const mfem::ND_FECollection *>(blf.FESpace()->FEColl()))
73  {
74  if (fec->GetClosedBasisType() != mfem::BasisType::GaussLobatto ||
75  fec->GetOpenBasisType() != mfem::BasisType::IntegratedGLL)
76  mooseError("Low-Order-Refined solver requires the FESpace closed-basis to be GaussLobatto "
77  "and the open-basis to be IntegratedGLL for ND elements.");
78  }
79  else if (auto fec = dynamic_cast<const mfem::RT_FECollection *>(blf.FESpace()->FEColl()))
80  {
81  if (fec->GetClosedBasisType() != mfem::BasisType::GaussLobatto ||
82  fec->GetOpenBasisType() != mfem::BasisType::IntegratedGLL)
83  mooseError("Low-Order-Refined solver requires the FESpace closed-basis to be GaussLobatto "
84  "and the open-basis to be IntegratedGLL for RT elements.");
85  }
86 }
87 
88 #endif
static InputParameters validParams()
virtual void checkSpectralEquivalence(mfem::ParBilinearForm &blf) const
Checks for the correct configuration of quadrature bases for LOR spectral equivalence.
MFEMSolverBase * _preconditioner
Preconditioner to be used for the problem.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static InputParameters validParams()
void setPreconditioner(T &solver)
Retrieves the preconditioner userobject if present, sets the member pointer to said object if still u...
void registerBase(const std::string &value)
This method must be called from every base "Moose System" to create linkage with the Action System...
This class adds a getMFEMProblem method.
MFEMSolverBase(const InputParameters &parameters)
Patch for mfem::HypreGMRES to reset preconditioning matrix at every nonlinear/time iteration...
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...
Definition: MooseBase.h:271
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
Patch for mfem::HyprePCG to reset preconditioning matrix at every nonlinear/time iteration.
mfem::Solver & getSolver()
Returns the wrapped MFEM solver.
bool isParamSetByUser(const std::string &name) const
Test if the supplied parameter is set by a user, as opposed to not set or set to default.
Definition: MooseBase.h:205