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 : #pragma once 13 : 14 : #include "libmesh/ignore_warnings.h" 15 : #include "mfem/miniapps/common/pfem_extras.hpp" 16 : #include "libmesh/restore_warnings.h" 17 : 18 : #include "MooseError.h" 19 : 20 : namespace Moose::MFEM 21 : { 22 : 23 : /// NonlinearFormIntegrator which scales its results by a constant value 24 : class NLScaleIntegrator : public mfem::NonlinearFormIntegrator 25 : { 26 : private: 27 : mfem::NonlinearFormIntegrator * _integrator{nullptr}; 28 : mfem::real_t _scale; 29 : bool _own_integrator; 30 : 31 : public: 32 : NLScaleIntegrator(mfem::NonlinearFormIntegrator * integ) 33 : : _integrator{integ}, _scale{1}, _own_integrator{true} 34 : { 35 : } 36 : NLScaleIntegrator(mfem::NonlinearFormIntegrator * integ, mfem::real_t scale) 37 : : _integrator{integ}, _scale{scale}, _own_integrator{true} 38 : { 39 : } 40 336 : NLScaleIntegrator(mfem::NonlinearFormIntegrator * integ, mfem::real_t scale, bool own) 41 336 : : _integrator{integ}, _scale{scale}, _own_integrator{own} 42 : { 43 336 : } 44 : 45 : void SetIntegrator(mfem::NonlinearFormIntegrator * integ) 46 : { 47 : if (_integrator && _own_integrator) 48 : { 49 : delete _integrator; 50 : } 51 : 52 : _integrator = integ; 53 : } 54 : 55 : void SetScale(mfem::real_t scale) { _scale = scale; } 56 : 57 : void SetOwn(bool own) { _own_integrator = own; } 58 : 59 77664 : void CheckIntegrator() 60 : { 61 77664 : if (!_integrator) 62 0 : mooseError("Integrator not set"); 63 77664 : } 64 : 65 : void SetIntegrationMode(mfem::NonlinearFormIntegrator::Mode m) 66 : { 67 : integrationMode = m; 68 : _integrator->SetIntegrationMode(m); 69 : } 70 : 71 : bool Patchwise() const { return integrationMode != Mode::ELEMENTWISE; } 72 : 73 : /// Set the memory type used for GeometricFactors and other large allocations 74 : /// in PA extensions. 75 : void SetPAMemoryType(mfem::MemoryType mt) 76 : { 77 : pa_mt = mt; 78 : _integrator->SetPAMemoryType(mt); 79 : } 80 : 81 : virtual void AssembleElementVector(const mfem::FiniteElement & el, 82 : mfem::ElementTransformation & Tr, 83 : const mfem::Vector & elfun, 84 : mfem::Vector & elvect) override; 85 : 86 : virtual void AssembleFaceVector(const mfem::FiniteElement & el1, 87 : const mfem::FiniteElement & el2, 88 : mfem::FaceElementTransformations & Tr, 89 : const mfem::Vector & elfun, 90 : mfem::Vector & elvect) override; 91 : 92 : virtual void AssembleElementGrad(const mfem::FiniteElement & el, 93 : mfem::ElementTransformation & Tr, 94 : const mfem::Vector & elfun, 95 : mfem::DenseMatrix & elmat) override; 96 : 97 : virtual void AssembleFaceGrad(const mfem::FiniteElement & el1, 98 : const mfem::FiniteElement & el2, 99 : mfem::FaceElementTransformations & Tr, 100 : const mfem::Vector & elfun, 101 : mfem::DenseMatrix & elmat) override; 102 : 103 : virtual mfem::real_t GetElementEnergy(const mfem::FiniteElement & el, 104 : mfem::ElementTransformation & Tr, 105 : const mfem::Vector & elfun) override; 106 : 107 672 : virtual ~NLScaleIntegrator() 108 336 : { 109 336 : if (_own_integrator) 110 336 : delete _integrator; 111 672 : }; 112 : }; 113 : } // namespace Moose::MFEM 114 : 115 : #endif