MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
lor_dg.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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#ifndef MFEM_LOR_DG
13#define MFEM_LOR_DG
14
15#include "lor_batched.hpp"
16
17namespace mfem
18{
19
20// BatchedLORKernel specialization for DG spaces. Not user facing. See the
21// classes BatchedLORAssembly and BatchedLORKernel .
23{
24 IntegrationRule ir_face; ///< Collocated Gauss-Lobatto face quadrature rule.
25 real_t kappa; ///< DG penalty parameter.
26public:
27 template <int ORDER, int SDIM> void Assemble2D();
28 template <int ORDER> void Assemble3D();
30 FiniteElementSpace &fes_ho_,
31 Vector &X_vert_,
32 Vector &sparse_ij_,
33 Array<int> &sparse_mapping_)
34 : BatchedLORKernel(fes_ho_, X_vert_, sparse_ij_, sparse_mapping_),
35 ir_face(GetLobattoIntRule(fes_ho_.GetMesh()->GetTypicalFaceGeometry(),
36 fes_ho_.GetMaxElementOrder() + 1))
37 {
40
41 auto *integ = GetInteriorFaceIntegrator<DGDiffusionIntegrator>(a);
42 if (integ)
43 {
44 kappa = integ->GetPenaltyParameter();
45 }
46 else
47 {
48 kappa = 0.0;
49 }
50 }
51
52 /// @brief Compute and return the face info array.
53 ///
54 /// The face info array has shape (6, nf), where @a nf is the number of
55 /// faces. For each face @a i, the column (:,i) has entries (e0, f0, o0, e1,
56 /// f1, o1), where @a e is adjacent element, @a f is the local face index,
57 /// and @a o is the orientation. For boundary and shared faces, (e1, f1, o1)
58 /// are all set to -1.
59 Array<int> GetFaceInfo() const;
60
61 /// @brief Compute and return the boundary penalty factor.
62 ///
63 /// The returned vector has shape (nq, nf), where @a nq is the number of
64 /// nodes per face, and @a nf is the number of faces.
65 ///
66 /// The boundary penalty factor is $J_f / h = J_f^2 / J_e$ (since $h = J_e /
67 /// J_f$), where $J_f$ is the face Jacobian determinant, and $J_e$ is the
68 /// element Jacobian determinant.
70
71 /// Assemble the face penalty terms in the matrix @a sparse_ij.
72 void AssembleFaceTerms();
73};
74
75}
76
77#include "lor_dg_impl.hpp"
78
79#endif
Abstract base class for the batched LOR assembly kernels.
CoefficientVector c2
Coefficient of second integrator.
CoefficientVector c1
Coefficient of first integrator.
BatchedLOR_DG(BilinearForm &a, FiniteElementSpace &fes_ho_, Vector &X_vert_, Vector &sparse_ij_, Array< int > &sparse_mapping_)
Definition lor_dg.hpp:29
Array< int > GetFaceInfo() const
Compute and return the face info array.
void AssembleFaceTerms()
Assemble the face penalty terms in the matrix sparse_ij.
Vector GetBdrPenaltyFactor() const
Compute and return the boundary penalty factor.
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:208
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Vector data type.
Definition vector.hpp:82
real_t kappa
Definition ex24.cpp:54
Mesh * GetMesh(int type)
Definition ex29.cpp:218
real_t a
Definition lissajous.cpp:41
IntegrationRule GetLobattoIntRule(Geometry::Type geom, int nd1d)
Return the Gauss-Lobatto rule for geometry geom with nd1d points per dimension.
void ProjectLORCoefficient(BilinearForm &a, CoefficientVector &coeff_vector)
float real_t
Definition config.hpp:46