Neumann Laplace Direct Method#
keys: homogeneous Neumann bvp, hypersingular operator, unknown Dirichlet 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 Neumann boundary value problem
The unique solution \(u\in H^1(\Omega)\) can be written in the following form (representation forumula)
Let’s carefully apply the Neumann trace \(\gamma_1\) to the representation formula and solve the resulting boundary integral equation for the Dirichlet trace \(u_0\) of \(u\) by discretisation of the following variational formulation:
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(2)
Define conforming finite element spaces for \(H^{\frac12}(\Gamma)\) and \(H^{-\frac12}(\Gamma)\) respectively:
fesL2 = SurfaceL2(mesh, order=2, dual_mapping=False) # TODO: check interpolation with option dual_mapping=True
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 = 1576
Project the given Neumann data \(u_1\) in finite element space and have a look at the boundary data:
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)
u1 = GridFunction(fesL2)
u1.Interpolate(graduexa*n, definedon=mesh.Boundaries(".*"))
Draw (u1, mesh, draw_vol=False, order=3);
The discretisation of the variational formulation leads to a system of linear equations, ie
where the linear operators are as follows
\(\mathrm{D}\) is the hypersingular operator and the stabilisation \((\mathrm D + \mathrm{S})\) is regular and symmetric.
\(\mathrm{M}\) is a mass matrix.
\(\mathrm{K}'\) is the adjoint double layer operator.
The stabilisation matrix \(S\) is needed to cope with the kernel of the hypersingular operator. The following stabilisation matrix is used here:
where \(v_i, v_j\) being \(H^{\frac12}(\Gamma)\)-conforming shape functions.
vH1m1 = LinearForm(vH1*1*ds(bonus_intorder=3)).Assemble()
S = (BaseMatrix(Matrix(vH1m1.vec.Reshape(1))))@(BaseMatrix(Matrix(vH1m1.vec.Reshape(fesH1.ndof))))
Now we assemble the stabilised system matrix and the right hand side and solve for the Dirichlet data \(u_0\) with an iterative solver:
u0 = GridFunction(fesH1)
pre = BilinearForm(uH1*vH1*ds).Assemble().mat.Inverse(freedofs=fesH1.FreeDofs())
with TaskManager(): # pajetrace=1000*1000*1000):
D = HypersingularOperator(fesH1, intorder=12, leafsize=40, eta=3., eps=1e-4, method="aca", testhmatrix=False)
M = BilinearForm( v.Trace() * uH1 * ds(bonus_intorder=3)).Assemble()
Kt = DoubleLayerPotentialOperator(fesL2, fesH1, intorder=12, leafsize=40, eta=3., eps=1e-4, method="aca")
rhs = ( (0.5 * M.mat.T - Kt.mat) * u1.vec ).Evaluate()
CG(mat=D.mat+S, pre=pre, rhs=rhs, sol=u0.vec, tol=1e-8, maxsteps=200, initialize=False, printrates=True)
CG iteration 1, residual = 0.6734119976551737
CG iteration 2, residual = 0.35707083594364425
CG iteration 3, residual = 0.2088885133580026
CG iteration 4, residual = 0.17063011718383092
CG iteration 5, residual = 0.12922502344210976
CG iteration 6, residual = 0.10837761567746897
CG iteration 7, residual = 0.06658980509444792
CG iteration 8, residual = 0.05390070177633284
CG iteration 9, residual = 0.036168391258839784
CG iteration 10, residual = 0.02374733108860615
CG iteration 11, residual = 0.01713597564632652
CG iteration 12, residual = 0.010529013246750856
CG iteration 13, residual = 0.006843245189197806
CG iteration 14, residual = 0.004377148034702332
CG iteration 15, residual = 0.0027921834233506377
CG iteration 16, residual = 0.0017274928252799912
CG iteration 17, residual = 0.0010047327067879371
CG iteration 18, residual = 0.0005839403720918436
CG iteration 19, residual = 0.0003410189817395494
CG iteration 20, residual = 0.00019007063675795846
CG iteration 21, residual = 0.00010675896994184973
CG iteration 22, residual = 5.8059772565693345e-05
CG iteration 23, residual = 3.1498114935420734e-05
CG iteration 24, residual = 1.7107557048835916e-05
CG iteration 25, residual = 9.03739609887231e-06
CG iteration 26, residual = 4.626967507501204e-06
CG iteration 27, residual = 2.4760648865551604e-06
CG iteration 28, residual = 1.3318926519814203e-06
CG iteration 29, residual = 8.230311534842697e-07
CG iteration 30, residual = 6.74521031615184e-07
CG iteration 31, residual = 7.082765813517587e-07
CG iteration 32, residual = 7.260965308834853e-07
CG iteration 33, residual = 5.508925208098296e-07
CG iteration 34, residual = 2.9975673881253844e-07
CG iteration 35, residual = 1.4416911434241806e-07
CG iteration 36, residual = 6.545008867522095e-08
CG iteration 37, residual = 3.0354147814622845e-08
CG iteration 38, residual = 1.5207593600053382e-08
CG iteration 39, residual = 9.546096972303662e-09
CG iteration 40, residual = 8.58026905893728e-09
CG iteration 41, residual = 8.208573615557572e-09
CG iteration 42, residual = 5.6662151788201356e-09
Let’s have a look at the numerical and exact Dirichlet data and compare them quantitatively:
Draw (u0, mesh, draw_vol=False, order=3);
u0exa = GridFunction(fesH1)
u0exa.Interpolate (uexa, definedon=mesh.Boundaries(".*"))
Draw (u0exa, mesh, draw_vol=False, order=3);
print ("L2 error of surface gradients =", sqrt (Integrate ( (grad(u0exa)-grad(u0))**2, mesh, BND)))
L2 error of surface gradients = 0.004797701047762403