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