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