This page was generated from unit-11.2-bem-Laplace/Laplace_DtN_direct.ipynb.

11.2.2. Dirichlet Laplace Direct Method

keys: homogeneous Dirichlet bvp, double and single layer potential, unknown Neumann data

[1]:
from netgen.occ import *
from ngsolve import *
from ngsolve.webgui import Draw
from ngsolve.bem import *
from ngsolve.solvers import CG

Consider the Dirichlet boundary value problem

\[\begin{split}\left\{ \begin{array}{rcl l} \Delta u &=& 0, \quad &\Omega \subset \mathbb R^3\,,\\ \gamma_0 u&=& u_0, \quad &\Gamma = \partial \Omega\,.\end{array} \right.\end{split}\]

Let us choose the following ansatz for the solution \(u\in H^1(\Omega)\) (direct ansatz)

\[u(x) = \underbrace{ \int\limits_\Gamma \displaystyle{\frac{1}{4\,\pi}\, \frac{1}{\| x-y\|} } \, u_1(y)\, \mathrm{d}\sigma_y}_{\displaystyle{ \mathrm{ SL}(u_1) }} - \underbrace{ \int\limits_\Gamma \displaystyle{\frac{1}{4\,\pi}\, \frac{\langle n(y) , x-y\rangle }{\| x-y\|^3} } \, u_0(y)\, \mathrm{d}\sigma_y}_{\displaystyle{ \mathrm{DL}(u_0) }}\]

and solve for the Neumann data \(u_1 \in H^{-\frac12}(\Gamma)\) by the boundary element method, i.e.,

\[\forall \, v\in H^{-\frac12}(\Gamma): \quad \left\langle \gamma_0 \left(\mathrm{SL}(u_1)\right), v \right\rangle_{-\frac12}= \left\langle u_0, v\right\rangle_{-\frac12} + \left\langle \gamma_0 \left(\mathrm{DL}(u_0)\right), v\right\rangle_{-\frac12}\,.\]

Define the geometry \(\Omega \subset \mathbb R^3\) and create a mesh:

[2]:
sp = Sphere( (0,0,0), 1)
mesh = Mesh( OCCGeometry(sp).GenerateMesh(maxh=0.2)).Curve(4)

Create the finite element spaces for \(H^{-\frac12}(\Gamma)\) and \(H^{\frac12}(\Gamma)\) according to the given mesh:

[3]:
fesL2 = SurfaceL2(mesh, order=3, dual_mapping=True)
u,v = fesL2.TnT()
fesH1 = H1(mesh, order=4)
uH1,vH1 = fesH1.TnT()
print ("ndofL2 = ", fesL2.ndof, "ndof H1 = ", fesH1.ndof)
ndofL2 =  7220 ndof H1 =  21241

Compute the interpolation of exact Dirichlet data \(u_0\) in finite element space:

[4]:
uexa = 1/ sqrt( (x-1)**2 + (y-1)**2 + (z-1)**2 )
u0 = GridFunction(fesH1)
u0.Interpolate (uexa)
Draw (u0, mesh, draw_vol=False, order=3);

The discretisation of the above variational formulation leads to a system of linear equations, ie

\[\mathrm{V} \, \mathrm{u}_1 = \left( \frac12 \,\mathrm{M} + \mathrm{K} \right) \, \mathrm{u}_0\,,\]

where the linear operators are as follows

  • \(\mathrm{V}\) is the single layer operator. \(\mathrm V\) is regular and symmetric.

  • \(\mathrm{M}\) is a mass matrix.

  • \(\mathrm{K}\) is the double layer operator.

We approximate the linear operators as hmatrices and solve for the Neumann data \(u_1\) with an iterative solver:

[5]:
u1 = GridFunction(fesL2)
pre = BilinearForm(u*v*ds, diagonal=True).Assemble().mat.Inverse()

with TaskManager():
    V = SingleLayerPotentialOperator(fesL2, intorder=12)

    M = BilinearForm( uH1 * v.Trace() * ds(bonus_intorder=3)).Assemble()
    K = DoubleLayerPotentialOperator(fesH1, fesL2, intorder=12)

    rhs = ( (0.5 * M.mat + K.mat)*u0.vec).Evaluate()
    CG(mat = V.mat, pre=pre, rhs = rhs, sol=u1.vec, tol=1e-8, maxsteps=50, initialize=False, printrates=True)
CG iteration 1, residual = 0.26882443508335313
CG iteration 2, residual = 0.0655351773816169
CG iteration 3, residual = 0.0155370066463854
CG iteration 4, residual = 0.0036775507940641212
CG iteration 5, residual = 0.0009473297864010314
CG iteration 6, residual = 0.0013803247805890552
CG iteration 7, residual = 0.00019659306334787909
CG iteration 8, residual = 5.781718854672164e-05
CG iteration 9, residual = 4.6435183675259584e-05
CG iteration 10, residual = 3.399771902764883e-05
CG iteration 11, residual = 1.7058052474862087e-05
CG iteration 12, residual = 9.988909715217277e-06
CG iteration 13, residual = 4.6690621552090545e-06
CG iteration 14, residual = 2.487873868095676e-06
CG iteration 15, residual = 5.520638665054005e-06
CG iteration 16, residual = 2.092351512206722e-06
CG iteration 17, residual = 1.9107086755868024e-06
CG iteration 18, residual = 1.7874571284593366e-06
CG iteration 19, residual = 9.145440297840876e-07
CG iteration 20, residual = 6.252906065869179e-07
CG iteration 21, residual = 5.046629419244661e-07
CG iteration 22, residual = 3.8856766260962595e-07
CG iteration 23, residual = 3.7945837805127596e-07
CG iteration 24, residual = 8.134055124988375e-07
CG iteration 25, residual = 2.09118902442277e-07
CG iteration 26, residual = 1.008877969611272e-07
CG iteration 27, residual = 8.446743073574429e-08
CG iteration 28, residual = 5.0345338319179604e-08
CG iteration 29, residual = 5.520760455913188e-08
CG iteration 30, residual = 2.678769121938202e-08
CG iteration 31, residual = 5.1129528647339e-08
CG iteration 32, residual = 3.460406667119261e-08
CG iteration 33, residual = 2.581844133111495e-08
CG iteration 34, residual = 2.5905766307499332e-08
CG iteration 35, residual = 1.3211032628548991e-08
CG iteration 36, residual = 7.520160045252428e-09
CG iteration 37, residual = 4.500493978520459e-09
CG iteration 38, residual = 4.859495648427841e-09
CG iteration 39, residual = 3.0818496614806833e-09
CG iteration 40, residual = 3.545134522555692e-09
CG iteration 41, residual = 5.557912254344217e-09
CG iteration 42, residual = 2.067215289960801e-09
[6]:
Draw (u1, mesh, draw_vol=False, order=3);

Let’s have a look at the exact Neumann data and compute the error of the numerical solution:

[7]:
graduexa = CF( (uexa.Diff(x), uexa.Diff(y), uexa.Diff(z)) )
n = specialcf.normal(3)
u1exa = graduexa*n
Draw (u1exa, mesh, draw_vol=False, order=3);
print ("L2-error =", sqrt (Integrate ( (u1exa-u1)**2, mesh, BND)))
L2-error = 0.00024938694903869495
[ ]:

[ ]: