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
Let us choose the following ansatz for the solution \(u\in H^1(\Omega)\) (indirect ansatz)
and solve for the density \(j\in H^{-\frac12}(\Gamma)\) by the boundary element method, i.e. the numerical solution of the variational formulation
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
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)
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.5510990298602576
CG iteration 3, residual = 0.1311274118731449
CG iteration 4, residual = 0.031078921260816843
CG iteration 5, residual = 0.007372481975125887
CG iteration 6, residual = 0.0017929371359744126
CG iteration 7, residual = 0.0005259800974229256
CG iteration 8, residual = 0.0002580777788551423
CG iteration 9, residual = 0.0002832379605424888
CG iteration 10, residual = 0.00012036673026200608
CG iteration 11, residual = 0.0002350713728555535
CG iteration 12, residual = 4.492368457430169e-05
CG iteration 13, residual = 4.801474714715425e-05
CG iteration 14, residual = 2.782522946264007e-05
CG iteration 15, residual = 1.4627805195051116e-05
CG iteration 16, residual = 1.4781005110121252e-05
CG iteration 17, residual = 9.771555083328643e-06
CG iteration 18, residual = 6.175962509923531e-06
CG iteration 19, residual = 1.5387353065062504e-05
CG iteration 20, residual = 3.2801496585677222e-06
CG iteration 21, residual = 2.970067358547831e-06
CG iteration 22, residual = 3.493169500416989e-06
CG iteration 23, residual = 2.5937980739747343e-06
CG iteration 24, residual = 1.2604332742721193e-06
CG iteration 25, residual = 8.37292446318865e-07
CG iteration 26, residual = 1.2440302529203402e-06
CG iteration 27, residual = 5.872586568701751e-07
CG iteration 28, residual = 4.4146187934116504e-07
CG iteration 29, residual = 5.781842795522604e-07
CG iteration 30, residual = 3.5267750820843127e-07
CG iteration 31, residual = 1.778225528323564e-07
CG iteration 32, residual = 1.5441963858753221e-07
CG iteration 33, residual = 2.1911311136998104e-07
CG iteration 34, residual = 1.2550251284310553e-07
CG iteration 35, residual = 8.723802307658936e-08
CG iteration 36, residual = 6.888879164602111e-08
CG iteration 37, residual = 1.0272143720383242e-07
CG iteration 38, residual = 3.437270110360088e-08
CG iteration 39, residual = 4.024009220704828e-08
CG iteration 40, residual = 5.7328691039674025e-08
CG iteration 41, residual = 1.8547249842788816e-08
[6]:
Draw (j, order=3);
Notes:
For details on the analysis of boundary integral equations derived from elliptic partial differential equations, see for instance Strongly Elliptic Systems and Boundary Integral Equations.
The integration of singular pairings is done as proposed in Randelementmethoden
The adaptive cross approximation is done as proposed in Hierarchical Matrices.
[ ]: