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.94329927058755
GMRes iteration 2, residual = 10.547738043200772
GMRes iteration 3, residual = 3.4483285857126784
GMRes iteration 4, residual = 3.392582032557781
GMRes iteration 5, residual = 0.5788213577922074
GMRes iteration 6, residual = 0.5480288625881221
GMRes iteration 7, residual = 0.3072544271327939
GMRes iteration 8, residual = 0.2841874645111465
GMRes iteration 9, residual = 0.1824812078400896
GMRes iteration 10, residual = 0.1507671235950272
GMRes iteration 11, residual = 0.10359016191343255
GMRes iteration 12, residual = 0.06769258397843109
GMRes iteration 13, residual = 0.04957676433661125
GMRes iteration 14, residual = 0.04077663438781509
GMRes iteration 15, residual = 0.029549248738200473
GMRes iteration 16, residual = 0.02177231567856068
GMRes iteration 17, residual = 0.009617154326378692
GMRes iteration 18, residual = 0.009543687101545366
GMRes iteration 19, residual = 0.007177973272878086
GMRes iteration 20, residual = 0.006301261570665648
GMRes iteration 21, residual = 0.005055304866041536
GMRes iteration 22, residual = 0.005054995717526628
GMRes iteration 23, residual = 0.0036880696907067327
GMRes iteration 24, residual = 0.0036768629048953465
GMRes iteration 25, residual = 0.0027411430596397572
GMRes iteration 26, residual = 0.0025845366098847626
GMRes iteration 27, residual = 0.002120620352569214
GMRes iteration 28, residual = 0.001867718076944775
GMRes iteration 29, residual = 0.0016414092270849768
GMRes iteration 30, residual = 0.0015190049603145993
GMRes iteration 31, residual = 0.0013581604963692126
GMRes iteration 32, residual = 0.0011735636453785354
GMRes iteration 33, residual = 0.0010714592425202756
GMRes iteration 34, residual = 0.0008878630459072692
GMRes iteration 35, residual = 0.000880931007966234
GMRes iteration 36, residual = 0.0007096725923288402
GMRes iteration 37, residual = 0.0007090869058464017
GMRes iteration 38, residual = 0.0005151521223399838
GMRes iteration 39, residual = 0.0005068076762784271
GMRes iteration 40, residual = 0.0004220967453020624
GMRes iteration 41, residual = 0.0004085282630930321
GMRes iteration 42, residual = 0.00034712202287610883
GMRes iteration 43, residual = 0.00030989598008640274
GMRes iteration 44, residual = 0.0002831356083425698
GMRes iteration 45, residual = 0.00026616071125240727
GMRes iteration 46, residual = 0.00025368602275096933
GMRes iteration 47, residual = 0.0002235450674676796
GMRes iteration 48, residual = 0.00021509094019227178
GMRes iteration 49, residual = 0.0001830513814265563
GMRes iteration 50, residual = 0.0001822255885720813
GMRes iteration 51, residual = 0.00015588335559792166
GMRes iteration 52, residual = 0.00015575319659043435
GMRes iteration 53, residual = 0.00012174296275482825
GMRes iteration 54, residual = 0.00011889742135243199
GMRes iteration 55, residual = 0.00010362889362512071
GMRes iteration 56, residual = 9.575243859721118e-05
GMRes iteration 57, residual = 8.64989862207691e-05
GMRes iteration 58, residual = 7.242180401820503e-05
GMRes iteration 59, residual = 7.157449491427172e-05
GMRes iteration 60, residual = 5.904494915375189e-05
GMRes iteration 61, residual = 5.90391791335325e-05
GMRes iteration 62, residual = 4.573585808971645e-05
GMRes iteration 63, residual = 4.379659941161777e-05
GMRes iteration 64, residual = 3.916789452076836e-05
GMRes iteration 65, residual = 3.384627488618706e-05
GMRes iteration 66, residual = 3.2776889914174414e-05
GMRes iteration 67, residual = 2.6732852686483722e-05
GMRes iteration 68, residual = 2.6715936662759847e-05
GMRes iteration 69, residual = 2.3007581833782443e-05
GMRes iteration 70, residual = 2.2408336914611477e-05
GMRes iteration 71, residual = 1.9611589648531623e-05
GMRes iteration 72, residual = 1.8027035437390105e-05
GMRes iteration 73, residual = 1.7047808275474947e-05
GMRes iteration 74, residual = 1.3001367120113446e-05
GMRes iteration 75, residual = 1.300107554913688e-05
GMRes iteration 76, residual = 1.1003848479401298e-05
GMRes iteration 77, residual = 1.0766573668627699e-05
GMRes iteration 78, residual = 9.105793241820545e-06
GMRes iteration 79, residual = 8.390871053348453e-06
GMRes iteration 80, residual = 7.6257430891138095e-06
GMRes iteration 81, residual = 6.677094359708035e-06
GMRes iteration 82, residual = 6.589380239273722e-06
GMRes iteration 83, residual = 5.143764617357629e-06
GMRes iteration 84, residual = 5.035756017874148e-06
GMRes iteration 85, residual = 4.27826337348564e-06
GMRes iteration 86, residual = 4.015115594158826e-06
GMRes iteration 87, residual = 3.803982487523403e-06
GMRes iteration 88, residual = 3.2408215790213194e-06
GMRes iteration 89, residual = 3.1863842852630077e-06
GMRes iteration 90, residual = 2.4850564517874906e-06
GMRes iteration 91, residual = 2.4128169897917557e-06
GMRes iteration 92, residual = 2.22705020098585e-06
GMRes iteration 93, residual = 1.94752053268344e-06
GMRes iteration 94, residual = 1.8639861554316754e-06
GMRes iteration 95, residual = 1.5663637497545076e-06
GMRes iteration 96, residual = 1.5310920905741415e-06
GMRes iteration 97, residual = 1.2919076340443824e-06
GMRes iteration 98, residual = 1.291289018753588e-06
GMRes iteration 99, residual = 9.998098831291678e-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:

[ ]: