https://mooseframework.inl.gov
MFEMHypreGMRES.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 "MFEMHypreGMRES.h"
13 #include "MFEMProblem.h"
14 
16 
19 {
21  params.addClassDescription("Hypre solver for the iterative solution of MFEM equation systems "
22  "using the generalized minimal residual method.");
23 
24  params.addParam<mfem::real_t>("l_tol", 1e-5, "Set the relative tolerance.");
25  params.addParam<mfem::real_t>("l_abs_tol", 1e-50, "Set the absolute tolerance.");
26  params.addParam<int>("l_max_its", 10000, "Set the maximum number of iterations.");
27  params.addParam<int>("kdim", 10, "Set the k-dimension.");
28  params.addParam<int>("print_level", 2, "Set the solver verbosity.");
29  params.addParam<UserObjectName>("preconditioner", "Optional choice of preconditioner to use.");
30 
31  return params;
32 }
33 
35 {
37 }
38 
39 void
41 {
42  auto solver =
43  std::make_unique<mfem::HypreGMRES>(getMFEMProblem().mesh().getMFEMParMesh().GetComm());
44  solver->SetTol(getParam<mfem::real_t>("l_tol"));
45  solver->SetAbsTol(getParam<mfem::real_t>("l_abs_tol"));
46  solver->SetMaxIter(getParam<int>("l_max_its"));
47  solver->SetKDim(getParam<int>("kdim"));
48  solver->SetPrintLevel(getParam<int>("print_level"));
49  setPreconditioner(*solver);
50  _solver = std::move(solver);
51 }
52 
53 void
54 MFEMHypreGMRES::updateSolver(mfem::ParBilinearForm & a, mfem::Array<int> & tdofs)
55 {
56  if (_lor && _preconditioner)
57  mooseError("LOR solver cannot take a preconditioner");
58 
59  if (_preconditioner)
60  {
61  _preconditioner->updateSolver(a, tdofs);
62  setPreconditioner(static_cast<mfem::HypreGMRES &>(*_solver));
63  }
64  else if (_lor)
65  {
67  mooseError("Low-Order-Refined solver requires the FESpace closed_basis to be GaussLobatto "
68  "and the open-basis to be IntegratedGLL for ND and RT elements.");
69 
70  mfem::ParLORDiscretization lor_disc(a, tdofs);
71  auto lor_solver = new mfem::LORSolver<mfem::HypreGMRES>(
72  lor_disc, getMFEMProblem().mesh().getMFEMParMesh().GetComm());
73  lor_solver->GetSolver().SetTol(getParam<mfem::real_t>("l_tol"));
74  lor_solver->GetSolver().SetAbsTol(getParam<mfem::real_t>("l_abs_tol"));
75  lor_solver->GetSolver().SetMaxIter(getParam<int>("l_max_its"));
76  lor_solver->GetSolver().SetKDim(getParam<int>("kdim"));
77  lor_solver->GetSolver().SetPrintLevel(getParam<int>("print_level"));
78 
79  _solver.reset(lor_solver);
80  }
81 }
82 
83 #endif
static InputParameters validParams()
virtual MFEMMesh & mesh() override
Overwritten mesh() method from base MooseMesh to retrieve the correct mesh type, in this case MFEMMes...
Definition: MFEMProblem.C:466
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:127
MFEMSolverBase * _preconditioner
Preconditioner to be used for the problem.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
registerMooseObject("MooseApp", MFEMHypreGMRES)
static InputParameters validParams()
void setPreconditioner(T &solver)
Retrieves the preconditioner userobject if present, sets the member pointer to said object if still u...
MFEMHypreGMRES(const InputParameters &)
Wrapper for mfem::HypreGMRES solver.
void constructSolver(const InputParameters &parameters) override
Override in derived classes to construct and set the solver options.
virtual bool checkSpectralEquivalence(mfem::ParBilinearForm &blf) const
Checks for the correct configuration of quadrature bases for LOR spectral equivalence.
mfem::ParMesh & getMFEMParMesh()
Accessors for the _mfem_par_mesh object.
Definition: MFEMMesh.h:45
bool _lor
Variable defining whether to use LOR solver.
MFEMProblem & getMFEMProblem()
Returns a reference to the MFEMProblem instance.
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:267
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...
virtual void updateSolver(mfem::ParBilinearForm &a, mfem::Array< int > &tdofs)=0
Updates the solver with the given bilinear form and essential dof list, in case an LOR or algebraic s...
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...
void updateSolver(mfem::ParBilinearForm &a, mfem::Array< int > &tdofs) override
Updates the solver with the bilinear form in case LOR solve is required.
std::unique_ptr< mfem::Solver > _solver
Solver to be used for the problem.