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 "MFEMHypreAMS.h" 13 : #include "MFEMProblem.h" 14 : 15 : registerMooseObject("MooseApp", MFEMHypreAMS); 16 : 17 : InputParameters 18 2312 : MFEMHypreAMS::validParams() 19 : { 20 2312 : InputParameters params = MFEMSolverBase::validParams(); 21 4624 : params.addClassDescription("Hypre auxiliary-space Maxwell solver and preconditioner for the " 22 : "iterative solution of MFEM equation systems."); 23 9248 : params.addParam<MFEMFESpaceName>("fespace", "H(curl) FESpace to use in HypreAMS setup."); 24 6936 : params.addParam<bool>("singular", 25 4624 : false, 26 : "Declare that the system is singular; use when solving curl-curl problem " 27 : "if mass term is zero"); 28 6936 : params.addParam<int>("print_level", 2, "Set the solver verbosity."); 29 : 30 2312 : return params; 31 0 : } 32 : 33 107 : MFEMHypreAMS::MFEMHypreAMS(const InputParameters & parameters) 34 : : MFEMSolverBase(parameters), 35 107 : _mfem_fespace(getMFEMProblem().getMFEMObject<MFEMFESpace>("MFEMFESpace", 36 321 : getParam<MFEMFESpaceName>("fespace"))) 37 : { 38 107 : constructSolver(); 39 107 : } 40 : 41 : void 42 107 : MFEMHypreAMS::constructSolver() 43 : { 44 107 : auto solver = std::make_unique<mfem::HypreAMS>(_mfem_fespace.getFESpace().get()); 45 321 : if (getParam<bool>("singular")) 46 27 : solver->SetSingularProblem(); 47 : 48 214 : solver->SetPrintLevel(getParam<int>("print_level")); 49 : 50 107 : _solver = std::move(solver); 51 107 : } 52 : 53 : void 54 128 : MFEMHypreAMS::updateSolver(mfem::ParBilinearForm & a, mfem::Array<int> & tdofs) 55 : { 56 128 : if (_lor) 57 : { 58 4 : checkSpectralEquivalence(a); 59 4 : if (_mfem_fespace.getFESpace()->GetMesh()->GetElement(0)->GetGeometryType() != 60 : mfem::Geometry::Type::CUBE) 61 0 : mooseError("LOR HypreAMS Solver only supports hex meshes."); 62 : 63 4 : auto lor_solver = new mfem::LORSolver<mfem::HypreAMS>(a, tdofs); 64 8 : lor_solver->GetSolver().SetPrintLevel(getParam<int>("print_level")); 65 12 : if (getParam<bool>("singular")) 66 0 : lor_solver->GetSolver().SetSingularProblem(); 67 : 68 4 : _solver.reset(lor_solver); 69 : } 70 128 : } 71 : 72 : #endif