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:

\[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.94329927058664
GMRes iteration 2, residual = 10.547738043196178
GMRes iteration 3, residual = 2.6122774434372293
GMRes iteration 4, residual = 1.919359854995171
GMRes iteration 5, residual = 0.41513121493734456
GMRes iteration 6, residual = 0.39044487452540866
GMRes iteration 7, residual = 0.17719618484152402
GMRes iteration 8, residual = 0.1379174980815782
GMRes iteration 9, residual = 0.08865764418873984
GMRes iteration 10, residual = 0.047706705693637816
GMRes iteration 11, residual = 0.0466210925552229
GMRes iteration 12, residual = 0.04497123877092446
GMRes iteration 13, residual = 0.020803872511064633
GMRes iteration 14, residual = 0.013029447152888268
GMRes iteration 15, residual = 0.012428849702764483
GMRes iteration 16, residual = 0.007375357812934093
GMRes iteration 17, residual = 0.007328667497114501
GMRes iteration 18, residual = 0.006116845107930565
GMRes iteration 19, residual = 0.005574751209948709
GMRes iteration 20, residual = 0.003669979651820781
GMRes iteration 21, residual = 0.003591976735924036
GMRes iteration 22, residual = 0.002994326308323427
GMRes iteration 23, residual = 0.002946930794913705
GMRes iteration 24, residual = 0.001979268009256806
GMRes iteration 25, residual = 0.0019408666098340928
GMRes iteration 26, residual = 0.0015575496265593243
GMRes iteration 27, residual = 0.0014246231844533505
GMRes iteration 28, residual = 0.0012204425765592055
GMRes iteration 29, residual = 0.0010847044870138376
GMRes iteration 30, residual = 0.0009837866451021156
GMRes iteration 31, residual = 0.0007490403198617734
GMRes iteration 32, residual = 0.0007488562559548595
GMRes iteration 33, residual = 0.0006134911148998155
GMRes iteration 34, residual = 0.0006124099757302618
GMRes iteration 35, residual = 0.0004939637916748199
GMRes iteration 36, residual = 0.00047858684476626115
GMRes iteration 37, residual = 0.00037174196502167306
GMRes iteration 38, residual = 0.0003709705652678626
GMRes iteration 39, residual = 0.00031421136012430077
GMRes iteration 40, residual = 0.0003104027602432238
GMRes iteration 41, residual = 0.00021837369378423065
GMRes iteration 42, residual = 0.00021449537656424627
GMRes iteration 43, residual = 0.00019217140387005894
GMRes iteration 44, residual = 0.0001788029673000178
GMRes iteration 45, residual = 0.0001623610281205856
GMRes iteration 46, residual = 0.00013935861769303966
GMRes iteration 47, residual = 0.00013877792368833944
GMRes iteration 48, residual = 0.00011238312296675299
GMRes iteration 49, residual = 0.00011226444042012831
GMRes iteration 50, residual = 8.692527619920691e-05
GMRes iteration 51, residual = 8.664936049783665e-05
GMRes iteration 52, residual = 7.291910106071894e-05
GMRes iteration 53, residual = 7.253481554861397e-05
GMRes iteration 54, residual = 5.70476438134183e-05
GMRes iteration 55, residual = 4.796453972807561e-05
GMRes iteration 56, residual = 4.672528342550095e-05
GMRes iteration 57, residual = 3.819921676798488e-05
GMRes iteration 58, residual = 3.807472571989822e-05
GMRes iteration 59, residual = 2.8152152672954432e-05
GMRes iteration 60, residual = 2.7724898193442326e-05
GMRes iteration 61, residual = 2.4337896037500556e-05
GMRes iteration 62, residual = 2.4273010198833463e-05
GMRes iteration 63, residual = 2.0368339155673966e-05
GMRes iteration 64, residual = 1.987311148855257e-05
GMRes iteration 65, residual = 1.6038071007721704e-05
GMRes iteration 66, residual = 1.4477033113576382e-05
GMRes iteration 67, residual = 1.4163714626545746e-05
GMRes iteration 68, residual = 1.040982180201345e-05
GMRes iteration 69, residual = 1.0328819439526262e-05
GMRes iteration 70, residual = 8.190115548042825e-06
GMRes iteration 71, residual = 7.817839190329928e-06
GMRes iteration 72, residual = 7.018534795812861e-06
GMRes iteration 73, residual = 5.897813610819529e-06
GMRes iteration 74, residual = 5.545469612552916e-06
GMRes iteration 75, residual = 4.159611582701412e-06
GMRes iteration 76, residual = 4.1288141343843395e-06
GMRes iteration 77, residual = 3.6120195126236136e-06
GMRes iteration 78, residual = 3.351576956733913e-06
GMRes iteration 79, residual = 2.984680855264741e-06
GMRes iteration 80, residual = 2.562703080702694e-06
GMRes iteration 81, residual = 2.4522388343613752e-06
GMRes iteration 82, residual = 2.1344289454240862e-06
GMRes iteration 83, residual = 2.0774997652487286e-06
GMRes iteration 84, residual = 1.6116740868794285e-06
GMRes iteration 85, residual = 1.5868599920985996e-06
GMRes iteration 86, residual = 1.348775717628033e-06
GMRes iteration 87, residual = 1.2296459395447935e-06
GMRes iteration 88, residual = 1.184705017169241e-06
GMRes iteration 89, residual = 9.075407480055652e-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:

[ ]: