https://mooseframework.inl.gov
MFEMHypreFGMRES.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 "MFEMHypreFGMRES.h"
13 #include "MFEMProblem.h"
14 
16 
19 {
21  params.addClassDescription("Hypre solver for the iterative solution of MFEM equation systems "
22  "using the flexible generalized minimal residual method.");
23  params.addParam<mfem::real_t>("l_tol", 1e-5, "Set the relative tolerance.");
24  params.addParam<int>("l_max_its", 10000, "Set the maximum number of iterations.");
25  params.addParam<int>("kdim", 10, "Set the k-dimension.");
26  params.addParam<int>("print_level", 2, "Set the solver verbosity.");
27  params.addParam<UserObjectName>("preconditioner", "Optional choice of preconditioner to use.");
28 
29  return params;
30 }
31 
33 {
35 }
36 
37 void
39 {
40  auto solver = std::make_unique<mfem::HypreFGMRES>(getMFEMProblem().getComm());
41  solver->SetTol(getParam<mfem::real_t>("l_tol"));
42  solver->SetMaxIter(getParam<int>("l_max_its"));
43  solver->SetKDim(getParam<int>("kdim"));
44  solver->SetPrintLevel(getParam<int>("print_level"));
45  setPreconditioner(*solver);
46  _solver = std::move(solver);
47 }
48 
49 void
50 MFEMHypreFGMRES::updateSolver(mfem::ParBilinearForm & a, mfem::Array<int> & tdofs)
51 {
52  if (_lor && _preconditioner)
53  mooseError("LOR solver cannot take a preconditioner");
54 
55  if (_preconditioner)
56  {
57  _preconditioner->updateSolver(a, tdofs);
58  setPreconditioner(static_cast<mfem::HypreFGMRES &>(*_solver));
59  }
60  else if (_lor)
61  {
63  mfem::ParLORDiscretization lor_disc(a, tdofs);
64  auto lor_solver = new mfem::LORSolver<mfem::HypreFGMRES>(lor_disc, getMFEMProblem().getComm());
65  lor_solver->GetSolver().SetTol(getParam<mfem::real_t>("l_tol"));
66  lor_solver->GetSolver().SetMaxIter(getParam<int>("l_max_its"));
67  lor_solver->GetSolver().SetKDim(getParam<int>("kdim"));
68  lor_solver->GetSolver().SetPrintLevel(getParam<int>("print_level"));
69 
70  _solver.reset(lor_solver);
71  }
72 }
73 
74 #endif
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
virtual void checkSpectralEquivalence(mfem::ParBilinearForm &blf) const
Checks for the correct configuration of quadrature bases for LOR spectral equivalence.
MFEMHypreFGMRES(const InputParameters &parameters)
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...
registerMooseObject("MooseApp", MFEMHypreFGMRES)
static InputParameters validParams()
MPI_Comm getComm()
Return the MPI communicator associated with this FE problem&#39;s mesh.
Definition: MFEMProblem.h:192
void constructSolver(const InputParameters &parameters) override
Override in derived classes to construct and set the solver options.
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:271
Wrapper for mfem::HypreFGMRES solver.
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.