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.94329927058823
GMRes iteration 2, residual = 10.547738043193325
GMRes iteration 3, residual = 2.6122774434371316
GMRes iteration 4, residual = 1.9193598549952613
GMRes iteration 5, residual = 0.4151312149374159
GMRes iteration 6, residual = 0.39044487452562565
GMRes iteration 7, residual = 0.17719618484132052
GMRes iteration 8, residual = 0.13791749808140097
GMRes iteration 9, residual = 0.0886576441870934
GMRes iteration 10, residual = 0.04770670569347106
GMRes iteration 11, residual = 0.04662109255457841
GMRes iteration 12, residual = 0.04497123877080145
GMRes iteration 13, residual = 0.020803872510229676
GMRes iteration 14, residual = 0.013029447152889687
GMRes iteration 15, residual = 0.012428849702681011
GMRes iteration 16, residual = 0.007375357812987344
GMRes iteration 17, residual = 0.007328667494713534
GMRes iteration 18, residual = 0.006116845104316371
GMRes iteration 19, residual = 0.005574751205831709
GMRes iteration 20, residual = 0.003669979651684257
GMRes iteration 21, residual = 0.0035919767417809086
GMRes iteration 22, residual = 0.002994326305410735
GMRes iteration 23, residual = 0.0029469307913325463
GMRes iteration 24, residual = 0.0019792680085865327
GMRes iteration 25, residual = 0.0019408666096342923
GMRes iteration 26, residual = 0.001557549634615168
GMRes iteration 27, residual = 0.0014246231875803454
GMRes iteration 28, residual = 0.0012204425808163266
GMRes iteration 29, residual = 0.0010847044876405555
GMRes iteration 30, residual = 0.0009837866441093782
GMRes iteration 31, residual = 0.000749040324429303
GMRes iteration 32, residual = 0.0007488562613133386
GMRes iteration 33, residual = 0.0006134912662315953
GMRes iteration 34, residual = 0.0006124100913278831
GMRes iteration 35, residual = 0.0004939639241300725
GMRes iteration 36, residual = 0.00047858689569459556
GMRes iteration 37, residual = 0.00037174196683985295
GMRes iteration 38, residual = 0.0003709705670467327
GMRes iteration 39, residual = 0.00031421136030416455
GMRes iteration 40, residual = 0.00031040276060670933
GMRes iteration 41, residual = 0.00021837369457490392
GMRes iteration 42, residual = 0.00021449537691047268
GMRes iteration 43, residual = 0.00019217140462721258
GMRes iteration 44, residual = 0.0001788029670156233
GMRes iteration 45, residual = 0.00016236102798872397
GMRes iteration 46, residual = 0.00013935861765601842
GMRes iteration 47, residual = 0.00013877792357912787
GMRes iteration 48, residual = 0.00011238312342499078
GMRes iteration 49, residual = 0.00011226444091701894
GMRes iteration 50, residual = 8.692527648485735e-05
GMRes iteration 51, residual = 8.664936081321003e-05
GMRes iteration 52, residual = 7.291910128468859e-05
GMRes iteration 53, residual = 7.253481584409563e-05
GMRes iteration 54, residual = 5.704764438893811e-05
GMRes iteration 55, residual = 4.796454046018698e-05
GMRes iteration 56, residual = 4.6725284164729884e-05
GMRes iteration 57, residual = 3.819921772393178e-05
GMRes iteration 58, residual = 3.807472657839828e-05
GMRes iteration 59, residual = 2.815215346523813e-05
GMRes iteration 60, residual = 2.772489960918793e-05
GMRes iteration 61, residual = 2.4337896756309134e-05
GMRes iteration 62, residual = 2.4273011266735696e-05
GMRes iteration 63, residual = 2.0368339473942554e-05
GMRes iteration 64, residual = 1.9873112253767385e-05
GMRes iteration 65, residual = 1.6038071022703806e-05
GMRes iteration 66, residual = 1.4477033143172106e-05
GMRes iteration 67, residual = 1.4163714566702878e-05
GMRes iteration 68, residual = 1.040982205395557e-05
GMRes iteration 69, residual = 1.0328819788063796e-05
GMRes iteration 70, residual = 8.19011556730384e-06
GMRes iteration 71, residual = 7.817839225577603e-06
GMRes iteration 72, residual = 7.018534787770111e-06
GMRes iteration 73, residual = 5.8978138052035984e-06
GMRes iteration 74, residual = 5.545469723485792e-06
GMRes iteration 75, residual = 4.15961170821016e-06
GMRes iteration 76, residual = 4.1288142552402566e-06
GMRes iteration 77, residual = 3.6120195945247696e-06
GMRes iteration 78, residual = 3.351577097604422e-06
GMRes iteration 79, residual = 2.9846808906892785e-06
GMRes iteration 80, residual = 2.562703131791921e-06
GMRes iteration 81, residual = 2.452238885374249e-06
GMRes iteration 82, residual = 2.134428960751776e-06
GMRes iteration 83, residual = 2.0774997639547917e-06
GMRes iteration 84, residual = 1.611674108765397e-06
GMRes iteration 85, residual = 1.5868600359149258e-06
GMRes iteration 86, residual = 1.348775697395756e-06
GMRes iteration 87, residual = 1.2296459389365097e-06
GMRes iteration 88, residual = 1.184705010387461e-06
GMRes iteration 89, residual = 9.075407298071629e-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
[ ]: