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 4th order problems, Stokes, …

We consider the diffusion equation

divλu=f in Ωu=uD on ΓDλun=g on ΓN

Primal variational formulation

Find uH1,u=uD on ΓD such that

Ωλuv=Ωfv+ΓNgvvH1,v=0 on ΓD

First order system

Find scalar u and the flux σ such that

λ1σ=u,divσ=f

with boundary conditions

σn=g on ΓN, and u=uD on ΓD

Mixed variational formulation

Find σH(div) and uL2 such that σn=g on ΓN and

Ωλ1στ+Ωdivτu=0Ωdivσv=Ωfv+ΓDuDτn

for all test-functions τH(div) and vL2 with τn=0.

Here σn is the normal trace operator σn|ΓN, which is meaningful in H(div).

A Compact notation is the single-liner

Find (σ,u)H(div)×L2 such that σn=g on ΓN and

Ωλ1στ+divσv+divτu=Ωfv+ΓDuDτn

for all test-functions (τ,v)H(div)×L2 with τn=0.

[1]:
from ngsolve import *
from ngsolve.webgui import Draw

mesh = Mesh(unit_square.GenerateMesh(maxh=0.1))
[2]:
source = sin(3.14*x)
ud = CF(5)
g = mesh.BoundaryCF( {"left" : y*(1-y)}, default=0)
lam = 10

Setup and solve primal problem:

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

ap = BilinearForm(lam*grad(up)*grad(vp)*dx).Assemble()
fp = LinearForm(source*vp*dx + g*vp * ds).Assemble()

gfup = GridFunction(fesp)
gfup.Set(ud, BND)

r = 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");

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 = 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((1/lam*sigma*tau + div(sigma)*v + div(tau)*u)*dx).Assemble()
fm = LinearForm(-source*v*dx + ud*(tau.Trace()*normal)*ds).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.

[8]:
gfsigma, gfu = gfm.components

Just to try something:

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

Now set the essential boundary conditions for the flux part:

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

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.0010272008852204048
err-flux: 0.0014920319417886373

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 u^Pk+1,dc and λ^P0,dc such that

λu^v^+λ^v^=σv^v^u^μ^=uμ^μ^
[13]:
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

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

[ ]: