This page was generated from wta/fembem.ipynb.
FEM-BEM Coupling¶
The ngbem boundary element addon project initiated by Lucy Weggeler (see https://weggler.github.io/docu-ngsbem/intro.html) is now partly integrated into core NGSolve. Find a short and sweet introduction to the boundary element method there.
In this demo we simulate a plate capacitor on an unbounded domain.
[1]:
from ngsolve import *
from netgen.occ import *
from ngsolve.solvers import GMRes
from ngsolve.webgui import Draw
from ngsolve.bem import *
[2]:
largebox = Box ((-2,-2,-2), (2,2,2) )
eltop = Box ( (-1,-1,0.5), (1,1,1) )
elbot = Box ( (-1,-1,-1), (1,1,-0.5))
largebox.faces.name = "outer" # coupling boundary
eltop.faces.name = "topface" # Dirichlet boundary
elbot.faces.name = "botface" # Dirichlet boundary
eltop.edges.hpref = 1
elbot.edges.hpref = 1
shell = largebox-eltop-elbot # FEM domain
shell.solids.name = "air"
mesh = shell.GenerateMesh(maxh=0.8)
mesh.RefineHP(2)
ea = { "euler_angles" : (-67, 0, 110) }
Draw (mesh, clipping={"x":1, "y":0, "z":0, "dist" : 1.1}, **ea);
On the exterior domain \(\Omega^c\), the solution can be expressed by the representation formula:
where \(\gamma_0 u = u\) and \(\gamma_1 u = \frac{\partial u}{\partial n}\) are Dirichlet and Neumann traces. These traces are related by the Calderon projector
.
The \(V\), \(K\) are the single layer and double layer potential operators, and \(D\) is the hypersingular operator.
On the FEM domain we have the variational formulation
We use Calderon’s represenataion formula for the Neumann trace:
To get a closed system, we use also the first equation of the Calderon equations. To see the structure of the discretized system, the dofs are split into degrees of freedom inside \(\Omega\), and those on the boundary \(\Gamma\). The FEM matrix \(A\) is split accordingly. We see, the coupled system is symmetric, but indefinite:
Generate the finite element space for \(H^1(\Omega)\) and set the given Dirichlet boundary conditions on the surfaces of the plates:
[3]:
order = 4
fesH1 = H1(mesh, order=order, dirichlet="topface|botface")
print ("H1-ndof = ", fesH1.ndof)
H1-ndof = 90703
The finite element space \(\verb-fesH1-\) provides \(H^{\frac12}(\Gamma)\) conforming element to discretize the Dirichlet trace on the coupling boundary \(\Gamma\). However we still need \(H^{-\frac12}(\Gamma)\) conforming elements to discretize the Neumann trace of \(u\) on the coupling boundary. Here it is:
[4]:
fesL2 = SurfaceL2(mesh, order=order-1, dual_mapping=True, definedon=mesh.Boundaries("outer"))
print ("L2-ndof = ", fesL2.ndof)
L2-ndof = 4035
[5]:
fes = fesH1 * fesL2
u,dudn = fes.TrialFunction()
v,dvdn = fes.TestFunction()
a = BilinearForm(grad(u)*grad(v)*dx, check_unused=False).Assemble()
gfudir = GridFunction(fes)
gfudir.components[0].Set ( mesh.BoundaryCF( { "topface" : 1, "botface" : -1 }), BND)
f = LinearForm(fes).Assemble()
res = (f.vec - a.mat * gfudir.vec).Evaluate()
Generate the the single layer potential \(V\), double layer potential \(K\) and hypersingular operator \(D\):
[6]:
n = specialcf.normal(3)
with TaskManager():
V = LaplaceSL(dudn*ds("outer"))*dvdn*ds("outer")
K = LaplaceDL(u*ds("outer"))*dvdn*ds("outer")
D = LaplaceSL(Cross(grad(u).Trace(),n)*ds("outer"))*Cross(grad(v).Trace(),n)*ds("outer")
M = BilinearForm(u*dvdn*ds("outer"), check_unused=False).Assemble()
Setup the coupled system matrix and the right hand side:
[7]:
sym = a.mat+D.mat - (0.5*M.mat+K.mat).T - (0.5*M.mat+K.mat) - V.mat
rhs = res
bfpre = BilinearForm(grad(u)*grad(v)*dx+1e-10*u*v*dx + dudn*dvdn*ds("outer") ).Assemble()
pre = bfpre.mat.Inverse(freedofs=fes.FreeDofs(), inverse="sparsecholesky")
Compute the solution of the coupled system:
[8]:
with TaskManager():
sol_sym = GMRes(A=sym, b=rhs, pre=pre, tol=1e-6, maxsteps=200, printrates=True)
GMRes iteration 1, residual = 47.943299270588106
GMRes iteration 2, residual = 10.547738043196349
GMRes iteration 3, residual = 2.612277443437126
GMRes iteration 4, residual = 1.9193598549949833
GMRes iteration 5, residual = 0.41513121493989497
GMRes iteration 6, residual = 0.39044487452891075
GMRes iteration 7, residual = 0.17719618484678076
GMRes iteration 8, residual = 0.13791749808247258
GMRes iteration 9, residual = 0.08865764419699401
GMRes iteration 10, residual = 0.04770670569335261
GMRes iteration 11, residual = 0.046621092554917545
GMRes iteration 12, residual = 0.04497123877047272
GMRes iteration 13, residual = 0.020803872508183916
GMRes iteration 14, residual = 0.01302944715267673
GMRes iteration 15, residual = 0.012428849701733793
GMRes iteration 16, residual = 0.007375357812441008
GMRes iteration 17, residual = 0.007328667494335652
GMRes iteration 18, residual = 0.006116845103964604
GMRes iteration 19, residual = 0.0055747512054370695
GMRes iteration 20, residual = 0.003669979651581257
GMRes iteration 21, residual = 0.0035919767415667644
GMRes iteration 22, residual = 0.0029943263052337786
GMRes iteration 23, residual = 0.0029469307910805404
GMRes iteration 24, residual = 0.0019792680085019853
GMRes iteration 25, residual = 0.001940866609625777
GMRes iteration 26, residual = 0.0015575496344463668
GMRes iteration 27, residual = 0.0014246231875360866
GMRes iteration 28, residual = 0.0012204425806869351
GMRes iteration 29, residual = 0.0010847044876622376
GMRes iteration 30, residual = 0.0009837866440103882
GMRes iteration 31, residual = 0.0007490403244034343
GMRes iteration 32, residual = 0.0007488562612908829
GMRes iteration 33, residual = 0.0006134912662545438
GMRes iteration 34, residual = 0.0006124100913657718
GMRes iteration 35, residual = 0.0004939639241532798
GMRes iteration 36, residual = 0.00047858689573854874
GMRes iteration 37, residual = 0.0003717419668628652
GMRes iteration 38, residual = 0.00037097056706578027
GMRes iteration 39, residual = 0.0003142113602324799
GMRes iteration 40, residual = 0.00031040276053255174
GMRes iteration 41, residual = 0.00021837369460081445
GMRes iteration 42, residual = 0.00021449537693121585
GMRes iteration 43, residual = 0.00019217140473600918
GMRes iteration 44, residual = 0.00017880296704259237
GMRes iteration 45, residual = 0.00016236102805599248
GMRes iteration 46, residual = 0.00013935861766658353
GMRes iteration 47, residual = 0.00013877792359174608
GMRes iteration 48, residual = 0.00011238312343492715
GMRes iteration 49, residual = 0.00011226444092507389
GMRes iteration 50, residual = 8.692527653387605e-05
GMRes iteration 51, residual = 8.664936085325821e-05
GMRes iteration 52, residual = 7.291910130359435e-05
GMRes iteration 53, residual = 7.25348158244282e-05
GMRes iteration 54, residual = 5.7047644138445705e-05
GMRes iteration 55, residual = 4.7964540200689665e-05
GMRes iteration 56, residual = 4.6725283752009727e-05
GMRes iteration 57, residual = 3.819921753944568e-05
GMRes iteration 58, residual = 3.8074726387271704e-05
GMRes iteration 59, residual = 2.815215347544157e-05
GMRes iteration 60, residual = 2.772489957500856e-05
GMRes iteration 61, residual = 2.4337896840966394e-05
GMRes iteration 62, residual = 2.4273011331130447e-05
GMRes iteration 63, residual = 2.036833959943495e-05
GMRes iteration 64, residual = 1.987311228801424e-05
GMRes iteration 65, residual = 1.603807105381812e-05
GMRes iteration 66, residual = 1.4477033159720984e-05
GMRes iteration 67, residual = 1.4163714666162307e-05
GMRes iteration 68, residual = 1.0409822052527678e-05
GMRes iteration 69, residual = 1.0328819767802272e-05
GMRes iteration 70, residual = 8.190115556817571e-06
GMRes iteration 71, residual = 7.817839199647946e-06
GMRes iteration 72, residual = 7.018534793077373e-06
GMRes iteration 73, residual = 5.897813784440728e-06
GMRes iteration 74, residual = 5.545469688078597e-06
GMRes iteration 75, residual = 4.159611606363843e-06
GMRes iteration 76, residual = 4.128814162142612e-06
GMRes iteration 77, residual = 3.61201954019791e-06
GMRes iteration 78, residual = 3.3515770490497482e-06
GMRes iteration 79, residual = 2.984680885947756e-06
GMRes iteration 80, residual = 2.562703125826963e-06
GMRes iteration 81, residual = 2.4522388781859235e-06
GMRes iteration 82, residual = 2.13442896458787e-06
GMRes iteration 83, residual = 2.077499764538598e-06
GMRes iteration 84, residual = 1.6116741008183059e-06
GMRes iteration 85, residual = 1.5868600309042541e-06
GMRes iteration 86, residual = 1.3487756791279234e-06
GMRes iteration 87, residual = 1.2296459301583446e-06
GMRes iteration 88, residual = 1.1847049978615616e-06
GMRes iteration 89, residual = 9.075407246279038e-07
[9]:
gfu = GridFunction(fes)
gfu.vec[:] = sol_sym + gfudir.vec
Draw(gfu.components[0], clipping={"x" : 1, "y":0, "z":0, "dist":0.0, "function" : True }, **ea, order=2);
The Neumann data:
[10]:
Draw (gfu.components[1], **ea);
References:
M. Costabel: Principles of boundary element methods
[ ]: