https://mooseframework.inl.gov
ScaleIntegrator.C
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 MFEM_ENABLED
11 
12 #include "ScaleIntegrator.h"
13 
14 namespace Moose::MFEM
15 {
16 
17 void
18 ScaleIntegrator::SetIntRule(const mfem::IntegrationRule * ir)
19 {
20  IntRule = ir;
21  _integrator->SetIntRule(ir);
22 }
23 
24 void
25 ScaleIntegrator::AssembleElementMatrix(const mfem::FiniteElement & el,
26  mfem::ElementTransformation & Trans,
27  mfem::DenseMatrix & elmat)
28 {
30  _integrator->AssembleElementMatrix(el, Trans, elmat);
31  elmat *= _scale;
32 }
33 
34 void
35 ScaleIntegrator::AssembleElementMatrix2(const mfem::FiniteElement & el1,
36  const mfem::FiniteElement & el2,
37  mfem::ElementTransformation & Trans,
38  mfem::DenseMatrix & elmat)
39 {
41  _integrator->AssembleElementMatrix2(el1, el2, Trans, elmat);
42  elmat *= _scale;
43 }
44 
45 void
46 ScaleIntegrator::AssembleFaceMatrix(const mfem::FiniteElement & el1,
47  const mfem::FiniteElement & el2,
48  mfem::FaceElementTransformations & Trans,
49  mfem::DenseMatrix & elmat)
50 {
52  _integrator->AssembleFaceMatrix(el1, el2, Trans, elmat);
53  elmat *= _scale;
54 }
55 
56 void
57 ScaleIntegrator::AssemblePA(const mfem::FiniteElementSpace & fes)
58 {
60  _integrator->AssemblePA(fes);
61 }
62 
63 void
65 {
66  _integrator->AssembleDiagonalPA(diag);
67  diag *= _scale;
68 }
69 
70 void
71 ScaleIntegrator::AssemblePAInteriorFaces(const mfem::FiniteElementSpace & fes)
72 {
73  _integrator->AssemblePAInteriorFaces(fes);
74 }
75 
76 void
77 ScaleIntegrator::AssemblePABoundaryFaces(const mfem::FiniteElementSpace & fes)
78 {
79  _integrator->AssemblePABoundaryFaces(fes);
80 }
81 
82 void
83 ScaleIntegrator::AddMultPA(const mfem::Vector & x, mfem::Vector & y) const
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 }
92 
93 void
94 ScaleIntegrator::AddMultTransposePA(const mfem::Vector & x, mfem::Vector & y) const
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 }
103 
104 void
105 ScaleIntegrator::AssembleEA(const mfem::FiniteElementSpace & fes,
106  mfem::Vector & emat,
107  const bool add)
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 }
123 
125 {
126  if (_own_integrator)
127  delete _integrator;
128 }
129 
130 } // namespace Moose::MFEM
131 
132 #endif
virtual void SetIntRule(const mfem::IntegrationRule *ir)
mfem::BilinearFormIntegrator * _integrator
virtual void AssemblePAInteriorFaces(const mfem::FiniteElementSpace &fes)
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 void AddMultTransposePA(const mfem::Vector &x, mfem::Vector &y) const
virtual void AssemblePABoundaryFaces(const mfem::FiniteElementSpace &fes)
virtual void AssembleDiagonalPA(mfem::Vector &diag)
virtual void AssembleElementMatrix(const mfem::FiniteElement &el, 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 AssembleElementMatrix2(const mfem::FiniteElement &trial_fe, const mfem::FiniteElement &test_fe, mfem::ElementTransformation &Trans, mfem::DenseMatrix &elmat)
virtual void AssemblePA(const mfem::FiniteElementSpace &fes)