In this tutorial, we apply a Zienkiewicz-Zhu type error estimator and run an adaptive loop with these steps: $$ \text{Solve}\rightarrow \text{Estimate}\rightarrow \text{Mark}\rightarrow \text{Refine}\rightarrow \text{Solve} \rightarrow \ldots $$
import netgen.gui
%gui tk
from ngsolve import *
from netgen.geom2d import SplineGeometry
The following geometry represents a heating chip embedded in another material that conducts away the heat.
# point numbers 0, 1, ... 11
# sub-domain numbers (1), (2), (3)
#
#
# 7-------------6
# | |
# | (2) |
# | |
# 3------4-------------5------2
# | |
# | 11 |
# | / \ |
# | 10 (3) 9 |
# | \ / (1) |
# | 8 |
# | |
# 0---------------------------1
#
def MakeGeometry():
geometry = SplineGeometry()
# point coordinates ...
pnts = [ (0,0), (1,0), (1,0.6), (0,0.6), \
(0.2,0.6), (0.8,0.6), (0.8,0.8), (0.2,0.8), \
(0.5,0.15), (0.65,0.3), (0.5,0.45), (0.35,0.3) ]
pnums = [geometry.AppendPoint(*p) for p in pnts]
# start-point, end-point, boundary-condition, domain on left side, domain on right side:
lines = [ (0,1,1,1,0), (1,2,2,1,0), (2,5,2,1,0), (5,4,2,1,2), (4,3,2,1,0), (3,0,2,1,0), \
(5,6,2,2,0), (6,7,2,2,0), (7,4,2,2,0), \
(8,9,2,3,1), (9,10,2,3,1), (10,11,2,3,1), (11,8,2,3,1) ]
for p1,p2,bc,left,right in lines:
geometry.Append(["line", pnums[p1], pnums[p2]], bc=bc, leftdomain=left, rightdomain=right)
geometry.SetMaterial(1,"base")
geometry.SetMaterial(2,"chip")
geometry.SetMaterial(3,"top")
return geometry
mesh = Mesh(MakeGeometry().GenerateMesh(maxh=0.2))
The problem is to find $u$ in $H_{0,D}^1$ satisfying $$ \int_\Omega \lambda \nabla u \cdot \nabla v = \int_\Omega f v $$ for all $v$ in $H_{0,D}^1$. We expect the solution to have singularities due to the nonconvex re-enrant angles and discontinuities in $\lambda$.
fes = H1(mesh, order=3, dirichlet=[1])
u = fes.TrialFunction()
v = fes.TestFunction()
# one heat conductivity coefficient per sub-domain
lam = CoefficientFunction([1, 1000, 10])
a = BilinearForm(fes)
a += SymbolicBFI(lam*grad(u)*grad(v))
# heat-source in inner subdomain
f = LinearForm(fes)
f += SymbolicLFI(CoefficientFunction([0, 0, 1])*v)
c = Preconditioner(a, type="multigrid", flags= {"inverse" : "sparsecholesky" })
gfu = GridFunction(fes)
Draw (gfu)
Note that the linear system is not yet assembled above.
Since we must solve multiple times, we define a function to solve the boundary value problem, where assembly, update, and solve occurs.
def SolveBVP():
fes.Update()
gfu.Update()
a.Assemble()
f.Assemble()
inv = CGSolver(a.mat, c.mat)
gfu.vec.data = inv * f.vec
Redraw (blocking=True)
SolveBVP()
We implement a gradient-recovery-type error estimator. For this, we need an H(div) space for flux recovery. We must compute the flux of the computed solution and interpolate it into this H(div) space.
space_flux = HDiv(mesh, order=2)
gf_flux = GridFunction(space_flux, "flux")
flux = lam * grad(gfu)
gf_flux.Set(flux)
err = 1/lam*(flux-gf_flux)*(flux-gf_flux)
Draw(err, mesh, 'error_representation')
Element-wise error estimator: On each element $T$, set
$$
\eta_T^2 = \int_T \frac{1}{\lambda}
|\lambda \nabla u_h - I_h(\lambda \nabla u_h) |^2
$$
where $u_h$ is the computed solution gfu
and $I_h$ is the interpolation performed by Set
in NGSolve.
eta2 = Integrate(err, mesh, VOL, element_wise=True)
print(eta2)
The above values, one per element, lead us to identify elements which might have large error.
We mark elements with large error estimator for refinement.
maxerr = max(eta2)
print ("maxerr = ", maxerr)
for el in mesh.Elements():
mesh.SetRefinementFlag(el, eta2[el.nr] > 0.25*maxerr)
Draw(gfu)
Refine marked elements:
mesh.Refine()
SolveBVP()
Redraw()
l = [] # l = list of estimated total error
def CalcError():
# compute the flux:
space_flux.Update()
gf_flux.Update()
flux = lam * grad(gfu)
gf_flux.Set(flux)
# compute estimator:
err = 1/lam*(flux-gf_flux)*(flux-gf_flux)
eta2 = Integrate(err, mesh, VOL, element_wise=True)
maxerr = max(eta2)
l.append ((fes.ndof, sqrt(sum(eta2))))
print("ndof =", fes.ndof, " maxerr =", maxerr)
# mark for refinement:
for el in mesh.Elements():
mesh.SetRefinementFlag(el, eta2[el.nr] > 0.25*maxerr)
CalcError()
mesh.Refine()
while fes.ndof < 50000:
SolveBVP()
CalcError()
mesh.Refine()
import matplotlib.pyplot as plt
plt.yscale('log')
plt.xscale('log')
plt.xlabel("ndof")
plt.ylabel("H1 error-estimate")
ndof,err = zip(*l)
plt.plot(ndof,err, "-*")
plt.ion()
plt.show()