This page was generated from unit-2.5-mixed/mixed.ipynb.

2.5 Mixed formulation for second order equations

Motivation: * exact flux conservation * useful for a posteriori error estimates * model problem for \(4^{th}\) order problems, Stokes, …

We consider the diffusion equation

\[\begin{split}\DeclareMathOperator{\Div}{div} \begin{array}{rcll} -\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}\end{split}\]

Primal variational formulation

Find \(u \in H^1, u = u_D\) on \(\Gamma_D\) such that

\[\int_\Omega \lambda \nabla u \nabla v = \int_\Omega f v + \int_{\Gamma_N} g v \quad \forall v \in H^1, v = 0 \text{ on } \Gamma_D\]

First order system

Find scalar \(u\) and the flux \(\sigma\) such that

\[\lambda^{-1} \sigma = \nabla u, \quad \Div \sigma = -f\]

with boundary conditions

\[\sigma \cdot n = g \text{ on } \Gamma_N, \quad \text{ and } \quad u = u_D \text{ on } \Gamma_D\]

Mixed variational formulation

Find \((\sigma, u) \in H(\Div) \times L_2\) such that \(\sigma_n = g\) on \(\Gamma_N\) and

\[\int_\Omega \lambda^{-1} \sigma \tau + \Div \sigma v + \Div \tau u = -\int_\Omega f v + \int_{\Gamma_D} u_D \tau_n\]

for all test-functions \((\tau, v) \in H(\Div) \times L_2\) with \(\tau_n = 0\).

Here \(\sigma_n\) is the normal trace operator \(\sigma \cdot n |_{\Gamma_N}\), which is meaningful in \(H(\Div)\).

[1]:
from netgen.geom2d import unit_square
from ngsolve import *
from ngsolve.webgui import Draw

mesh = Mesh(unit_square.GenerateMesh(maxh=0.1))
[2]:
source = sin(3.14*x)
ud = CoefficientFunction(5)
g = CoefficientFunction([y*(1-y) if bc=="left" else 0 for bc in mesh.GetBoundaries()])
lam = 10

Setup and solve primal problem:

[3]:
fesp = H1(mesh, order=4, dirichlet="bottom")
up, vp = fesp.TnT()

ap = BilinearForm(fesp)
ap += lam*grad(up)*grad(vp)*dx
ap.Assemble()

fp = LinearForm(fesp)
fp += source*vp*dx + g*vp * ds
fp.Assemble()

gfup = GridFunction(fesp, "u-primal")
gfup.Set(ud, BND)

r = fp.vec.CreateVector()
r.data = fp.vec - ap.mat * gfup.vec
gfup.vec.data += ap.mat.Inverse(freedofs=fesp.FreeDofs()) * r
Draw (gfup)
Draw (lam * grad(gfup), mesh, "flux-primal")
[3]:
BaseWebGuiScene

Solving the mixed problem

Define spaces for mixed problem. Note that essential boundary conditions are set to the \(H(\Div)\)-component on the opposite boundary. Creating a space from a list of spaces generates a product space:

[4]:
order_flux=1
V = HDiv(mesh, order=order_flux, dirichlet="right|top|left")
Q = L2(mesh, order=order_flux-1)
fesm = FESpace([V,Q])

The space provides now a list of trial and test-functions:

[5]:
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:

[6]:
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.

[7]:
am = BilinearForm(fesm)
am += (1/lam*sigma*tau + div(sigma)*v + div(tau)*u)*dx
am.Assemble()
fm = LinearForm(fesm)
fm += -source*v*dx +  ud*(tau.Trace()*normal)*ds
fm.Assemble()

gfm = GridFunction(fesm, name="gfmixed")

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.

[8]:
gfsigma, gfu = gfm.components

Just to try something:

[9]:
gfsigma.Set(CoefficientFunction( (sin(10*x), sin(10*y)) ))
gfu.Set (x)
Draw (gfsigma, mesh, "flux-mixed")
Draw (gfu, mesh, "u-mixed")
[9]:
BaseWebGuiScene

Now set the essential boundary conditions for the flux part:

[10]:
gfsigma.Set(g*normal, BND)
[11]:
# fm.vec.data -= am.mat * gfm.vec
# gfm.vec.data += am.mat.Inverse(freedofs=fesm.FreeDofs(), inverse="umfpack") * fm.vec
solvers.BVP(bf=am, lf=fm, gf=gfm)
Draw (gfsigma, mesh, "flux-mixed")
Draw (gfu, mesh, "u-mixed")
[11]:
BaseWebGuiScene

Calculate the difference:

[12]:
print ("err-u:   ", sqrt(Integrate( (gfup-gfu)**2, mesh)))
errflux = lam * grad(gfup) - gfsigma
print ("err-flux:", sqrt(Integrate(errflux*errflux, mesh)))
err-u:    0.0010183129055822298
err-flux: 0.001453535299447794

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

\[\begin{split}\begin{array}{ccccl} \int \lambda \nabla \widehat u \nabla \widehat v & + & \int \widehat \lambda \widehat v & = & \int \sigma \nabla \widehat v & \forall \, \widehat v\\ \int \widehat u \widehat \mu & & & = & \int u \widehat \mu & \forall \, \widehat \mu \end{array}\end{split}\]
[13]:
fespost_u = L2(mesh, order=order_flux+1)
fespost_lam = L2(mesh, order=0)
fes_post = FESpace([fespost_u,fespost_lam])

u,la = fes_post.TrialFunction()
v,mu = fes_post.TestFunction()

a = BilinearForm(fes_post)
a += (lam*grad(u)*grad(v)+la*v+mu*u)*dx
a.Assemble()
f = LinearForm(fes_post)
f += (gfsigma*grad(v)+gfu*mu)*dx
f.Assemble()

gfpost = GridFunction(fes_post)
gfpost.vec.data = a.mat.Inverse() * f.vec

Draw (gfpost.components[0], mesh, "upost")
print ("err-upost:   ", sqrt(Integrate( (gfup-gfpost.components[0])**2, mesh)))
err-upost:    1.1258001649640946e-05
[ ]: