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/miniapps/common/pfem_extras.hpp" 16 : #include "libmesh/restore_warnings.h" 17 : 18 : namespace Moose::MFEM 19 : { 20 : 21 : /// Integrator which scales its results by a constant value 22 : class ScaleIntegrator : public mfem::BilinearFormIntegrator 23 : { 24 : private: 25 : mfem::BilinearFormIntegrator * _integrator{nullptr}; 26 : double _scale; 27 : bool _own_integrator; 28 : 29 : public: 30 : ScaleIntegrator(mfem::BilinearFormIntegrator * integ) 31 : : _integrator{integ}, _scale{1}, _own_integrator{true} 32 : { 33 : } 34 2 : ScaleIntegrator(mfem::BilinearFormIntegrator * integ, double scale) 35 2 : : _integrator{integ}, _scale{scale}, _own_integrator{true} 36 : { 37 2 : } 38 109 : ScaleIntegrator(mfem::BilinearFormIntegrator * integ, double scale, bool own) 39 109 : : _integrator{integ}, _scale{scale}, _own_integrator{own} 40 : { 41 109 : } 42 : 43 : void SetIntegrator(mfem::BilinearFormIntegrator * integ) 44 : { 45 : if (_integrator && _own_integrator) 46 : { 47 : delete _integrator; 48 : } 49 : 50 : _integrator = integ; 51 : } 52 : 53 : void SetScale(double scale) { _scale = scale; } 54 : 55 : void SetOwn(bool own) { _own_integrator = own; } 56 : 57 161730 : void CheckIntegrator() 58 : { 59 161730 : if (!_integrator) 60 0 : mooseError("Integrator not set"); 61 161730 : } 62 : 63 : virtual void SetIntRule(const mfem::IntegrationRule * ir); 64 : 65 : virtual void AssembleElementMatrix(const mfem::FiniteElement & el, 66 : mfem::ElementTransformation & Trans, 67 : mfem::DenseMatrix & elmat); 68 : virtual void AssembleElementMatrix2(const mfem::FiniteElement & trial_fe, 69 : const mfem::FiniteElement & test_fe, 70 : mfem::ElementTransformation & Trans, 71 : mfem::DenseMatrix & elmat); 72 : 73 : using mfem::BilinearFormIntegrator::AssembleFaceMatrix; 74 : virtual void AssembleFaceMatrix(const mfem::FiniteElement & el1, 75 : const mfem::FiniteElement & el2, 76 : mfem::FaceElementTransformations & Trans, 77 : mfem::DenseMatrix & elmat); 78 : 79 : using mfem::BilinearFormIntegrator::AssemblePA; 80 : virtual void AssemblePA(const mfem::FiniteElementSpace & fes); 81 : 82 : virtual void AssembleDiagonalPA(mfem::Vector & diag); 83 : 84 : virtual void AssemblePAInteriorFaces(const mfem::FiniteElementSpace & fes); 85 : 86 : virtual void AssemblePABoundaryFaces(const mfem::FiniteElementSpace & fes); 87 : 88 : virtual void AddMultTransposePA(const mfem::Vector & x, mfem::Vector & y) const; 89 : 90 : virtual void AddMultPA(const mfem::Vector & x, mfem::Vector & y) const; 91 : 92 : virtual void 93 : AssembleEA(const mfem::FiniteElementSpace & fes, mfem::Vector & emat, const bool add); 94 : 95 : virtual ~ScaleIntegrator(); 96 : }; 97 : } // namespace Moose::MFEM 98 : 99 : #endif