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 "MFEMHypreFGMRES.h" 13 : #include "MFEMProblem.h" 14 : 15 : registerMooseObject("MooseApp", MFEMHypreFGMRES); 16 : 17 : InputParameters 18 8634 : MFEMHypreFGMRES::validParams() 19 : { 20 8634 : InputParameters params = MFEMSolverBase::validParams(); 21 8634 : params.addClassDescription("Hypre solver for the iterative solution of MFEM equation systems " 22 : "using the flexible generalized minimal residual method."); 23 8634 : params.addParam<mfem::real_t>("l_tol", 1e-5, "Set the relative tolerance."); 24 8634 : params.addParam<int>("l_max_its", 10000, "Set the maximum number of iterations."); 25 8634 : params.addParam<int>("kdim", 10, "Set the k-dimension."); 26 8634 : params.addParam<int>("print_level", 2, "Set the solver verbosity."); 27 8634 : params.addParam<UserObjectName>("preconditioner", "Optional choice of preconditioner to use."); 28 : 29 8634 : return params; 30 0 : } 31 : 32 2 : MFEMHypreFGMRES::MFEMHypreFGMRES(const InputParameters & parameters) : MFEMSolverBase(parameters) 33 : { 34 2 : constructSolver(parameters); 35 2 : } 36 : 37 : void 38 2 : MFEMHypreFGMRES::constructSolver(const InputParameters &) 39 : { 40 : auto solver = 41 2 : std::make_unique<mfem::HypreFGMRES>(getMFEMProblem().mesh().getMFEMParMesh().GetComm()); 42 2 : solver->SetTol(getParam<mfem::real_t>("l_tol")); 43 2 : solver->SetMaxIter(getParam<int>("l_max_its")); 44 2 : solver->SetKDim(getParam<int>("kdim")); 45 2 : solver->SetPrintLevel(getParam<int>("print_level")); 46 2 : setPreconditioner(*solver); 47 2 : _solver = std::move(solver); 48 2 : } 49 : 50 : void 51 2 : MFEMHypreFGMRES::updateSolver(mfem::ParBilinearForm & a, mfem::Array<int> & tdofs) 52 : { 53 2 : if (_lor && _preconditioner) 54 0 : mooseError("LOR solver cannot take a preconditioner"); 55 : 56 2 : if (_preconditioner) 57 : { 58 0 : _preconditioner->updateSolver(a, tdofs); 59 0 : setPreconditioner(static_cast<mfem::HypreFGMRES &>(*_solver)); 60 : } 61 2 : else if (_lor) 62 : { 63 1 : if (!checkSpectralEquivalence(a)) 64 0 : mooseError("Low-Order-Refined solver requires the FESpace closed_basis to be GaussLobatto " 65 : "and the open-basis to be IntegratedGLL for ND and RT elements."); 66 : 67 1 : mfem::ParLORDiscretization lor_disc(a, tdofs); 68 : auto lor_solver = new mfem::LORSolver<mfem::HypreFGMRES>( 69 1 : lor_disc, getMFEMProblem().mesh().getMFEMParMesh().GetComm()); 70 1 : lor_solver->GetSolver().SetTol(getParam<mfem::real_t>("l_tol")); 71 1 : lor_solver->GetSolver().SetMaxIter(getParam<int>("l_max_its")); 72 1 : lor_solver->GetSolver().SetKDim(getParam<int>("kdim")); 73 1 : lor_solver->GetSolver().SetPrintLevel(getParam<int>("print_level")); 74 : 75 1 : _solver.reset(lor_solver); 76 1 : } 77 2 : } 78 : 79 : #endif