https://mooseframework.inl.gov
ScaleIntegrator.h
Go to the documentation of this file.
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 #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 
22 class ScaleIntegrator : public mfem::BilinearFormIntegrator
23 {
24 private:
25  mfem::BilinearFormIntegrator * _integrator{nullptr};
26  mfem::real_t _scale;
28 
29 public:
30  ScaleIntegrator(mfem::BilinearFormIntegrator * integ)
31  : _integrator{integ}, _scale{1}, _own_integrator{true}
32  {
33  }
34  ScaleIntegrator(mfem::BilinearFormIntegrator * integ, mfem::real_t scale)
35  : _integrator{integ}, _scale{scale}, _own_integrator{true}
36  {
37  }
38  ScaleIntegrator(mfem::BilinearFormIntegrator * integ, mfem::real_t scale, bool own)
39  : _integrator{integ}, _scale{scale}, _own_integrator{own}
40  {
41  }
42 
43  void SetIntegrator(mfem::BilinearFormIntegrator * integ)
44  {
46  {
47  delete _integrator;
48  }
49 
50  _integrator = integ;
51  }
52 
53  void SetScale(mfem::real_t scale) { _scale = scale; }
54 
55  void SetOwn(bool own) { _own_integrator = own; }
56 
58  {
59  if (!_integrator)
60  mooseError("Integrator not set");
61  }
62 
63  virtual void SetIntRule(const mfem::IntegrationRule * ir) override;
64 
65  virtual void AssembleElementMatrix(const mfem::FiniteElement & el,
66  mfem::ElementTransformation & Trans,
67  mfem::DenseMatrix & elmat) override;
68  virtual void AssembleElementMatrix2(const mfem::FiniteElement & trial_fe,
69  const mfem::FiniteElement & test_fe,
70  mfem::ElementTransformation & Trans,
71  mfem::DenseMatrix & elmat) override;
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) override;
78 
79  using mfem::BilinearFormIntegrator::AssemblePA;
80  virtual void AssemblePA(const mfem::FiniteElementSpace & fes) override;
81 
82  virtual void AssembleDiagonalPA(mfem::Vector & diag) override;
83 
84  virtual void AssemblePAInteriorFaces(const mfem::FiniteElementSpace & fes) override;
85 
86  virtual void AssemblePABoundaryFaces(const mfem::FiniteElementSpace & fes) override;
87 
88  virtual void AddMultTransposePA(const mfem::Vector & x, mfem::Vector & y) const override;
89 
90  virtual void AddMultPA(const mfem::Vector & x, mfem::Vector & y) const override;
91 
92  virtual void
93  AssembleEA(const mfem::FiniteElementSpace & fes, mfem::Vector & emat, const bool add) override;
94 
95  virtual void AssembleEABoundary(const mfem::FiniteElementSpace & fes,
96  mfem::Vector & emat,
97  const bool add) override;
98 
99  virtual void AssembleMF(const mfem::FiniteElementSpace & fes) override;
100 
101  virtual void AddMultMF(const mfem::Vector & x, mfem::Vector & y) const override;
102 
103  virtual void AssembleDiagonalMF(mfem::Vector & diag) override;
104 
105  virtual ~ScaleIntegrator();
106 };
107 } // namespace Moose::MFEM
108 
109 #endif
virtual void AssembleEA(const mfem::FiniteElementSpace &fes, mfem::Vector &emat, const bool add) override
virtual void AssembleDiagonalMF(mfem::Vector &diag) override
virtual void AssembleElementMatrix(const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, mfem::DenseMatrix &elmat) override
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
void scale(MeshBase &mesh, const Real xs, const Real ys=0., const Real zs=0.)
virtual void AssemblePA(const mfem::FiniteElementSpace &fes) override
virtual void SetIntRule(const mfem::IntegrationRule *ir) override
mfem::BilinearFormIntegrator * _integrator
virtual void AssembleMF(const mfem::FiniteElementSpace &fes) override
void SetScale(mfem::real_t scale)
virtual void AssembleEABoundary(const mfem::FiniteElementSpace &fes, mfem::Vector &emat, const bool add) override
virtual void AddMultPA(const mfem::Vector &x, mfem::Vector &y) const override
virtual void AssemblePABoundaryFaces(const mfem::FiniteElementSpace &fes) override
virtual void AssemblePAInteriorFaces(const mfem::FiniteElementSpace &fes) override
ScaleIntegrator(mfem::BilinearFormIntegrator *integ)
ScaleIntegrator(mfem::BilinearFormIntegrator *integ, mfem::real_t scale, bool own)
virtual void AssembleDiagonalPA(mfem::Vector &diag) override
Integrator which scales its results by a constant value.
virtual void AssembleFaceMatrix(const mfem::FiniteElement &el1, const mfem::FiniteElement &el2, mfem::FaceElementTransformations &Trans, mfem::DenseMatrix &elmat) override
virtual void AddMultTransposePA(const mfem::Vector &x, mfem::Vector &y) const override
void SetIntegrator(mfem::BilinearFormIntegrator *integ)
virtual void AddMultMF(const mfem::Vector &x, mfem::Vector &y) const override
ScaleIntegrator(mfem::BilinearFormIntegrator *integ, mfem::real_t scale)
virtual void AssembleElementMatrix2(const mfem::FiniteElement &trial_fe, const mfem::FiniteElement &test_fe, mfem::ElementTransformation &Trans, mfem::DenseMatrix &elmat) override