Neumann Laplace Direct Method

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

keys: homogeneous Neumann bvp, hypersingular operator, unknown Dirichlet data

Neumann Laplace Direct Method#

Consider the Neumann boundary value problem

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

The solution \(u\in H^1(\Omega)\) is given by

\[ (1) \quad \quad \quad u(x) = \underbrace{ \int\limits_\Gamma \displaystyle{\frac{1}{4\,\pi}\, \frac{1}{\| x-y\|} } \, j(y)\, \mathrm{d}s_y}_{\displaystyle{ \mathrm{LaplaceSL}(j) }} - \underbrace{ \int\limits_\Gamma \displaystyle{\frac{1}{4\,\pi}\, \frac{\langle n(y), x-y \rangle}{\| x-y\|^3} } \, m( y)\, \mathrm{d}s_y}_{\displaystyle{ \mathrm{LaplaceDL}(m) }}\,.\]

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

sp = Sphere( (0,0,0), 1)
mesh = Mesh( OCCGeometry(sp).GenerateMesh(maxh=0.2)).Curve(2)
fesL2 = SurfaceL2(mesh, order=2, dual_mapping=False) 
u,v = fesL2.TnT()
fesH1 = H1(mesh, order=2, definedon=mesh.Boundaries(".*"))
uH1,vH1 = fesH1.TnT()
print ("ndofL2 = ", fesL2.ndof, "ndof H1 = ", fesH1.ndof)
ndofL2 =  4332 ndof H1 =  1580

Consider as Neumann data \(m = n\cdot \nabla \displaystyle \frac{1}{\|x-x_0\|}\) and compute its interpolation in \(H^{-\frac12}(\Gamma)\):

uexa = 1/ sqrt( (x-1)**2 + (y-1)**2 + (z-1)**2 )
graduexa = CF( (uexa.Diff(x), uexa.Diff(y), uexa.Diff(z)) )
n = specialcf.normal(3)
j = GridFunction(fesL2)
j.Interpolate(graduexa*n, definedon=mesh.Boundaries(".*"))
Draw (j, mesh, draw_vol=False, order=3);

Boundary Integral Equation

Let’s carefully apply the Neumann trace to (1) and derive a boundary integral equation for \(m\):

\[\begin{split}\begin{array}{r c ccc} \forall \, v\in H^{\frac12}(\Gamma): && \displaystyle \int\limits_\Gamma v \cdot \gamma_1 \left(\mathrm{DL}(m) \right) \, \mathrm d s &=& \displaystyle \int\limits_\Gamma j \cdot v \,\mathrm d s &-& \displaystyle \int\limits_\Gamma v \cdot \gamma_1 \left(\mathrm{SL}(j)\right) \, \mathrm d s \\[2ex] & \Rightarrow & \mathrm{D}\,\mathrm{m} &=& \dfrac12 \, \mathrm{M} \, \mathrm j &-& \mathrm{K}' \, \mathrm j\,. \end{array} \end{split}\]

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

\[ \left( \mathrm{D} + \mathrm{S}\right) \, \mathrm{m} = \left( \frac12 \,\mathrm{M} - \mathrm{K}' \right) \, \mathrm{j} \,, \]

where the linear operators are as follows

  • \(\mathrm{D} + \mathrm S\) is regular and symmetric, it is the stabilized hypersingular operator (see Laplace NtD indirect )

  • \(\mathrm{K}'\) is the adjoint double layer potential operator, it is the transposed of \(\mathrm K\) as discussed in Laplace DtN direct).

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

vH1m1 = LinearForm(vH1*1*ds(bonus_intorder=3)).Assemble()
S = (BaseMatrix(Matrix(vH1m1.vec.Reshape(1))))@(BaseMatrix(Matrix(vH1m1.vec.Reshape(fesH1.ndof))))

The assembly of the hypersingular operator is declarative and close to their actual implementations (explicit formula are given below).

with TaskManager(): 
    D = LaplaceSL(curl(uH1)*ds)*curl(vH1)*ds
    K = LaplaceDL(uH1*ds)*v*ds
    M = BilinearForm( uH1 * v * ds(bonus_intorder=3)).Assemble()

Left and right hand side of linear system are computed and the solution is the Dirichlet data \(u_0\):

m = GridFunction(fesH1)
pre = BilinearForm(uH1*vH1*ds).Assemble().mat.Inverse(freedofs=fesH1.FreeDofs())
with TaskManager(): 
    lhs = D.mat + S
    rhs = ( (0.5 * M.mat.T - K.mat.T) * j.vec ).Evaluate()
    CG(mat=lhs, pre=pre, rhs=rhs, sol=m.vec, tol=1e-8, 
       maxsteps=200, initialize=False, printrates=False)

Let’s have a look at the numerical and exact Dirichlet data and compare them quantitatively:

Draw (m, mesh, draw_vol=False, order=3);
mexa = GridFunction(fesH1)
mexa.Interpolate (uexa, definedon=mesh.Boundaries(".*"))
Draw (mexa, mesh, draw_vol=False, order=3);
print ("L2 error of surface gradients =", 
       sqrt (Integrate ( (grad(mexa)-grad(m))**2, mesh, BND)))
L2 error of surface gradients = 0.004568783275440319

Evaluation of the Solution

screen = WorkPlane(Axes( (0,0,0), Z, X)).RectangleC(0.5,0.5).Face()
screen.faces.name="screen"
vismesh = screen.GenerateMesh(maxh=0.05)
LSPotential = LaplaceSL( u*ds )
LDPotential = LaplaceDL( uH1*ds )
repformula = LSPotential(j) - LDPotential(m)
fes_screen = H1(vismesh, order=5)
gf_screen = GridFunction(fes_screen)
with TaskManager():
    gf_screen.Set(repformula, definedon=vismesh.Boundaries("screen"))
Draw(gf_screen)
 
BaseWebGuiScene

Have a look at the error on the screen (it must be a constant):

Draw(repformula - uexa, vismesh) # constant
BaseWebGuiScene

References:

For details on the explicit representation of the hypersingular operator have a look into Numerische Näherungsverfahren für elliptische Randwertprobleme, p.127, p.259 (1st edition).