Processing math: 100%

Discontinuous Galerkin Methods

  • Use discontinuous finite element spaces to solve PDEs.
  • Allows upwind-stabilization for convection-dominated problems
  • Requires additional jump terms for consistency

Interior penalty DG form for Δu:

A(u,v)=TTuv12FF{nu}[v]12FF{nu}[v]+αp2hFF[u][v]

with jump-term over facets: [u]=ulefturight

and averaging operator {nu}=12(nleftuleft+nlefturight)

DG form for divbu, where b is the given wind:

B(u,v)=Tbuv+FFbnuupwindv
In [1]:
import netgen.gui
%gui tk
from netgen.geom2d import unit_square
from ngsolve import *
mesh = Mesh(unit_square.GenerateMesh(maxh=0.3))

The space is responsible for allocating the matrix graph. Tell it that it should reserve entries for the coupling terms:

In [2]:
order=4
fes = L2(mesh, order=order, flags = { "dgjumps" : True })
u = fes.TrialFunction()
v = fes.TestFunction()

Every facet has a master element. The value from the other element is referred to via the Other() operator:

In [3]:
jump_u = u-u.Other()
jump_v = v-v.Other()
n = specialcf.normal(2)
mean_dudn = 0.5*n * (grad(u)+grad(u.Other()))
mean_dvdn = 0.5*n * (grad(v)+grad(v.Other()))

Integrals on facets are computed by setting skeleton=True. This iterates over all internal facets. Additionally setting BND iterates only over boundary facets:

In [4]:
alpha = 4
h = specialcf.mesh_size
a = BilinearForm(fes)
a += SymbolicBFI(grad(u)*grad(v))
a += SymbolicBFI(alpha*order**2/h*jump_u*jump_v, skeleton=True)
a += SymbolicBFI(alpha*order**2/h*u*v, BND, skeleton=True)
a += SymbolicBFI(-mean_dudn*jump_v -mean_dvdn*jump_u, skeleton=True)
a += SymbolicBFI(-n*grad(u)*v-n*grad(v)*u, BND, skeleton=True)
a.Assemble()
In [5]:
f = LinearForm(fes)
f += SymbolicLFI(1*v)
f.Assemble()
In [6]:
gfu = GridFunction(fes, name="uDG")
gfu.vec.data = a.mat.Inverse() * f.vec
Draw (gfu)

Next we are solving a convection-diffusion problem:

In [7]:
alpha = 4
h = specialcf.mesh_size
acd = BilinearForm(fes)
acd += SymbolicBFI(grad(u)*grad(v))
acd += SymbolicBFI(alpha*order**2/h*jump_u*jump_v, skeleton=True)
acd += SymbolicBFI(alpha*order**2/h*u*v, BND, skeleton=True)
acd += SymbolicBFI(-mean_dudn*jump_v -mean_dvdn*jump_u, skeleton=True)
acd += SymbolicBFI(-n*grad(u)*v-n*grad(v)*u, BND, skeleton=True)

The IfPos checks whether the first argument is positive. Then it returns the second one, else the third one. This is used to define the upwind flux. The check is performed in every integration-point on the skeleton:

In [8]:
b = CoefficientFunction( (20,1) )
acd += SymbolicBFI(-b * u * grad(v))
uup = IfPos(b*n, u, u.Other())
acd += SymbolicBFI(b*n*uup*jump_v, skeleton=True)
In [9]:
acd.Assemble()
In [10]:
gfu = GridFunction(fes)
gfu.vec.data = acd.mat.Inverse(freedofs=fes.FreeDofs(),inverse="umfpack") * f.vec
Draw (gfu)

Hybrid Discontinuous Galerkin methods

use additionally the hybrid facet variable on the skeleton:

A(u,ˆu;v,ˆv)=TTuvTTnu(vˆv)TTnu(uˆu)+αp2hFF(uˆu)(vˆv)

the jump-term is now replaced by the difference uˆu.

No additional matrix entries across elements are produced. Dirichlet boundary conditions are set as usual to the facet variable:

In [11]:
order=4
V = L2(mesh, order=order)
F = FacetFESpace(mesh, order=order, dirichlet="bottom|left|right|top")
fes = FESpace([V,F])
u,uhat = fes.TrialFunction()
v,vhat = fes.TestFunction()

Now, the jump is the difference between element-term and facet-term:

In [12]:
jump_u = u-uhat
jump_v = v-vhat
In [13]:
alpha = 4
h = specialcf.mesh_size
n = specialcf.normal(mesh.dim)

a = BilinearForm(fes)
a += SymbolicBFI(grad(u)*grad(v))
a += SymbolicBFI(alpha*order**2/h*jump_u*jump_v, element_boundary=True)
a += SymbolicBFI(-grad(u)*n*jump_v - grad(v)*n*jump_u, element_boundary=True)


b = CoefficientFunction( (20,1) )
a += SymbolicBFI(-b * u * grad(v))
uup = IfPos(b*n, u, uhat)
a += SymbolicBFI(b*n*uup*jump_v, element_boundary=True)
a.Assemble()

f = LinearForm(fes)
f += SymbolicLFI(1*v)
f.Assemble()
In [14]:
gfu = GridFunction(fes)
inv = a.mat.Inverse(fes.FreeDofs(), "umfpack")
gfu.vec.data = inv * f.vec

Draw (gfu.components[0], mesh, "u-HDG")