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/ngbem/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.94329927058765
GMRes iteration 2, residual = 10.547738043191677
GMRes iteration 3, residual = 2.6122774434371627
GMRes iteration 4, residual = 1.9193598549952564
GMRes iteration 5, residual = 0.41513121493724003
GMRes iteration 6, residual = 0.39044487452539667
GMRes iteration 7, residual = 0.17719618484093516
GMRes iteration 8, residual = 0.13791749808112977
GMRes iteration 9, residual = 0.08865764418669783
GMRes iteration 10, residual = 0.04770670569325261
GMRes iteration 11, residual = 0.046621092672894056
GMRes iteration 12, residual = 0.04497123843408976
GMRes iteration 13, residual = 0.02080387343973401
GMRes iteration 14, residual = 0.013029447085518918
GMRes iteration 15, residual = 0.01242885010088187
GMRes iteration 16, residual = 0.007375357480150988
GMRes iteration 17, residual = 0.007328667020748155
GMRes iteration 18, residual = 0.006116845094996015
GMRes iteration 19, residual = 0.005574750978768432
GMRes iteration 20, residual = 0.0036699794848263794
GMRes iteration 21, residual = 0.0035919764871824174
GMRes iteration 22, residual = 0.0029943262420837045
GMRes iteration 23, residual = 0.0029469307130900136
GMRes iteration 24, residual = 0.001979267988589037
GMRes iteration 25, residual = 0.0019408665966375875
GMRes iteration 26, residual = 0.0015575497467427476
GMRes iteration 27, residual = 0.001424623290267094
GMRes iteration 28, residual = 0.0012204426815134583
GMRes iteration 29, residual = 0.0010847044638458076
GMRes iteration 30, residual = 0.0009837866109222084
GMRes iteration 31, residual = 0.0007490403633009601
GMRes iteration 32, residual = 0.0007488563027355072
GMRes iteration 33, residual = 0.0006134913389615955
GMRes iteration 34, residual = 0.0006124101462970083
GMRes iteration 35, residual = 0.0004939639686464623
GMRes iteration 36, residual = 0.00047858687844772377
GMRes iteration 37, residual = 0.0003717419629269158
GMRes iteration 38, residual = 0.000370970559501477
GMRes iteration 39, residual = 0.00031421136916909
GMRes iteration 40, residual = 0.0003104027688915253
GMRes iteration 41, residual = 0.0002183736880826148
GMRes iteration 42, residual = 0.00021449536290087237
GMRes iteration 43, residual = 0.00019217139232859728
GMRes iteration 44, residual = 0.00017880288009988998
GMRes iteration 45, residual = 0.0001623610100438626
GMRes iteration 46, residual = 0.00013935860106523016
GMRes iteration 47, residual = 0.00013877790838410383
GMRes iteration 48, residual = 0.00011238309526630581
GMRes iteration 49, residual = 0.00011226441055210988
GMRes iteration 50, residual = 8.692525397286349e-05
GMRes iteration 51, residual = 8.664933665470679e-05
GMRes iteration 52, residual = 7.291906629814632e-05
GMRes iteration 53, residual = 7.25347733044331e-05
GMRes iteration 54, residual = 5.7047649102582595e-05
GMRes iteration 55, residual = 4.796457860195715e-05
GMRes iteration 56, residual = 4.6725296273216014e-05
GMRes iteration 57, residual = 3.8199326068255294e-05
GMRes iteration 58, residual = 3.807481560304819e-05
GMRes iteration 59, residual = 2.815216539063003e-05
GMRes iteration 60, residual = 2.7724918391572095e-05
GMRes iteration 61, residual = 2.4337904809523378e-05
GMRes iteration 62, residual = 2.4273020164859047e-05
GMRes iteration 63, residual = 2.036836419381598e-05
GMRes iteration 64, residual = 1.9873160190189243e-05
GMRes iteration 65, residual = 1.60380764691404e-05
GMRes iteration 66, residual = 1.4477033504297908e-05
GMRes iteration 67, residual = 1.4163716004744369e-05
GMRes iteration 68, residual = 1.0409828566087441e-05
GMRes iteration 69, residual = 1.0328826494164312e-05
GMRes iteration 70, residual = 8.190116204104962e-06
GMRes iteration 71, residual = 7.817838665845138e-06
GMRes iteration 72, residual = 7.018536558068308e-06
GMRes iteration 73, residual = 5.897828734523412e-06
GMRes iteration 74, residual = 5.545477398623388e-06
GMRes iteration 75, residual = 4.159616756242356e-06
GMRes iteration 76, residual = 4.128821672753138e-06
GMRes iteration 77, residual = 3.6120362662786823e-06
GMRes iteration 78, residual = 3.3516478987688176e-06
GMRes iteration 79, residual = 2.9846879774381086e-06
GMRes iteration 80, residual = 2.5627300490045576e-06
GMRes iteration 81, residual = 2.45224540996611e-06
GMRes iteration 82, residual = 2.1344375226514643e-06
GMRes iteration 83, residual = 2.0775045595766064e-06
GMRes iteration 84, residual = 1.611676389134043e-06
GMRes iteration 85, residual = 1.5868640825028031e-06
GMRes iteration 86, residual = 1.348775825689277e-06
GMRes iteration 87, residual = 1.229647490554399e-06
GMRes iteration 88, residual = 1.1847056226992243e-06
GMRes iteration 89, residual = 9.075413356006354e-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
[ ]: