4.2. Mixed formulation for second order equations#
We have used \(H(div)\) elements in the adaptivity tutorial. Here is used to introduce a mixed formulation for second order equations.
Why mixed formulations?:
exact flux conservation
useful for a posteriori error estimates
model problem for \(4^{th}\) order problems, Stokes, …
We consider the diffusion equation $$
\begin{array}{rcll} -\text{div} \lambda \nabla u & = & f & \text{ in } \Omega \ u & = & u_D & \text{ on } \Gamma_D \ \lambda \frac{\partial u}{\partial n} & = & g & \text{ on } \Gamma_N \end{array} $$
Primal variational formulation
Find \(u \in H^1, u = u_D\) on \(\Gamma_D\) such that
We can recast the fisrt equation as a set of two equations introducing the flux \(\sigma = \lambda \nabla u\):
Find scalar \(u\) and the flux \(\sigma\) such that
with boundary conditions
Mixed variational formulation
Find \(\sigma \in H(\text{div} )\) and \(u \in L_2\) such that \(\sigma_n = g\) on \(\Gamma_N\) and
\begin{array}{rlll} \int_\Omega \lambda^{-1} \sigma \tau &+ \int_\Omega \text{div} \tau u & = & 0 \[0.4 cm] \int_\Omega \text{div} \sigma v & & = & -\int_\Omega f v + \int_{\Gamma_D} u_D \tau_n \end{array}
for all test-functions \(\tau \in H(\text{div} )\) and \(v \in L_2\) with \(\tau_n = 0\).
Here \(\sigma_n\) is the normal trace operator \(\sigma \cdot n |_{\Gamma_N}\), which is meaningful in \(H(\text{div} )\).
The big-Bilinear Form
A Compact notation (and meaningful) notation in the context of FEM is the single-liner
Find \((\sigma, u) \in H(\text{div} ) \times L_2\) such that \(\sigma_n = g\) on \(\Gamma_N\) and
for all test-functions \((\tau, v) \in H(\text{div} ) \times L_2\) with \(\tau_n = 0\).
the above can be seen as
\begin{array}{rlll} a((\sigma, u), (\tau, v)) & = & L((\tau, v)) \end{array}
import netgen.occ
from ngsolve import *
from ngsolve.webgui import Draw
from netgen.occ import *
def MakeGeometryOCC():
# create a rectangle with boundary conditions
square = Circle((0, 0), 2).Face()
square.edges.name = "outer"
square.faces.name = "air"
# create a rectangle with boundary conditions
el1 = MoveTo(-0.4, 0.2).Rectangle(0.8, 0.1).Face()
el1.edges.name = "el1"
el1.vertices.name = "el1"
# create a rectangle with boundary conditions
el2 = MoveTo(-0.4, -0.2).Rectangle(0.8, 0.1).Face()
el2.edges.name = "el2"
el2.vertices.name = "el2"
dielec = MoveTo(-0.9, -0.1).Rectangle(1.8, 0.3).Face()
dielec.faces.name = "dielec"
air = square # subtract the rectangles from the air rectangle
shape = Glue([air - dielec, dielec])
shape = shape - el1 - el2
# adding extra specifications on the shape
# predefined mesh size for the shape
Draw(shape)
return OCCGeometry(shape, dim=2)
mesh = Mesh(MakeGeometryOCC().GenerateMesh(maxh=0.1))
mesh.Curve(3)
<ngsolve.comp.Mesh at 0x7f151eaa1130>
lam = CoefficientFunction([1, 2])
Setup and solve primal problem:
4.2.1. Solving the mixed problem#
Define spaces for mixed problem. Note that essential boundary conditions are set to the \(H(\text{div})\)-component on the opposite boundary. Creating a space from a list of spaces generates a product space:
order_flux=1
V = HDiv(mesh, order=order_flux, dirichlet = "outer")
Q = L2(mesh, order=order_flux-1)
fesm = V*Q
The space provides now a list of trial and test-functions:
sigma, u = fesm.TrialFunction()
tau, v = fesm.TestFunction()
The normal vector is provided as a special coefficient function (which may be used only at the boundary). The orientation depends on the definition of the geometry. In 2D, it is the tangential vector rotated to the right, and is the outer vector in our case. Since every CoefficientFunction must know its vector-dimension, we have to provide the spatial dimension:
normal = specialcf.normal(mesh.dim)
print (normal)
coef normal vector, real, dim=2
Define the forms on the product space. For the boundary term, we have to use the Trace operator, which provides the projection to normal direction.
# define the boundary conditions
electrode = mesh.BoundaryCF({"el1": 1, "el2": -1}, default=0)
ud = electrode
am = BilinearForm((1/lam*sigma*tau + div(sigma)*v + div(tau)*u)*dx).Assemble()
fm = LinearForm( ud*(tau.Trace()*normal)*ds("el.*")).Assemble()
gfm = GridFunction(fesm)
The proxy-functions used for the forms know to which component of the product space they belong to. To access components of the solution, we have to unpack its components. They don’t have own coefficient vectors, but refer to ranges of the big coefficient vector.
gfsigma, gfu = gfm.components
Just to try something:
Now set the essential boundary conditions for the flux part:
res = fm.vec.data - am.mat * gfm.vec
gfm.vec.data += am.mat.Inverse(freedofs=fesm.FreeDofs(), inverse="pardiso") * res
# solvers.BVP(bf=am, lf=fm, gf=gfm)
Draw (gfsigma, mesh, "flux-mixed")
Draw (gfu, mesh, "u-mixed");
Calculate the difference:
4.2.2. Post-processing for the scalar variable#
The scalar variable is approximated one order lower than the vector variable, what is its gradient. Knowing the gradient of a function more accurate, and knowing its mean value, one can recover the function itself. For this post-processing trick we refer to [Arnold+Brezzi 85]
find \(\widehat u \in P^{k+1, dc}\) and \(\widehat \lambda \in P^{0, dc}\) such that
fespost_u = L2(mesh, order=order_flux+1)
fespost_lam = L2(mesh, order=0)
fes_post = fespost_u*fespost_lam
(u,la),(v,mu) = fes_post.TnT()
a = BilinearForm( (lam*grad(u)*grad(v)+la*v+mu*u)*dx).Assemble()
f = LinearForm((gfsigma*grad(v)+gfu*mu)*dx).Assemble()
gfpost = GridFunction(fes_post)
gfpost.vec.data = a.mat.Inverse() * f.vec
# if we zoom in the deformation we can observe theat the elements are not conforming
Draw (gfpost.components[0], mesh, "upost", deformation=True)
BaseWebGuiScene