3.5 Scalar PDE on surfaces (with HDG)

We are solving the scalar linear transport problem

\[\partial_t u - \varepsilon \Delta_\Gamma u + \operatorname{div}_\Gamma ( \mathbf{w} u ) = 0 \quad \text{ on } \Gamma,\]

where \(\Gamma\) is a closed oriented surface in \(\mathbb{R}^3\).

In [1]:
import netgen.gui
%gui tk
from math import pi
from ngsolve import *

Mesh of (only) the surface:

We consider the unit sphere (only surface mesh!)

In [2]:
from netgen.csg import *
geo = CSGeometry()
geo.Add(Sphere(Pnt(0,0,0),1).bc("sphere"))
Out[2]:
0
In [3]:
from netgen.meshing import MeshingParameters
from netgen.meshing import MeshingStep
mp = MeshingParameters(maxh=0.5, perfstepsend = MeshingStep.MESHSURFACE)

mesh = Mesh(geo.GenerateMesh(mp=mp))
In [4]:
order = 4
if order > 0:
    mesh.Curve(order)
Draw(mesh)

Tangential velocity field:

In [5]:
b = CoefficientFunction((y,-x,0))
Draw (b, mesh, "b")
eps = 5e-5

some local quantities

normal, tangential, co-normal and mesh size

In [6]:
n = specialcf.normal(mesh.dim)
h = specialcf.mesh_size
local_tang = specialcf.tangential(mesh.dim)
con = Cross(n,local_tang) #co-normal pointing outside of a surface element
bn = InnerProduct(b,con)

Discretization

Surface FESpaces:

In [7]:
fesl2 = SurfaceL2(mesh, order=order)
fesedge = FacetSurface(mesh, order=order)

fes = FESpace([fesl2,fesedge])
u,uE = fes.TrialFunction()
v,vE = fes.TestFunction()

surface gradients, co-normal derivatives and jumps:

In [8]:
gradu = u.Trace().Deriv()
gradv = v.Trace().Deriv()
jumpu = u - uE
jumpv = v - vE
dudn = InnerProduct(gradu,con)
dvdn = InnerProduct(gradv,con)

Diffusion formulation (Hybrid Interior Penalty on Surface):

\[ \sum_T \int_T \nabla u \nabla v - \sum_T \int_{\partial T} n \nabla u (v-\widehat v) - \sum_T \int_{\partial T} n \nabla u (u-\widehat u) + \frac{\alpha p^2}{h} \sum_F \int_F (u-\widehat u)(v-\widehat v)\]
In [9]:
a = BilinearForm(fes)
# diffusion
alpha = 10*(order+1)**2
a += SymbolicBFI(eps*gradu * gradv, VOL_or_BND = BND)
a += SymbolicBFI(- eps*dudn * jumpv,element_boundary=True, VOL_or_BND = BND)
a += SymbolicBFI(- eps*dvdn * jumpu,element_boundary=True, VOL_or_BND = BND)
a += SymbolicBFI(eps*alpha/h * jumpu * jumpv,element_boundary=True, VOL_or_BND = BND)

Convection formulation (Hybrid Upwinding):

\[ \sum_T \int_T - u \mathbf{w} \cdot \nabla v + \sum_T \int_{\partial T} (\mathbf{w} \cdot \mathbf{n}) u^{hyb-upw} v + \sum_T \int_{\partial T_{out}} (\mathbf{w} \cdot \mathbf{n}) (\hat{u} - u) \hat{v}\]

with

\[\begin{split}u^{hyb-upw} = \left\{ \begin{array}{cll} u & \text{if } \mathbf{w} \cdot \mathbf{n} & (outflow), \\ \hat{u} & \text{else} & (inflow) . \end{array} \right.\end{split}\]
alt alt
In [10]:

# convection (surface version of Egger+Schöberl formulation:
a += SymbolicBFI(- u * InnerProduct(b,gradv), VOL_or_BND = BND)
a += SymbolicBFI(IfPos(bn, bn*u, bn*uE) * v, element_boundary = True, VOL_or_BND = BND)
a += SymbolicBFI(IfPos(bn, bn, 0) * (uE - u) * vE, element_boundary = True, VOL_or_BND = BND)

a.Assemble()

Mass matrix:

In [11]:
m = BilinearForm(fes)
m += SymbolicBFI(u*v, VOL_or_BND = BND)
m.Assemble()

\(M^\ast\)-matrix for implicit Euler

In [12]:
dt = 0.02

mstar = m.mat.CreateMatrix()
mstar.AsVector().data = m.mat.AsVector() + dt * a.mat.AsVector()
mstarinv = mstar.Inverse()

from math import pi
T = 8*pi

Time stepping

Visualization of solution (Deformation in normal direction)

In [13]:
gfu = GridFunction(fes)
gfu.components[0].Set(1.5*exp(-20*(x*x+(y-1)**2+(z)**2)),definedon=mesh.Boundaries("sphere"))
Draw (gfu.components[0], mesh, "u")
Draw (gfu.components[0]*n, mesh, "un", sd=4, autoscale=False)
res = gfu.vec.CreateVector()
SetVisualization(deformation=True)

The time loop

In [14]:
for i in range(int(T/dt)):
    res.data = m.mat * gfu.vec
    gfu.vec.data = mstarinv * res
    Redraw(blocking=True)

Remark:

  • We can apply static condensation as in the volume