We consider the diffusion equation $$ \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} $$
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 $$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 $$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)$.
from netgen.geom2d import unit_square
from ngsolve import *
import netgen.gui
%gui tk
mesh = Mesh(unit_square.GenerateMesh(maxh=0.1))
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:
fesp = H1(mesh, order=4, dirichlet="bottom")
up = fesp.TrialFunction()
vp = fesp.TestFunction()
ap = BilinearForm(fesp)
ap += SymbolicBFI(lam * grad(up) * grad(vp))
ap.Assemble()
fp = LinearForm(fesp)
fp += SymbolicLFI(source * vp)
fp += SymbolicLFI(g * vp, BND)
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")
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:
V = HDiv(mesh, order=4, dirichlet="right|top|left")
Q = L2(mesh, order=3)
fesm = FESpace([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)
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.
am = BilinearForm(fesm)
am += SymbolicBFI(1/lam * sigma * tau + div(sigma)*v + div(tau)* u)
am.Assemble()
fm = LinearForm(fesm)
fm += SymbolicLFI(-source*v)
fm += SymbolicLFI(ud * (tau.Trace()*normal), BND)
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.
gfsigma = gfm.components[0]
gfu = gfm.components[1]
Just to try something:
gfsigma.Set(CoefficientFunction( (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:
gfsigma.Set(g*normal, BND)
fm.vec.data -= am.mat * gfm.vec
gfm.vec.data += am.mat.Inverse(freedofs=fesm.FreeDofs(), inverse="umfpack") * fm.vec
Redraw()
Calculate the difference:
print ("err-u: ", sqrt(Integrate( (gfup-gfu)*(gfup-gfu), mesh)))
errflux = lam * grad(gfup) - gfsigma
print ("err-flux:", sqrt(Integrate(errflux*errflux, mesh)))