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

11.2.1 Dirichlet Laplace Indirect MethodΒΆ

keys: homogeneous Dirichlet bvp, single layer potential, ACA

[1]:
from netgen.occ import *
from ngsolve import *
from ngsolve.webgui import Draw
from ngsolve.bem import *
from ngsolve import Projector, Preconditioner
from ngsolve.krylovspace import CG

We 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)\) (indirect ansatz)

\[u(x) = \underbrace{ \int\limits_\Gamma \displaystyle{\frac{1}{4\,\pi}\, \frac{1}{\| x-y\|} } \, j(y)\, \mathrm{d}\sigma_y }_{\displaystyle{ \mathrm{SL}(j) } }\]

and solve for the density \(j\in H^{-\frac12}(\Gamma)\) by the boundary element method, i.e. the numerical solution of the variational formulation

\[\forall \, v\in H^{-\frac12}(\Gamma): \quad \left\langle \gamma_0 \left(\mathrm{SL}(j)\right), v \right\rangle_{-\frac12} = \left\langle u_0, 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.3)).Curve(4)

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

[3]:
fesL2 = SurfaceL2(mesh, order=4, dual_mapping=True)
u,v = fesL2.TnT()

Define Dirichlet data \(u_0\) and compute the right hand side vector:

[4]:
u0 = 1/ sqrt( (x-1)**2 + (y-1)**2 + (z-1)**2 )
rhs = LinearForm (u0*v.Trace()*ds(bonus_intorder=3)).Assemble()

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

\[\mathrm{V} \, \mathrm{j} = \mathrm{rhs} \,,\]

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

Demo 1: Assemble the single layer operator \(V\) as dense matrix and solve for unknwon density \(j\):

[5]:
j = GridFunction(fesL2)
pre = BilinearForm(u*v*ds, diagonal=True).Assemble().mat.Inverse()
with TaskManager():
    # V = SingleLayerPotentialOperator(fesL2, intorder=10)
    V = LaplaceSL(u*ds)*v*ds
    CG(mat = V.mat, pre=pre, rhs = rhs.vec, sol=j.vec, tol=1e-8, maxsteps=200, initialize=False, printrates=True)
CG iteration 1, residual = 2.185606954606291
CG iteration 2, residual = 0.5510990298605162
CG iteration 3, residual = 0.1311274118731144
CG iteration 4, residual = 0.03107892126080938
CG iteration 5, residual = 0.00737248197512428
CG iteration 6, residual = 0.0017929371359740698
CG iteration 7, residual = 0.0005259800974227988
CG iteration 8, residual = 0.00025807777882633784
CG iteration 9, residual = 0.0002832379521244805
CG iteration 10, residual = 0.00012036653419387316
CG iteration 11, residual = 0.0002350720197571702
CG iteration 12, residual = 4.492369025523945e-05
CG iteration 13, residual = 4.801474718476479e-05
CG iteration 14, residual = 2.782522946270083e-05
CG iteration 15, residual = 1.462780519504024e-05
CG iteration 16, residual = 1.4781005108292774e-05
CG iteration 17, residual = 9.771555034947063e-06
CG iteration 18, residual = 6.175958029417542e-06
CG iteration 19, residual = 1.5387382836454795e-05
CG iteration 20, residual = 3.2801500430824397e-06
CG iteration 21, residual = 2.970067358691875e-06
CG iteration 22, residual = 3.493169500410065e-06
CG iteration 23, residual = 2.593798073962397e-06
CG iteration 24, residual = 1.260433272980252e-06
CG iteration 25, residual = 8.372904985691209e-07
CG iteration 26, residual = 1.2440365951762135e-06
CG iteration 27, residual = 5.872586615825374e-07
CG iteration 28, residual = 4.4146186934114627e-07
CG iteration 29, residual = 5.781842321683408e-07
CG iteration 30, residual = 3.526775208042795e-07
CG iteration 31, residual = 1.778225531414204e-07
CG iteration 32, residual = 1.544195099290323e-07
CG iteration 33, residual = 2.1911322969691716e-07
CG iteration 34, residual = 1.2550255964749062e-07
CG iteration 35, residual = 8.72380230193785e-08
CG iteration 36, residual = 6.888878583827084e-08
CG iteration 37, residual = 1.0272145305453892e-07
CG iteration 38, residual = 3.4372701170999706e-08
CG iteration 39, residual = 4.024003368369137e-08
CG iteration 40, residual = 5.732884935783174e-08
CG iteration 41, residual = 1.854725022228517e-08
[6]:
Draw (j, order=3);

Notes:

[ ]: