https://mooseframework.inl.gov
Public Member Functions | Private Attributes | List of all members
Moose::MFEM::ScaleIntegrator Class Reference

Integrator which scales its results by a constant value. More...

#include <ScaleIntegrator.h>

Inheritance diagram for Moose::MFEM::ScaleIntegrator:
[legend]

Public Member Functions

 ScaleIntegrator (mfem::BilinearFormIntegrator *integ)
 
 ScaleIntegrator (mfem::BilinearFormIntegrator *integ, double scale)
 
 ScaleIntegrator (mfem::BilinearFormIntegrator *integ, double scale, bool own)
 
void SetIntegrator (mfem::BilinearFormIntegrator *integ)
 
void SetScale (double scale)
 
void SetOwn (bool own)
 
void CheckIntegrator ()
 
virtual void SetIntRule (const mfem::IntegrationRule *ir)
 
virtual void AssembleElementMatrix (const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, mfem::DenseMatrix &elmat)
 
virtual void AssembleElementMatrix2 (const mfem::FiniteElement &trial_fe, const mfem::FiniteElement &test_fe, mfem::ElementTransformation &Trans, mfem::DenseMatrix &elmat)
 
virtual void AssembleFaceMatrix (const mfem::FiniteElement &el1, const mfem::FiniteElement &el2, mfem::FaceElementTransformations &Trans, mfem::DenseMatrix &elmat)
 
virtual void AssemblePA (const mfem::FiniteElementSpace &fes)
 
virtual void AssembleDiagonalPA (mfem::Vector &diag)
 
virtual void AssemblePAInteriorFaces (const mfem::FiniteElementSpace &fes)
 
virtual void AssemblePABoundaryFaces (const mfem::FiniteElementSpace &fes)
 
virtual void AddMultTransposePA (const mfem::Vector &x, mfem::Vector &y) const
 
virtual void AddMultPA (const mfem::Vector &x, mfem::Vector &y) const
 
virtual void AssembleEA (const mfem::FiniteElementSpace &fes, mfem::Vector &emat, const bool add)
 
virtual ~ScaleIntegrator ()
 

Private Attributes

mfem::BilinearFormIntegrator * _integrator {nullptr}
 
double _scale
 
bool _own_integrator
 

Detailed Description

Integrator which scales its results by a constant value.

Definition at line 22 of file ScaleIntegrator.h.

Constructor & Destructor Documentation

◆ ScaleIntegrator() [1/3]

Moose::MFEM::ScaleIntegrator::ScaleIntegrator ( mfem::BilinearFormIntegrator *  integ)
inline

Definition at line 30 of file ScaleIntegrator.h.

31  : _integrator{integ}, _scale{1}, _own_integrator{true}
32  {
33  }
mfem::BilinearFormIntegrator * _integrator

◆ ScaleIntegrator() [2/3]

Moose::MFEM::ScaleIntegrator::ScaleIntegrator ( mfem::BilinearFormIntegrator *  integ,
double  scale 
)
inline

Definition at line 34 of file ScaleIntegrator.h.

35  : _integrator{integ}, _scale{scale}, _own_integrator{true}
36  {
37  }
void scale(MeshBase &mesh, const Real xs, const Real ys=0., const Real zs=0.)
mfem::BilinearFormIntegrator * _integrator

◆ ScaleIntegrator() [3/3]

Moose::MFEM::ScaleIntegrator::ScaleIntegrator ( mfem::BilinearFormIntegrator *  integ,
double  scale,
bool  own 
)
inline

Definition at line 38 of file ScaleIntegrator.h.

39  : _integrator{integ}, _scale{scale}, _own_integrator{own}
40  {
41  }
void scale(MeshBase &mesh, const Real xs, const Real ys=0., const Real zs=0.)
mfem::BilinearFormIntegrator * _integrator

◆ ~ScaleIntegrator()

Moose::MFEM::ScaleIntegrator::~ScaleIntegrator ( )
virtual

Definition at line 124 of file ScaleIntegrator.C.

125 {
126  if (_own_integrator)
127  delete _integrator;
128 }
mfem::BilinearFormIntegrator * _integrator

Member Function Documentation

◆ AddMultPA()

void Moose::MFEM::ScaleIntegrator::AddMultPA ( const mfem::Vector &  x,
mfem::Vector &  y 
) const
virtual

Definition at line 83 of file ScaleIntegrator.C.

84 {
85  // y += Mx*scale
86  mfem::Vector Mx(y.Size());
87  Mx = 0.0;
88  _integrator->AddMultPA(x, Mx);
89  Mx *= _scale;
90  y += Mx;
91 }
mfem::BilinearFormIntegrator * _integrator

◆ AddMultTransposePA()

void Moose::MFEM::ScaleIntegrator::AddMultTransposePA ( const mfem::Vector &  x,
mfem::Vector &  y 
) const
virtual

Definition at line 94 of file ScaleIntegrator.C.

95 {
96  // y += M^T x*scale
97  mfem::Vector MTx(y.Size());
98  MTx = 0.0;
99  _integrator->AddMultTransposePA(x, MTx);
100  MTx *= _scale;
101  y += MTx;
102 }
mfem::BilinearFormIntegrator * _integrator

◆ AssembleDiagonalPA()

void Moose::MFEM::ScaleIntegrator::AssembleDiagonalPA ( mfem::Vector &  diag)
virtual

Definition at line 64 of file ScaleIntegrator.C.

65 {
66  _integrator->AssembleDiagonalPA(diag);
67  diag *= _scale;
68 }
mfem::BilinearFormIntegrator * _integrator

◆ AssembleEA()

void Moose::MFEM::ScaleIntegrator::AssembleEA ( const mfem::FiniteElementSpace &  fes,
mfem::Vector &  emat,
const bool  add 
)
virtual

Definition at line 105 of file ScaleIntegrator.C.

108 {
109  CheckIntegrator();
110  if (add)
111  {
112  mfem::Vector emat_scale(emat.Size());
113  _integrator->AssembleEA(fes, emat_scale, false);
114  emat_scale *= _scale;
115  emat += emat_scale;
116  }
117  else
118  {
119  _integrator->AssembleEA(fes, emat, add);
120  emat *= _scale;
121  }
122 }
mfem::BilinearFormIntegrator * _integrator

◆ AssembleElementMatrix()

void Moose::MFEM::ScaleIntegrator::AssembleElementMatrix ( const mfem::FiniteElement &  el,
mfem::ElementTransformation &  Trans,
mfem::DenseMatrix &  elmat 
)
virtual

Definition at line 25 of file ScaleIntegrator.C.

28 {
30  _integrator->AssembleElementMatrix(el, Trans, elmat);
31  elmat *= _scale;
32 }
mfem::BilinearFormIntegrator * _integrator

◆ AssembleElementMatrix2()

void Moose::MFEM::ScaleIntegrator::AssembleElementMatrix2 ( const mfem::FiniteElement &  trial_fe,
const mfem::FiniteElement &  test_fe,
mfem::ElementTransformation &  Trans,
mfem::DenseMatrix &  elmat 
)
virtual

Definition at line 35 of file ScaleIntegrator.C.

39 {
41  _integrator->AssembleElementMatrix2(el1, el2, Trans, elmat);
42  elmat *= _scale;
43 }
mfem::BilinearFormIntegrator * _integrator

◆ AssembleFaceMatrix()

void Moose::MFEM::ScaleIntegrator::AssembleFaceMatrix ( const mfem::FiniteElement &  el1,
const mfem::FiniteElement &  el2,
mfem::FaceElementTransformations &  Trans,
mfem::DenseMatrix &  elmat 
)
virtual

Definition at line 46 of file ScaleIntegrator.C.

50 {
52  _integrator->AssembleFaceMatrix(el1, el2, Trans, elmat);
53  elmat *= _scale;
54 }
mfem::BilinearFormIntegrator * _integrator

◆ AssemblePA()

void Moose::MFEM::ScaleIntegrator::AssemblePA ( const mfem::FiniteElementSpace &  fes)
virtual

Definition at line 57 of file ScaleIntegrator.C.

58 {
60  _integrator->AssemblePA(fes);
61 }
mfem::BilinearFormIntegrator * _integrator

◆ AssemblePABoundaryFaces()

void Moose::MFEM::ScaleIntegrator::AssemblePABoundaryFaces ( const mfem::FiniteElementSpace &  fes)
virtual

Definition at line 77 of file ScaleIntegrator.C.

78 {
79  _integrator->AssemblePABoundaryFaces(fes);
80 }
mfem::BilinearFormIntegrator * _integrator

◆ AssemblePAInteriorFaces()

void Moose::MFEM::ScaleIntegrator::AssemblePAInteriorFaces ( const mfem::FiniteElementSpace &  fes)
virtual

Definition at line 71 of file ScaleIntegrator.C.

72 {
73  _integrator->AssemblePAInteriorFaces(fes);
74 }
mfem::BilinearFormIntegrator * _integrator

◆ CheckIntegrator()

void Moose::MFEM::ScaleIntegrator::CheckIntegrator ( )
inline

Definition at line 57 of file ScaleIntegrator.h.

Referenced by AssembleEA(), AssembleElementMatrix(), AssembleElementMatrix2(), AssembleFaceMatrix(), and AssemblePA().

58  {
59  if (!_integrator)
60  mooseError("Integrator not set");
61  }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
mfem::BilinearFormIntegrator * _integrator

◆ SetIntegrator()

void Moose::MFEM::ScaleIntegrator::SetIntegrator ( mfem::BilinearFormIntegrator *  integ)
inline

Definition at line 43 of file ScaleIntegrator.h.

44  {
46  {
47  delete _integrator;
48  }
49 
50  _integrator = integ;
51  }
mfem::BilinearFormIntegrator * _integrator

◆ SetIntRule()

void Moose::MFEM::ScaleIntegrator::SetIntRule ( const mfem::IntegrationRule *  ir)
virtual

Definition at line 18 of file ScaleIntegrator.C.

19 {
20  IntRule = ir;
21  _integrator->SetIntRule(ir);
22 }
mfem::BilinearFormIntegrator * _integrator

◆ SetOwn()

void Moose::MFEM::ScaleIntegrator::SetOwn ( bool  own)
inline

Definition at line 55 of file ScaleIntegrator.h.

◆ SetScale()

void Moose::MFEM::ScaleIntegrator::SetScale ( double  scale)
inline

Definition at line 53 of file ScaleIntegrator.h.

53 { _scale = scale; }
void scale(MeshBase &mesh, const Real xs, const Real ys=0., const Real zs=0.)

Member Data Documentation

◆ _integrator

mfem::BilinearFormIntegrator* Moose::MFEM::ScaleIntegrator::_integrator {nullptr}
private

◆ _own_integrator

bool Moose::MFEM::ScaleIntegrator::_own_integrator
private

Definition at line 27 of file ScaleIntegrator.h.

Referenced by SetIntegrator(), SetOwn(), and ~ScaleIntegrator().

◆ _scale

double Moose::MFEM::ScaleIntegrator::_scale
private

The documentation for this class was generated from the following files: