MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilininteg_br2.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "bilininteg.hpp"
13 #include "pfespace.hpp"
14 #include <algorithm>
15 
16 namespace mfem
17 {
18 
20  FiniteElementSpace &fes, double e) : eta(e), Q(NULL)
21 {
23 }
24 
26  FiniteElementSpace &fes, Coefficient &Q_, double e) : eta(e), Q(&Q_)
27 {
29 }
30 
32  FiniteElementSpace *fes, double e) : eta(e), Q(NULL)
33 {
35 }
36 
38 {
39  MFEM_VERIFY(fes.IsDGSpace(),
40  "The BR2 integrator is only defined for DG spaces.");
41  // Precompute local mass matrix inverses needed for the lifting operators
42  // First compute offsets and total size needed (e.g. for mixed meshes or
43  // p-refinement)
44  int nel = fes.GetNE();
45  Minv_offsets.SetSize(nel+1);
46  ipiv_offsets.SetSize(nel+1);
47  ipiv_offsets[0] = 0;
48  Minv_offsets[0] = 0;
49  for (int i=0; i<nel; ++i)
50  {
51  int dof = fes.GetFE(i)->GetDof();
52  ipiv_offsets[i+1] = ipiv_offsets[i] + dof;
53  Minv_offsets[i+1] = Minv_offsets[i] + dof*dof;
54  }
55 
56 #ifdef MFEM_USE_MPI
57  // When running in parallel, we also need to compute the local mass matrices
58  // of face neighbor elements
59  ParFiniteElementSpace *pfes = dynamic_cast<ParFiniteElementSpace *>(&fes);
60  if (pfes != NULL)
61  {
62  ParMesh *pmesh = pfes->GetParMesh();
63  pfes->ExchangeFaceNbrData();
64  int nel_nbr = pmesh->GetNFaceNeighborElements();
65  Minv_offsets.SetSize(nel+nel_nbr+1);
66  ipiv_offsets.SetSize(nel+nel_nbr+1);
67  for (int i=0; i<nel_nbr; ++i)
68  {
69  int dof = pfes->GetFaceNbrFE(i)->GetDof();
70  ipiv_offsets[nel+i+1] = ipiv_offsets[nel+i] + dof;
71  Minv_offsets[nel+i+1] = Minv_offsets[nel+i] + dof*dof;
72  }
73  nel += nel_nbr;
74  }
75 #endif
76  // The final "offset" is the total size of all the blocks
79 
80  // Assemble the local mass matrices and compute LU factorization
81  MassIntegrator mi;
82  for (int i=0; i<nel; ++i)
83  {
84  const FiniteElement *fe = NULL;
85  ElementTransformation *tr = NULL;
86  if (i < fes.GetNE())
87  {
88  fe = fes.GetFE(i);
89  tr = fes.GetElementTransformation(i);
90  }
91  else
92  {
93 #ifdef MFEM_USE_MPI
94  int inbr = i - fes.GetNE();
95  fe = pfes->GetFaceNbrFE(inbr);
96  tr = pfes->GetParMesh()->GetFaceNbrElementTransformation(inbr);
97 #endif
98  }
99  int dof = fe->GetDof();
100  double *Minv_el = &Minv[Minv_offsets[i]];
101  int *ipiv_el = &ipiv[ipiv_offsets[i]];
102  DenseMatrix Me(Minv_el, dof, dof);
103  mi.AssembleElementMatrix(*fe, *tr, Me);
104  LUFactors lu(Minv_el, ipiv_el);
105  lu.Factor(dof);
106  }
107 }
108 
110  const FiniteElement &el1, const FiniteElement &el2,
112 {
113  int ndof1 = el1.GetDof();
114  shape1.SetSize(ndof1);
115 
116  R11.SetSize(ndof1, ndof1);
117  R11 = 0.0;
118  LUFactors M1inv(&Minv[Minv_offsets[Trans.Elem1No]],
119  &ipiv[ipiv_offsets[Trans.Elem1No]]);
120  LUFactors M2inv;
121 
122  double factor = Geometries.NumBdr(Trans.Elem1->GetGeometryType());
123 
124  int ndof2;
125  if (Trans.Elem2No >= 0)
126  {
127  ndof2 = el2.GetDof();
128  shape2.SetSize(ndof2);
129  R12.SetSize(ndof1, ndof2);
130  R21.SetSize(ndof2, ndof1);
131  R22.SetSize(ndof2, ndof2);
132  M2inv.data = &Minv[Minv_offsets[Trans.Elem2No]];
133  M2inv.ipiv = &ipiv[ipiv_offsets[Trans.Elem2No]];
134 
135  R12 = 0.0;
136  R21 = 0.0;
137  R22 = 0.0;
138 
139  Geometry::Type geom2 = Trans.Elem2->GetGeometryType();
140  factor = std::max(factor, double(Geometries.NumBdr(geom2)));
141  }
142  else
143  {
144  ndof2 = 0;
145  }
146 
147  int ndofs = ndof1 + ndof2;
148 
149  Re.SetSize(ndofs, ndofs);
150  MinvRe.SetSize(ndofs, ndofs);
151 
152  elmat.SetSize(ndofs);
153  elmat = 0.0;
154 
155  const IntegrationRule *ir = IntRule;
156  if (ir == NULL)
157  {
158  int order;
159  if (ndof2)
160  {
161  order = 2*std::max(el1.GetOrder(), el2.GetOrder());
162  }
163  else
164  {
165  order = 2*el1.GetOrder();
166  }
167  ir = &IntRules.Get(Trans.FaceGeom, order);
168  }
169 
170  for (int p = 0; p < ir->GetNPoints(); p++)
171  {
172  const IntegrationPoint &ip = ir->IntPoint(p);
173  Trans.SetAllIntPoints(&ip);
174 
175  const IntegrationPoint &eip1 = Trans.Elem1->GetIntPoint();
176  el1.CalcShape(eip1, shape1);
177  double q = Q ? Q->Eval(*Trans.Elem1, eip1) : 1.0;
178  if (ndof2)
179  {
180  const IntegrationPoint &eip2 = Trans.Elem2->GetIntPoint();
181  el2.CalcShape(eip2, shape2);
182  // Set coefficient value q to the average of the values on either side
183  if (Q) { q = 0.5*(q + Q->Eval(*Trans.Elem2, eip2)); }
184  }
185  // Take sqrt here because
186  // eta (r_e([u]), r_e([v])) = (sqrt(eta) r_e([u]), sqrt(eta) r_e([v]))
187  double w = sqrt((factor + 1)*eta*q)*ip.weight*Trans.Face->Weight();
188  // r_e is defined by, (r_e([u]), tau) = <[u], {tau}>, so we pick up a
189  // factor of 0.5 on interior faces from the average term.
190  if (ndof2) { w *= 0.5; }
191 
192  for (int i = 0; i < ndof1; i++)
193  {
194  const double wsi = w*shape1(i);
195  for (int j = 0; j < ndof1; j++)
196  {
197  R11(i, j) += wsi*shape1(j);
198  }
199  }
200 
201  if (ndof2)
202  {
203  for (int i = 0; i < ndof2; i++)
204  {
205  const double wsi = w*shape2(i);
206  for (int j = 0; j < ndof1; j++)
207  {
208  R21(i, j) += wsi*shape1(j);
209  R12(j, i) -= wsi*shape1(j);
210  }
211  for (int j = 0; j < ndof2; j++)
212  {
213  R22(i, j) -= wsi*shape2(j);
214  }
215  }
216  }
217  }
218 
219  MinvR11 = R11;
220  M1inv.Solve(ndof1, ndof1, MinvR11.Data());
221  for (int i = 0; i < ndof1; i++)
222  {
223  for (int j = 0; j < ndof1; j++)
224  {
225  Re(i, j) = R11(i, j);
226  MinvRe(i, j) = MinvR11(i, j);
227  }
228  }
229 
230  if (ndof2)
231  {
232  MinvR12 = R12;
233  MinvR21 = R21;
234  MinvR22 = R22;
235  M1inv.Solve(ndof1, ndof2, MinvR12.Data());
236  M2inv.Solve(ndof2, ndof1, MinvR21.Data());
237  M2inv.Solve(ndof2, ndof2, MinvR22.Data());
238 
239  for (int i = 0; i < ndof2; i++)
240  {
241  for (int j = 0; j < ndof1; j++)
242  {
243  Re(ndof1 + i, j) = R21(i, j);
244  MinvRe(ndof1 + i, j) = MinvR21(i, j);
245 
246  Re(j, ndof1 + i) = R12(j, i);
247  MinvRe(j, ndof1 + i) = MinvR12(j, i);
248  }
249  for (int j = 0; j < ndof2; j++)
250  {
251  Re(ndof1 + i, ndof1 + j) = R22(i, j);
252  MinvRe(ndof1 + i, ndof1 + j) = MinvR22(i, j);
253  }
254  }
255  }
256 
257  // Compute the matrix associated with (r_e([u]), r_e([u])).
258  // The matrix for r_e([u]) is `MinvRe`, and so we need to form the product
259  // `(MinvRe)^T M MinvRe`. Using `Minv^T M = Minv M = I`, we obtain
260  // `Re^T MinvRe`.
261  MultAtB(Re, MinvRe, elmat);
262 }
263 
264 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
ParMesh * GetParMesh() const
Definition: pfespace.hpp:277
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition: eltrans.hpp:162
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:327
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Abstract parallel finite element space.
Definition: pfespace.hpp:28
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition: eltrans.hpp:98
DGDiffusionBR2Integrator(class FiniteElementSpace &fes, double e=1.0)
ElementTransformation * Elem2
Definition: eltrans.hpp:522
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:614
Geometry Geometries
Definition: fe.cpp:49
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1517
void SetAllIntPoints(const IntegrationPoint *face_ip)
Set the integration point in the Face and the two neighboring elements, if present.
Definition: eltrans.hpp:566
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Definition: fespace.hpp:925
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
virtual void Solve(int m, int n, double *X) const
Definition: densemat.cpp:3381
double * data
Definition: densemat.hpp:599
int GetNFaceNeighborElements() const
Definition: pmesh.hpp:463
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:647
double * Data() const
Returns the matrix data array.
Definition: densemat.hpp:111
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
void PrecomputeMassInverse(class FiniteElementSpace &fes)
Precomputes the inverses (LU factorizations) of the local mass matrices.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
int NumBdr(int GeomType)
Return the number of boundary &quot;faces&quot; of a given Geometry::Type.
Definition: geom.hpp:126
const IntegrationRule * IntRule
Definition: nonlininteg.hpp:30
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: pmesh.cpp:1976
Class for integration point with weight.
Definition: intrules.hpp:25
virtual bool Factor(int m, double TOL=0.0)
Compute the LU factorization of the current matrix.
Definition: densemat.cpp:3230
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2781
ElementTransformation * Elem1
Definition: eltrans.hpp:522
ElementTransformation * Face
Definition: eltrans.hpp:523
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
Class for parallel meshes.
Definition: pmesh.hpp:32
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379