HDG Kernels
HDG kernels are an advanced systems that should only be developed by users with a fair amount of finite element experience. For background on hybridization, we encourage the user to read (Cockburn et al., 2009) which presents a unified framework for considering hybridization of discontinuous Galerkin, mixed, and continuous Galerkin methods for elliptic problems. (Cockburn et al., 2008) presents a single-face hybridizable discontinuous Galerkin (HDG) method for an elliptic problem, in which a non-zero stabilization term is added to only one face of a given element. (Nguyen et al., 2010) presents an HDG method for Stokes flow. (Nguyen et al., 2011) extends HDG to Navier-Stokes. More HDG literature may be found by looking at the research of Bernardo Cockburn, his former postdoc Sander Rhebergen, and Rhebergen's former postdoc Tamas Horvath. Work by Tan Bui-Thanh on upwind HDG methods, like in (Bui-Thanh, 2015) is also worth noting.
A hybridized finite element formulation starts with some primal finite element discretization. Then some continuity property of the finite element space is broken. For instance Raviart-Thomas finite elements may be used to solve a mixed formulation description of a Poisson problem. The Raviart-Thomas elements ensure continuity of the normal component of the vector field across element faces. We break that continuity in the finite element space used in the hybridized method and instead introduce degrees of freedom, that live only on the mesh skeleton (the faces of the mesh), that are responsible for ensuring the continuity that was lost by breaking the finite element space. In libMesh/MOOSE implementation terms, when hybridizing the Raviart-Thomas description of the Poisson problem, we change from using a RAVIART_THOMAS
basis to an L2_RAVIART_THOMAS
basis and introduce a SIDE_HIERARCHIC
variable whose degrees of freedom live on the mesh skeleton. We will refer to the variables that exist "before" the hybridization as primal variables and the variable(s) that live on the mesh skeleton as Lagrange multipliers (LMs) or dual variable(s).
We note that some classes of HDG methods, such as the LDG method in (Cockburn et al., 2008), have the gradient as an independent primal variable. With these methods, for diffusion or diffusion-dominated problems, the primal gradient and primal scalar variable fields can be used to postprocess a scalar field that converges with order in the norm, where is the polynomial order of the primal scalar variable. However, as advection becomes dominant, the higher order convergence is lost and consequently so is the value of having the gradient as an independent variable. In advection-dominated cases, interior penalty HDG methods, such as that outlined in (Rhebergen and Wells, 2017), may be a good choice.
Implementation in MOOSE
HDG kernels derive from Kernels. However, they add additional interfaces: computeResidualOnSide
and computeJacobianOnSide
which must be overridden and computeResidualAndJacobianOnSide
which may be optionally overridden if the HDG kernel developer wishes to enable the ability to compute the residual and Jacobian together. These interfaces will be called on internal faces on a per-element basis. This means that a given internal face will be visited twice, once from each element side. External boundary condition integration occurs with standard boundary condition classes, see BCs System.
There are currently two HDG implementations in MOOSE: L-HDG and IP-HDG. Both L-HDG and IP-HDG kernel classes inherit from HDGKernel
but that is where their similarity ends. L-HDG currently implements physics monolithically, e.g. the L-HDG discretization of the Navier-Stokes equations, both mass and momentum, is contained entirely within a single kernel. However, the MOOSE IP-HDG implementation is modular; like is typically the case in MOOSE, there is a single kernel object per PDE term. So for a 2D setup of the Navier-Stokes equations, there are five kernels. Two advection kernels for the x- and y-momentum component equations, two stress kernel (which contains both viscous stress and pressure) for the x- and y-momentum component equations, and one advection kernel for the mass equation.
This difference in kernel design naturally has consequences for boundary conditions. The monolithic kernel design for L-HDG leads to monolithic boundary conditions. Modularity for IP-HDG kernels means modular boundary conditions, e.g. a user may end up specifying multiple integrated boundary conditions on a single boundary like is done in
[BCs<<<{"href": "../BCs/index.html"}>>>]
[dirichlet_diff]
type = DiffusionIPHDGDirichletBC<<<{"description": "Weakly imposes Dirichlet boundary conditions for a hybridized discretization of a diffusion equation", "href": "../../source/bcs/DiffusionIPHDGDirichletBC.html"}>>>
functor<<<{"description": "The Dirichlet value for the primal variable. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = exact
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left right top bottom'
[]
[dirichlet_adv]
type = AdvectionIPHDGDirichletBC<<<{"description": "Weakly imposes Dirichlet boundary conditions for a hybridized discretization of an advection equation", "href": "../../source/bcs/AdvectionIPHDGDirichletBC.html"}>>>
functor<<<{"description": "The Dirichlet value for the primal variable. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = exact
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left right top bottom'
[]
[]
(moose/test/tests/hdgkernels/ip-advection-diffusion/mms-advection-diffusion.i)In the future we may make L-HDG design more modular, although the monolithic approach does have the advantage of less required user input. However, less user input may also be achieved in the future by leveraging the PhysicsPhysics system.
References
- Tan Bui-Thanh.
From godunov to a unified hybridized discontinuous galerkin framework for partial differential equations.
Journal of Computational Physics, 295:114–146, 2015.[BibTeX]
- Bernardo Cockburn, Bo Dong, and Johnny Guzmán.
A superconvergent ldg-hybridizable galerkin method for second-order elliptic problems.
Mathematics of Computation, 77(264):1887–1916, 2008.[BibTeX]
- Bernardo Cockburn, Jayadeep Gopalakrishnan, and Raytcho Lazarov.
Unified hybridization of discontinuous galerkin, mixed, and continuous galerkin methods for second order elliptic problems.
SIAM Journal on Numerical Analysis, 47(2):1319–1365, 2009.[BibTeX]
- Ngoc Cuong Nguyen, Jaime Peraire, and Bernardo Cockburn.
A hybridizable discontinuous galerkin method for stokes flow.
Computer Methods in Applied Mechanics and Engineering, 199(9-12):582–597, 2010.[BibTeX]
- Ngoc Cuong Nguyen, Jaume Peraire, and Bernardo Cockburn.
An implicit high-order hybridizable discontinuous galerkin method for the incompressible navier–stokes equations.
Journal of Computational Physics, 230(4):1147–1170, 2011.[BibTeX]
- Sander Rhebergen and Garth N Wells.
Analysis of a hybridized/interface stabilized finite element method for the stokes equations.
SIAM Journal on Numerical Analysis, 55(4):1982–2003, 2017.[BibTeX]