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