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
Primal variational formulation¶
Find \(u \in H^1, u = u_D\) on \(\Gamma_D\) such that
First order system¶
Find scalar \(u\) and the flux \(\sigma\) such that
with boundary conditions
Mixed variational formulation¶
Find \(\sigma \in H(\Div)\) and \(u \in L_2\) such that \(\sigma_n = g\) on \(\Gamma_N\) and
\begin{eqnarray*} \int_\Omega \lambda^{-1} \sigma \tau + \int_\Omega \Div \tau u & = & 0 \\ \int_\Omega \Div \sigma v & = & -\int_\Omega f v + \int_{\Gamma_D} u_D \tau_n \end{eqnarray*}
for all test-functions \(\tau \in H(\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(\Div)\).
A Compact notation is the single-liner
Find \((\sigma, u) \in H(\Div) \times L_2\) such that \(\sigma_n = g\) on \(\Gamma_N\) and
for all test-functions \((\tau, v) \in H(\Div) \times L_2\) with \(\tau_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.0010366712251541224
err-flux: 0.0014966050558590056
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
[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.1735501667459742e-05
[ ]:
[ ]: