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.94329927058755
GMRes iteration 2, residual = 10.539369746932744
GMRes iteration 3, residual = 2.611587156531044
GMRes iteration 4, residual = 1.9195757191367275
GMRes iteration 5, residual = 0.41669165805695724
GMRes iteration 6, residual = 0.39297352574755506
GMRes iteration 7, residual = 0.17795761889219466
GMRes iteration 8, residual = 0.13816526015257474
GMRes iteration 9, residual = 0.08671805759869106
GMRes iteration 10, residual = 0.047060025937660004
GMRes iteration 11, residual = 0.045344795301519754
GMRes iteration 12, residual = 0.04442191913718222
GMRes iteration 13, residual = 0.01670150254333788
GMRes iteration 14, residual = 0.012683670900526295
GMRes iteration 15, residual = 0.011115768703967474
GMRes iteration 16, residual = 0.007051225983157458
GMRes iteration 17, residual = 0.006714580260522305
GMRes iteration 18, residual = 0.005947766472135921
GMRes iteration 19, residual = 0.004793409918704949
GMRes iteration 20, residual = 0.0035486000831008047
GMRes iteration 21, residual = 0.0031764124771485307
GMRes iteration 22, residual = 0.002826216820531523
GMRes iteration 23, residual = 0.0024657041487294233
GMRes iteration 24, residual = 0.001863204293047011
GMRes iteration 25, residual = 0.0018471765310798427
GMRes iteration 26, residual = 0.001347896908018178
GMRes iteration 27, residual = 0.0013444294388166634
GMRes iteration 28, residual = 0.0010720995982663655
GMRes iteration 29, residual = 0.0010545592010335753
GMRes iteration 30, residual = 0.0007620541689166671
GMRes iteration 31, residual = 0.0006873897895759025
GMRes iteration 32, residual = 0.0006388026392287238
GMRes iteration 33, residual = 0.000552268269270044
GMRes iteration 34, residual = 0.0004835371772131719
GMRes iteration 35, residual = 0.00044792987828947914
GMRes iteration 36, residual = 0.00036474070614771095
GMRes iteration 37, residual = 0.0003467884081759419
GMRes iteration 38, residual = 0.00030219561951236665
GMRes iteration 39, residual = 0.00028858193078453256
GMRes iteration 40, residual = 0.00023060160538938005
GMRes iteration 41, residual = 0.0001999519113150148
GMRes iteration 42, residual = 0.00019138778612530627
GMRes iteration 43, residual = 0.00015891304799895453
GMRes iteration 44, residual = 0.00015494149238606477
GMRes iteration 45, residual = 0.00012383284911598037
GMRes iteration 46, residual = 0.00012381191857020163
GMRes iteration 47, residual = 9.62295311267093e-05
GMRes iteration 48, residual = 9.349714177633496e-05
GMRes iteration 49, residual = 7.68946957148184e-05
GMRes iteration 50, residual = 7.458496971614618e-05
GMRes iteration 51, residual = 6.382337099328349e-05
GMRes iteration 52, residual = 6.236933012711591e-05
GMRes iteration 53, residual = 4.9562817259132064e-05
GMRes iteration 54, residual = 4.124180682815443e-05
GMRes iteration 55, residual = 4.018343557303746e-05
GMRes iteration 56, residual = 3.2049203964697885e-05
GMRes iteration 57, residual = 3.18471220101913e-05
GMRes iteration 58, residual = 2.4094099808058865e-05
GMRes iteration 59, residual = 2.4009483362785104e-05
GMRes iteration 60, residual = 2.1199138645430218e-05
GMRes iteration 61, residual = 2.1157704554433656e-05
GMRes iteration 62, residual = 1.7273855200544128e-05
GMRes iteration 63, residual = 1.7055528279620253e-05
GMRes iteration 64, residual = 1.3764553642645894e-05
GMRes iteration 65, residual = 1.2844393430422604e-05
GMRes iteration 66, residual = 1.1069310052997508e-05
GMRes iteration 67, residual = 8.716928503148732e-06
GMRes iteration 68, residual = 8.571982306800147e-06
GMRes iteration 69, residual = 6.935384868113545e-06
GMRes iteration 70, residual = 6.9345940180380984e-06
GMRes iteration 71, residual = 5.211230013587087e-06
GMRes iteration 72, residual = 5.074424879102161e-06
GMRes iteration 73, residual = 4.2298351196508175e-06
GMRes iteration 74, residual = 3.6631646507713715e-06
GMRes iteration 75, residual = 3.589351241281644e-06
GMRes iteration 76, residual = 2.8148014188681046e-06
GMRes iteration 77, residual = 2.814783921071022e-06
GMRes iteration 78, residual = 2.265435391000516e-06
GMRes iteration 79, residual = 2.2419092720416823e-06
GMRes iteration 80, residual = 1.9689483431552116e-06
GMRes iteration 81, residual = 1.9162861173351253e-06
GMRes iteration 82, residual = 1.7652217267979612e-06
GMRes iteration 83, residual = 1.3683867901488392e-06
GMRes iteration 84, residual = 1.3585483908849019e-06
GMRes iteration 85, residual = 1.0769021910072119e-06
GMRes iteration 86, residual = 1.0765310007919902e-06
GMRes iteration 87, residual = 8.756346076518212e-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
[ ]: