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.547738043195999
GMRes iteration 3, residual = 2.612277474514308
GMRes iteration 4, residual = 1.9193596826204513
GMRes iteration 5, residual = 0.4151313306740289
GMRes iteration 6, residual = 0.39044492013918564
GMRes iteration 7, residual = 0.1771959584688743
GMRes iteration 8, residual = 0.13791736363646095
GMRes iteration 9, residual = 0.08865704985751684
GMRes iteration 10, residual = 0.047706546998883374
GMRes iteration 11, residual = 0.04662084838180726
GMRes iteration 12, residual = 0.04497110881914466
GMRes iteration 13, residual = 0.020803556528803415
GMRes iteration 14, residual = 0.013029444779844227
GMRes iteration 15, residual = 0.012428805616745218
GMRes iteration 16, residual = 0.0073753808217687745
GMRes iteration 17, residual = 0.007328681353587161
GMRes iteration 18, residual = 0.006116865233195386
GMRes iteration 19, residual = 0.005574735701212088
GMRes iteration 20, residual = 0.003669996591597239
GMRes iteration 21, residual = 0.0035919779200622136
GMRes iteration 22, residual = 0.0029943373277455596
GMRes iteration 23, residual = 0.0029469284845469973
GMRes iteration 24, residual = 0.0019792742948276155
GMRes iteration 25, residual = 0.0019408798214113644
GMRes iteration 26, residual = 0.0015575474869656079
GMRes iteration 27, residual = 0.0014246264757687376
GMRes iteration 28, residual = 0.001220435555438055
GMRes iteration 29, residual = 0.0010847085199403726
GMRes iteration 30, residual = 0.0009837759801017115
GMRes iteration 31, residual = 0.0007490407681332909
GMRes iteration 32, residual = 0.0007488572698788029
GMRes iteration 33, residual = 0.0006134945187180244
GMRes iteration 34, residual = 0.0006124123794552847
GMRes iteration 35, residual = 0.0004939701568627521
GMRes iteration 36, residual = 0.00047859098529587216
GMRes iteration 37, residual = 0.00037174511788517416
GMRes iteration 38, residual = 0.00037097383005258003
GMRes iteration 39, residual = 0.00031421434544862863
GMRes iteration 40, residual = 0.0003104065593116147
GMRes iteration 41, residual = 0.0002183771130508382
GMRes iteration 42, residual = 0.0002144971882939912
GMRes iteration 43, residual = 0.00019217529055723625
GMRes iteration 44, residual = 0.00017880441950033073
GMRes iteration 45, residual = 0.00016236399028432856
GMRes iteration 46, residual = 0.00013935962163312288
GMRes iteration 47, residual = 0.00013877943555283708
GMRes iteration 48, residual = 0.00011238424531434909
GMRes iteration 49, residual = 0.00011226538072709904
GMRes iteration 50, residual = 8.692583780380783e-05
GMRes iteration 51, residual = 8.66496764198899e-05
GMRes iteration 52, residual = 7.291955504901833e-05
GMRes iteration 53, residual = 7.253496453604106e-05
GMRes iteration 54, residual = 5.704810330246421e-05
GMRes iteration 55, residual = 4.796419002043011e-05
GMRes iteration 56, residual = 4.672499559647277e-05
GMRes iteration 57, residual = 3.819894343197312e-05
GMRes iteration 58, residual = 3.8074462418167446e-05
GMRes iteration 59, residual = 2.8152083531700356e-05
GMRes iteration 60, residual = 2.772480737481193e-05
GMRes iteration 61, residual = 2.4337921910361945e-05
GMRes iteration 62, residual = 2.4273028317065255e-05
GMRes iteration 63, residual = 2.0368421505686022e-05
GMRes iteration 64, residual = 1.9873153887986e-05
GMRes iteration 65, residual = 1.6038181652319842e-05
GMRes iteration 66, residual = 1.4477108124854297e-05
GMRes iteration 67, residual = 1.4163794996176848e-05
GMRes iteration 68, residual = 1.0409852351551227e-05
GMRes iteration 69, residual = 1.0328843174170287e-05
GMRes iteration 70, residual = 8.190109957482001e-06
GMRes iteration 71, residual = 7.81784073213391e-06
GMRes iteration 72, residual = 7.018519646148358e-06
GMRes iteration 73, residual = 5.8978042004634055e-06
GMRes iteration 74, residual = 5.545462739456978e-06
GMRes iteration 75, residual = 4.159597043106153e-06
GMRes iteration 76, residual = 4.128792693749889e-06
GMRes iteration 77, residual = 3.612000332161141e-06
GMRes iteration 78, residual = 3.351495780061599e-06
GMRes iteration 79, residual = 2.9846865727569175e-06
GMRes iteration 80, residual = 2.5626832076157546e-06
GMRes iteration 81, residual = 2.452254796139947e-06
GMRes iteration 82, residual = 2.1344356011103567e-06
GMRes iteration 83, residual = 2.0775153781217196e-06
GMRes iteration 84, residual = 1.6116773846755133e-06
GMRes iteration 85, residual = 1.5868612393630693e-06
GMRes iteration 86, residual = 1.3487745930763676e-06
GMRes iteration 87, residual = 1.2296483081939284e-06
GMRes iteration 88, residual = 1.184702020396127e-06
GMRes iteration 89, residual = 9.07539284347495e-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
[ ]: