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 : #pragma once 13 : #include "MFEMGeneralUserObject.h" 14 : #include "libmesh/ignore_warnings.h" 15 : #include <mfem.hpp> 16 : #include "libmesh/restore_warnings.h" 17 : #include <memory> 18 : 19 : /** 20 : * Base class for wrapping mfem::Solver-derived classes. 21 : */ 22 : class MFEMSolverBase : public MFEMGeneralUserObject 23 : { 24 : public: 25 : static InputParameters validParams(); 26 : 27 : MFEMSolverBase(const InputParameters & parameters); 28 : 29 : /// Retrieves the preconditioner userobject if present, sets the member pointer to 30 : /// said object if still unset, and sets the solver to use this preconditioner. 31 : template <typename T> 32 : void setPreconditioner(T & solver); 33 : 34 : /// Returns the wrapped MFEM solver 35 : mfem::Solver & getSolver(); 36 : 37 : /// Updates the solver with the given bilinear form and essential dof list, in case an LOR or algebraic solver is needed. 38 : virtual void updateSolver(mfem::ParBilinearForm & a, mfem::Array<int> & tdofs) = 0; 39 : 40 : /// Returns whether or not this solver (or its preconditioner) uses LOR 41 288 : bool isLOR() const { return _lor || (_preconditioner && _preconditioner->isLOR()); } 42 : 43 : protected: 44 : /// Override in derived classes to construct and set the solver options. 45 : virtual void constructSolver(const InputParameters & parameters) = 0; 46 : 47 : /// Checks for the correct configuration of quadrature bases for LOR spectral equivalence 48 : virtual bool checkSpectralEquivalence(mfem::ParBilinearForm & blf) const; 49 : 50 : // Variable defining whether to use LOR solver 51 : bool _lor; 52 : 53 : // Solver and preconditioner to be used for the problem 54 : std::unique_ptr<mfem::Solver> _solver; 55 : MFEMSolverBase * _preconditioner; 56 : }; 57 : 58 : inline mfem::Solver & 59 402 : MFEMSolverBase::getSolver() 60 : { 61 : mooseAssert(_solver, "Attempting to retrieve solver before it's been constructed"); 62 402 : return *_solver; 63 : } 64 : 65 : #endif