This page was generated from jupyter-files/unit-2.8-DG/DG.ipynb.
2.8 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 \(-\Delta u\):
with jump-term over facets:
and averaging operator
DG form for \(\Div (b u)\), where \(b\) is the given wind:
[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:
[2]:
order=4
fes = L2(mesh, order=order, dgjumps=True)
u,v = fes.TnT()
Every facet has a master element. The value from the other element is referred to via the Other()
operator:
[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:
[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()
[5]:
f = LinearForm(fes)
f += SymbolicLFI(1*v)
f.Assemble()
[6]:
gfu = GridFunction(fes, name="uDG")
gfu.vec.data = a.mat.Inverse() * f.vec
Draw (gfu)
DG requires a lot of additional matrix entries:
[7]:
print ("a nze:", a.mat.nze)
fes2 = L2(mesh, order=order)
a2 = BilinearForm(fes2)
a2 += SymbolicBFI(u*v)
a2.Assemble()
print ("a2 nze:", a2.mat.nze)
a nze: 23400
a2 nze: 6750
Next we are solving a convection-diffusion problem:
[8]:
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:
[9]:
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)
[10]:
acd.Assemble()
[11]:
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:
the jump-term is now replaced by the difference \(u - \widehat u\).
No additional matrix entries across elements are produced. Dirichlet boundary conditions are set as usual to the facet variable:
[12]:
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:
[13]:
jump_u = u-uhat
jump_v = v-vhat
[14]:
alpha = 4
condense = True
h = specialcf.mesh_size
n = specialcf.normal(mesh.dim)
a = BilinearForm(fes, eliminate_internal=condense)
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()
gfu = GridFunction(fes)
print ("A non-zero elements:", a.mat.nze)
A non-zero elements: 5825
[15]:
if not condense:
inv = a.mat.Inverse(fes.FreeDofs(), "umfpack")
gfu.vec.data = inv * f.vec
else:
f.vec.data += a.harmonic_extension_trans * f.vec
inv = a.mat.Inverse(fes.FreeDofs(True), "umfpack")
gfu.vec.data = inv * f.vec
gfu.vec.data += a.harmonic_extension * gfu.vec
gfu.vec.data += a.inner_solve * f.vec
Draw (gfu.components[0], mesh, "u-HDG")
[ ]: