Dirichlet Laplace Direct Method

Dirichlet Laplace Direct Method#

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

from netgen.occ import *
from ngsolve import *
from ngsolve.webgui import Draw
from ngbem import *
from ngsolve import Projector, Preconditioner
from ngsolve.krylovspace import CG
Loading ngbem library

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:

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:

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 =  20965

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

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)
BaseWebGuiScene

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:

u1 = GridFunction(fesL2)
pre = BilinearForm(u*v*ds, diagonal=True).Assemble().mat.Inverse()
with TaskManager(): 
    V = SingleLayerPotentialOperator(fesL2, intorder=12, leafsize=40, eta=3., eps=1e-4, method="aca", testhmatrix=False)

    M = BilinearForm( uH1 * v.Trace() * ds(bonus_intorder=3)).Assemble()
    K = DoubleLayerPotentialOperator(fesH1, fesL2, intorder=12, leafsize=40, eta=3., eps=1e-4, method="aca")
    
    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=200, initialize=False, printrates=True)

Draw (u1, mesh, draw_vol=False, order=3);
CG iteration 1, residual = 0.2688253801037298     
CG iteration 2, residual = 0.06553558458073683     
CG iteration 3, residual = 0.015544599077948492     
CG iteration 4, residual = 0.004632058606485141     
CG iteration 5, residual = 0.004593036294323392     
CG iteration 6, residual = 0.0008022189950998139     
CG iteration 7, residual = 0.0001976604628421484     
CG iteration 8, residual = 7.158924638153832e-05     
CG iteration 9, residual = 0.0001030320559234788     
CG iteration 10, residual = 7.699909626520429e-05     
CG iteration 11, residual = 3.150255444331728e-05     
CG iteration 12, residual = 1.8228831474636923e-05     
CG iteration 13, residual = 1.5343515082211592e-05     
CG iteration 14, residual = 1.7402177351023794e-05     
CG iteration 15, residual = 6.589075458617313e-06     
CG iteration 16, residual = 5.802384354063657e-06     
CG iteration 17, residual = 1.1118866472429553e-05     
CG iteration 18, residual = 4.963671995981006e-06     
CG iteration 19, residual = 6.692547041312336e-06     
CG iteration 20, residual = 4.667024702442362e-06     
CG iteration 21, residual = 2.183834188669911e-06     
CG iteration 22, residual = 1.4625200657010096e-06     
CG iteration 23, residual = 1.7967828880618047e-06     
CG iteration 24, residual = 1.7350325927622153e-06     
CG iteration 25, residual = 1.5891232841504342e-06     
CG iteration 26, residual = 1.0526815232099118e-06     
CG iteration 27, residual = 8.878557830179139e-07     
CG iteration 28, residual = 4.259435052307448e-07     
CG iteration 29, residual = 3.804753677114846e-07     
CG iteration 30, residual = 3.831943458012857e-07     
CG iteration 31, residual = 3.5379552308248864e-07     
CG iteration 32, residual = 3.4974492682269254e-07     
CG iteration 33, residual = 1.950002222287894e-07     
CG iteration 34, residual = 2.095693642989858e-07     
CG iteration 35, residual = 3.8584789038963966e-07     
CG iteration 36, residual = 1.4045168299806663e-07     
CG iteration 37, residual = 1.1182182402974015e-07     
CG iteration 38, residual = 5.517831213747081e-08     
CG iteration 39, residual = 9.929492891981467e-08     
CG iteration 40, residual = 6.762790390076364e-08     
CG iteration 41, residual = 5.366845377493877e-08     
CG iteration 42, residual = 3.1798784296967365e-08     
CG iteration 43, residual = 2.574529250934904e-08     
CG iteration 44, residual = 2.2837643129812035e-08     
CG iteration 45, residual = 3.6218388823949796e-08     
CG iteration 46, residual = 2.893687540266911e-08     
CG iteration 47, residual = 2.1352335236873292e-08     
CG iteration 48, residual = 1.0405930747743728e-08     
CG iteration 49, residual = 9.194248692395103e-09     
CG iteration 50, residual = 8.885071352642902e-09     
CG iteration 51, residual = 1.2971281650331245e-08     
CG iteration 52, residual = 5.514226106365831e-09     
CG iteration 53, residual = 6.345431815568637e-09     
CG iteration 54, residual = 3.1571079297916193e-09     
CG iteration 55, residual = 3.3915831150085376e-09     
CG iteration 56, residual = 4.196536341885275e-09     
CG iteration 57, residual = 2.36465715985326e-09     
for t in Timers():
    if "ngbem" in t["name"]:
        print (t)
{'name': 'ngbem - HMatrix::MultAdd', 'time': 0.13537381696893894, 'counts': 58, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem disjoint LaplaceDL', 'time': 2.7198759684103235, 'counts': 35954, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem common edge LaplaceDL', 'time': 10.260364542893987, 'counts': 1505, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem common vertex LaplaceDL', 'time': 12.188614705405996, 'counts': 4516, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem identic panel LaplaceDL', 'time': 4.067725101261571, 'counts': 501, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem FarFieldBlock - matvec2 LaplaceDL', 'time': 0.04467998658107885, 'counts': 5210, 'flops': 73607065.0, 'Gflop/s': 1.647428091018524}
{'name': 'ngbem FarFieldBlock - matvec1 LaplaceDL', 'time': 0.02620549503317173, 'counts': 5210, 'flops': 143971849.0, 'Gflop/s': 5.493956470494298}
{'name': 'ngbem FarFieldBlock - norm LaplaceDL', 'time': 0.2461732197505629, 'counts': 5210, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem - loops LaplaceDL', 'time': 29.319149350539437, 'counts': 739, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem FarFieldBlock - kernel LaplaceDL', 'time': 0.06524628595859225, 'counts': 10420, 'flops': 38852590.0, 'Gflop/s': 0.5954758869287566}
{'name': 'ngbem FarFieldBlock - ACA LaplaceDL', 'time': 0.5918590755041725, 'counts': 447, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem FarFieldBlock LaplaceDL', 'time': 1.1528825223453434, 'counts': 536, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem - all LaplaceDL', 'time': 30.27700374277697, 'counts': 739, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem disjoint LaplaceSL', 'time': 0.5176561936534021, 'counts': 6925, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem common edge LaplaceSL', 'time': 3.0531155271977166, 'counts': 532, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem common vertex LaplaceSL', 'time': 3.6576357577933054, 'counts': 1558, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem identic panel LaplaceSL', 'time': 1.2044378126764137, 'counts': 175, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem FarFieldBlock - matvec2 LaplaceSL', 'time': 0.03841715348942381, 'counts': 5557, 'flops': 60723789.0, 'Gflop/s': 1.5806425902095318}
{'name': 'ngbem FarFieldBlock - matvec1 LaplaceSL', 'time': 0.010158820965625905, 'counts': 5557, 'flops': 60619076.0, 'Gflop/s': 5.9671369546835145}
{'name': 'ngbem FarFieldBlock - norm LaplaceSL', 'time': 0.13202063303891864, 'counts': 5557, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem FarFieldBlock - kernel LaplaceSL', 'time': 0.029408960295076454, 'counts': 11114, 'flops': 25833388.0, 'Gflop/s': 0.8784189492181721}
{'name': 'ngbem FarFieldBlock - ACA LaplaceSL', 'time': 0.34680290399932384, 'counts': 548, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem FarFieldBlock LaplaceSL', 'time': 0.7455792471814339, 'counts': 635, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem - loops LaplaceSL', 'time': 8.445176257368969, 'counts': 669, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem - all LaplaceSL', 'time': 9.262878026834052, 'counts': 669, 'flops': 0.0, 'Gflop/s': 0.0}
{'name': 'ngbem - BaseClass::CalcHMatrix', 'time': 11.906686595789036, 'counts': 2, 'flops': 0.0, 'Gflop/s': 0.0}

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

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.0021893572396633047