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:

\[x \in \Omega^c: \quad u(x) = - \int\limits_\Gamma \displaystyle{\frac{1}{4\,\pi}\, \frac{1}{\| x-y\|} } \, \gamma_1 u (y)\, \mathrm{d}\sigma_y + \int\limits_\Gamma \displaystyle{\frac{1}{4\,\pi}\, \frac{\langle n(y) , x-y\rangle }{\| x-y\|^3} } \, \gamma_0 u (y)\, \mathrm{d}\sigma_y\,,\]

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

\[\begin{split}\left( \begin{array}{c} \gamma_0 u \\ \gamma_1 u \end{array}\right) = \left( \begin{array}{cc} -V & \frac12 + K \\ \frac12 - K^\intercal & -D \end{array} \right) \left( \begin{array}{c} \gamma_1 u \\ \gamma_0 u \end{array}\right)\end{split}\]

.

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

\[\int_{\Omega_\text{FEM}} \nabla u \nabla v \, dx - \int_\Gamma \gamma_1 u v \, ds = 0 \qquad \forall \, v \in H^1(\Omega_\text{FEM})\]

We use Calderon’s represenataion formula for the Neumann trace:

\[\int_{\Omega_\text{FEM}} \nabla u \nabla v \, dx - \int_\Gamma \left( \left( \tfrac{1}{2} - K^\intercal\right) \,\gamma_1 u - D \, \gamma_0 u\right) v = 0 \qquad \forall \, v \in H^1(\Omega_\text{FEM})\]

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:

\[\begin{split}\left( \begin{array}{ccc } A_{\Omega\Omega} & A_{\Omega\Gamma} & 0 \\ A_{\Gamma\Omega} & A_{\Gamma\Gamma } + D & -\frac12 M^\intercal + K^\intercal \\ 0 & -\frac12 M + K & -V \end{array}\right) \left( \begin{array}{c} u \\ \gamma_0 u \\ \gamma_1 u \end{array}\right) = \left( \begin{array}{c} F_{\Omega} \\ F_{\Gamma}\\ 0 \end{array}\right) \,.\end{split}\]

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.94329927058705
GMRes iteration 2, residual = 10.547738043196185
GMRes iteration 3, residual = 3.4483285857124577
GMRes iteration 4, residual = 3.3925820325575398
GMRes iteration 5, residual = 0.5788213577922388
GMRes iteration 6, residual = 0.5480288625880863
GMRes iteration 7, residual = 0.3072544271329453
GMRes iteration 8, residual = 0.28418746451133714
GMRes iteration 9, residual = 0.18248120784000013
GMRes iteration 10, residual = 0.15076712396703834
GMRes iteration 11, residual = 0.10359016193302874
GMRes iteration 12, residual = 0.06769258395152024
GMRes iteration 13, residual = 0.04957676436248499
GMRes iteration 14, residual = 0.04077663433402211
GMRes iteration 15, residual = 0.029549247836700876
GMRes iteration 16, residual = 0.021772314044945167
GMRes iteration 17, residual = 0.009617158157111678
GMRes iteration 18, residual = 0.009543689912971574
GMRes iteration 19, residual = 0.007177978274388733
GMRes iteration 20, residual = 0.006301263636977168
GMRes iteration 21, residual = 0.005055306875327309
GMRes iteration 22, residual = 0.005054997706598148
GMRes iteration 23, residual = 0.0036880715436810804
GMRes iteration 24, residual = 0.003676864390851893
GMRes iteration 25, residual = 0.0027411432192132002
GMRes iteration 26, residual = 0.0025845367519143855
GMRes iteration 27, residual = 0.00212062062278903
GMRes iteration 28, residual = 0.0018677181260951911
GMRes iteration 29, residual = 0.0016414089939346971
GMRes iteration 30, residual = 0.0015190049514389006
GMRes iteration 31, residual = 0.0013581603782113604
GMRes iteration 32, residual = 0.0011735638752092593
GMRes iteration 33, residual = 0.0010714595699004789
GMRes iteration 34, residual = 0.0008878631699477931
GMRes iteration 35, residual = 0.0008809310919198435
GMRes iteration 36, residual = 0.0007096726488685535
GMRes iteration 37, residual = 0.0007090869341951916
GMRes iteration 38, residual = 0.0005151520364418145
GMRes iteration 39, residual = 0.0005068076929562785
GMRes iteration 40, residual = 0.0004220967159799718
GMRes iteration 41, residual = 0.000408528272105588
GMRes iteration 42, residual = 0.00034712193170832755
GMRes iteration 43, residual = 0.0003098960190515055
GMRes iteration 44, residual = 0.0002831356098673024
GMRes iteration 45, residual = 0.0002661608017260366
GMRes iteration 46, residual = 0.0002536860599653408
GMRes iteration 47, residual = 0.00022354514201723552
GMRes iteration 48, residual = 0.0002150908894815591
GMRes iteration 49, residual = 0.00018305145151611648
GMRes iteration 50, residual = 0.00018222563238710427
GMRes iteration 51, residual = 0.00015588339944943322
GMRes iteration 52, residual = 0.000155753234543838
GMRes iteration 53, residual = 0.00012174289331266637
GMRes iteration 54, residual = 0.00011889738165237885
GMRes iteration 55, residual = 0.00010362881748169164
GMRes iteration 56, residual = 9.575240949729111e-05
GMRes iteration 57, residual = 8.649894113614913e-05
GMRes iteration 58, residual = 7.242179504890145e-05
GMRes iteration 59, residual = 7.157448059659974e-05
GMRes iteration 60, residual = 5.904497319868405e-05
GMRes iteration 61, residual = 5.9039202133481436e-05
GMRes iteration 62, residual = 4.573590257515461e-05
GMRes iteration 63, residual = 4.379661969759127e-05
GMRes iteration 64, residual = 3.916792012985113e-05
GMRes iteration 65, residual = 3.384628170360831e-05
GMRes iteration 66, residual = 3.277690634785984e-05
GMRes iteration 67, residual = 2.6732862616437977e-05
GMRes iteration 68, residual = 2.6715944323492467e-05
GMRes iteration 69, residual = 2.3007588887598507e-05
GMRes iteration 70, residual = 2.2408332209782483e-05
GMRes iteration 71, residual = 1.961158759220571e-05
GMRes iteration 72, residual = 1.802702750604061e-05
GMRes iteration 73, residual = 1.704780166752383e-05
GMRes iteration 74, residual = 1.3001365366132553e-05
GMRes iteration 75, residual = 1.3001073738674487e-05
GMRes iteration 76, residual = 1.1003849200572383e-05
GMRes iteration 77, residual = 1.07665726306825e-05
GMRes iteration 78, residual = 9.105792854290656e-06
GMRes iteration 79, residual = 8.390870230079759e-06
GMRes iteration 80, residual = 7.625742650695371e-06
GMRes iteration 81, residual = 6.677096635811354e-06
GMRes iteration 82, residual = 6.589382999816681e-06
GMRes iteration 83, residual = 5.143769267279606e-06
GMRes iteration 84, residual = 5.035758651288849e-06
GMRes iteration 85, residual = 4.278268557934796e-06
GMRes iteration 86, residual = 4.015118049271159e-06
GMRes iteration 87, residual = 3.8039881357510457e-06
GMRes iteration 88, residual = 3.240824445942112e-06
GMRes iteration 89, residual = 3.186389400559058e-06
GMRes iteration 90, residual = 2.4850598392615407e-06
GMRes iteration 91, residual = 2.4128190350842354e-06
GMRes iteration 92, residual = 2.2270512078106948e-06
GMRes iteration 93, residual = 1.9475218560445574e-06
GMRes iteration 94, residual = 1.8639864828153503e-06
GMRes iteration 95, residual = 1.5663646255384106e-06
GMRes iteration 96, residual = 1.531092683931934e-06
GMRes iteration 97, residual = 1.2919079791723415e-06
GMRes iteration 98, residual = 1.2912893525413357e-06
GMRes iteration 99, residual = 9.998097061280072e-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:

[ ]: