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.94329927058774
GMRes iteration 2, residual = 10.547738066512279
GMRes iteration 3, residual = 2.6122774561838593
GMRes iteration 4, residual = 1.9193598745774196
GMRes iteration 5, residual = 0.4151312183469935
GMRes iteration 6, residual = 0.3904448759672413
GMRes iteration 7, residual = 0.17719617347188396
GMRes iteration 8, residual = 0.13791748320195735
GMRes iteration 9, residual = 0.08865758874505879
GMRes iteration 10, residual = 0.04770667025278976
GMRes iteration 11, residual = 0.04662105056009847
GMRes iteration 12, residual = 0.04497121005324146
GMRes iteration 13, residual = 0.020803795188892086
GMRes iteration 14, residual = 0.013029442692585769
GMRes iteration 15, residual = 0.012428826178702933
GMRes iteration 16, residual = 0.007375353267668423
GMRes iteration 17, residual = 0.007328660477955353
GMRes iteration 18, residual = 0.006116844305180282
GMRes iteration 19, residual = 0.00557474073529345
GMRes iteration 20, residual = 0.003669981046576089
GMRes iteration 21, residual = 0.0035919739535645295
GMRes iteration 22, residual = 0.0029943263110651934
GMRes iteration 23, residual = 0.0029469272197431306
GMRes iteration 24, residual = 0.0019792658964475243
GMRes iteration 25, residual = 0.0019408671542832233
GMRes iteration 26, residual = 0.0015575454887884002
GMRes iteration 27, residual = 0.0014246225478179583
GMRes iteration 28, residual = 0.0012204383855294655
GMRes iteration 29, residual = 0.0010847049141660327
GMRes iteration 30, residual = 0.0009837817135757993
GMRes iteration 31, residual = 0.0007490390345336671
GMRes iteration 32, residual = 0.0007488551551952959
GMRes iteration 33, residual = 0.0006134908290596645
GMRes iteration 34, residual = 0.0006124094630833038
GMRes iteration 35, residual = 0.0004939647433187051
GMRes iteration 36, residual = 0.00047858769696751847
GMRes iteration 37, residual = 0.000371742951719414
GMRes iteration 38, residual = 0.00037097165102282077
GMRes iteration 39, residual = 0.0003142122701058231
GMRes iteration 40, residual = 0.00031040388522289456
GMRes iteration 41, residual = 0.00021837440921558586
GMRes iteration 42, residual = 0.0002144958223692561
GMRes iteration 43, residual = 0.0001921722093089571
GMRes iteration 44, residual = 0.00017880328469099637
GMRes iteration 45, residual = 0.00016236176272092158
GMRes iteration 46, residual = 0.00013935888729775054
GMRes iteration 47, residual = 0.00013877833727816004
GMRes iteration 48, residual = 0.00011238342651392257
GMRes iteration 49, residual = 0.00011226469228306348
GMRes iteration 50, residual = 8.692538273698849e-05
GMRes iteration 51, residual = 8.664938785098507e-05
GMRes iteration 52, residual = 7.291921147722146e-05
GMRes iteration 53, residual = 7.253483209598842e-05
GMRes iteration 54, residual = 5.704785240657879e-05
GMRes iteration 55, residual = 4.7964478625381804e-05
GMRes iteration 56, residual = 4.672528369133715e-05
GMRes iteration 57, residual = 3.819918644044966e-05
GMRes iteration 58, residual = 3.807470953304845e-05
GMRes iteration 59, residual = 2.8152174055890493e-05
GMRes iteration 60, residual = 2.7724896952249298e-05
GMRes iteration 61, residual = 2.433791253472575e-05
GMRes iteration 62, residual = 2.4273022230816905e-05
GMRes iteration 63, residual = 2.036834132007951e-05
GMRes iteration 64, residual = 1.9873104646594464e-05
GMRes iteration 65, residual = 1.603804673832704e-05
GMRes iteration 66, residual = 1.4477030905721276e-05
GMRes iteration 67, residual = 1.416369841157829e-05
GMRes iteration 68, residual = 1.0409818001997474e-05
GMRes iteration 69, residual = 1.0328820158165103e-05
GMRes iteration 70, residual = 8.190112235838925e-06
GMRes iteration 71, residual = 7.817841977489779e-06
GMRes iteration 72, residual = 7.018534324588462e-06
GMRes iteration 73, residual = 5.89781774707014e-06
GMRes iteration 74, residual = 5.545474210252832e-06
GMRes iteration 75, residual = 4.159615098864856e-06
GMRes iteration 76, residual = 4.128816004241347e-06
GMRes iteration 77, residual = 3.612024109547303e-06
GMRes iteration 78, residual = 3.3515738417511477e-06
GMRes iteration 79, residual = 2.9846858915626448e-06
GMRes iteration 80, residual = 2.5627047881882396e-06
GMRes iteration 81, residual = 2.4522443845632054e-06
GMRes iteration 82, residual = 2.134430987061428e-06
GMRes iteration 83, residual = 2.0775029275906883e-06
GMRes iteration 84, residual = 1.6116737372138804e-06
GMRes iteration 85, residual = 1.586859613194713e-06
GMRes iteration 86, residual = 1.3487753405374731e-06
GMRes iteration 87, residual = 1.2296466211567543e-06
GMRes iteration 88, residual = 1.1847314161527659e-06
GMRes iteration 89, residual = 9.075571370960842e-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:

[ ]: