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.94329927058897
GMRes iteration 2, residual = 10.547738043195517
GMRes iteration 3, residual = 2.6122774434372857
GMRes iteration 4, residual = 1.9193598549952162
GMRes iteration 5, residual = 0.4151312149373696
GMRes iteration 6, residual = 0.39044487452556176
GMRes iteration 7, residual = 0.17719618484116878
GMRes iteration 8, residual = 0.13791749808128562
GMRes iteration 9, residual = 0.0886576441869601
GMRes iteration 10, residual = 0.04770670569340013
GMRes iteration 11, residual = 0.04662109255459429
GMRes iteration 12, residual = 0.04497123877062714
GMRes iteration 13, residual = 0.02080387251300668
GMRes iteration 14, residual = 0.013029447152683206
GMRes iteration 15, residual = 0.012428849702072415
GMRes iteration 16, residual = 0.007375357812442446
GMRes iteration 17, residual = 0.007328667494073906
GMRes iteration 18, residual = 0.0061168451041507576
GMRes iteration 19, residual = 0.0055747512053160335
GMRes iteration 20, residual = 0.003669979651737646
GMRes iteration 21, residual = 0.003591976741851272
GMRes iteration 22, residual = 0.0029943263054992945
GMRes iteration 23, residual = 0.0029469307914403533
GMRes iteration 24, residual = 0.001979268008635592
GMRes iteration 25, residual = 0.001940866609689945
GMRes iteration 26, residual = 0.001557549634432896
GMRes iteration 27, residual = 0.0014246231874214305
GMRes iteration 28, residual = 0.001220442580724697
GMRes iteration 29, residual = 0.001084704487555052
GMRes iteration 30, residual = 0.0009837866441411156
GMRes iteration 31, residual = 0.0007490403244690277
GMRes iteration 32, residual = 0.0007488562613459386
GMRes iteration 33, residual = 0.0006134912663569409
GMRes iteration 34, residual = 0.0006124100914887338
GMRes iteration 35, residual = 0.0004939639241873984
GMRes iteration 36, residual = 0.000478586895811197
GMRes iteration 37, residual = 0.0003717419668553123
GMRes iteration 38, residual = 0.00037097056706618934
GMRes iteration 39, residual = 0.00031421136026626794
GMRes iteration 40, residual = 0.00031040276057089537
GMRes iteration 41, residual = 0.0002183736945636929
GMRes iteration 42, residual = 0.0002144953769022915
GMRes iteration 43, residual = 0.0001921714046511171
GMRes iteration 44, residual = 0.00017880296698568235
GMRes iteration 45, residual = 0.0001623610279743289
GMRes iteration 46, residual = 0.00013935861762850427
GMRes iteration 47, residual = 0.00013877792354042952
GMRes iteration 48, residual = 0.00011238312338471435
GMRes iteration 49, residual = 0.00011226444088038565
GMRes iteration 50, residual = 8.69252764752689e-05
GMRes iteration 51, residual = 8.664936080910505e-05
GMRes iteration 52, residual = 7.291910129353608e-05
GMRes iteration 53, residual = 7.253481582355446e-05
GMRes iteration 54, residual = 5.704764412365262e-05
GMRes iteration 55, residual = 4.796454019736064e-05
GMRes iteration 56, residual = 4.672528374965699e-05
GMRes iteration 57, residual = 3.819921756294961e-05
GMRes iteration 58, residual = 3.80747264129586e-05
GMRes iteration 59, residual = 2.815215348182639e-05
GMRes iteration 60, residual = 2.7724899588278344e-05
GMRes iteration 61, residual = 2.4337896877127376e-05
GMRes iteration 62, residual = 2.4273011365721145e-05
GMRes iteration 63, residual = 2.0368339655558123e-05
GMRes iteration 64, residual = 1.9873112326157423e-05
GMRes iteration 65, residual = 1.603807117724699e-05
GMRes iteration 66, residual = 1.4477033212955939e-05
GMRes iteration 67, residual = 1.4163714744644644e-05
GMRes iteration 68, residual = 1.0409822069017602e-05
GMRes iteration 69, residual = 1.0328819776430524e-05
GMRes iteration 70, residual = 8.190115576600966e-06
GMRes iteration 71, residual = 7.8178392083996e-06
GMRes iteration 72, residual = 7.0185348057855195e-06
GMRes iteration 73, residual = 5.897813784724025e-06
GMRes iteration 74, residual = 5.545469691992986e-06
GMRes iteration 75, residual = 4.1596115999647055e-06
GMRes iteration 76, residual = 4.128814156338762e-06
GMRes iteration 77, residual = 3.6120195409046164e-06
GMRes iteration 78, residual = 3.3515770498487764e-06
GMRes iteration 79, residual = 2.9846808707118993e-06
GMRes iteration 80, residual = 2.5627031190500025e-06
GMRes iteration 81, residual = 2.4522388601623355e-06
GMRes iteration 82, residual = 2.134428968536966e-06
GMRes iteration 83, residual = 2.077499782281205e-06
GMRes iteration 84, residual = 1.6116741250398705e-06
GMRes iteration 85, residual = 1.5868600438146168e-06
GMRes iteration 86, residual = 1.3487757094664411e-06
GMRes iteration 87, residual = 1.2296459407180102e-06
GMRes iteration 88, residual = 1.1847050133225304e-06
GMRes iteration 89, residual = 9.075407361217008e-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:

[ ]: