This page was generated from unit-11.3-bem-Helmholtz/BrakhageWerner.ipynb.
11.3.1 Helmholtz solver using Brakhage-Werner formulation¶
Combined field integral equations combine single and double layer integral operators, one simple option is the Brakhage-Werner formulation.
The solution is represented as
\[u = (i \kappa S - D) \phi,\]
where \(\phi\) solve the boundary integral equation
\[\big( \tfrac{1}{2} + K + i \kappa V \big) \phi = u_{in} \qquad \text{on} \, \Gamma\]
[1]:
from netgen.occ import *
from ngsolve import *
from ngsolve.webgui import Draw
from ngsolve.bem import *
[2]:
kappa=20
order=4
[3]:
screen = WorkPlane(Axes( (0,0,0), Z, X)).RectangleC(15,15).Face()
sp = Fuse(Sphere( (0,0,0), pi).faces)
screen.faces.name="screen"
sp.faces.name="sphere"
shape = Compound([screen,sp])
mesh = shape.GenerateMesh(maxh=5/kappa).Curve(order)
Draw (mesh);
[4]:
fes_sphere = Compress(SurfaceL2(mesh, order=order, complex=True, definedon=mesh.Boundaries("sphere")))
u,v = fes_sphere.TnT()
fes_screen = Compress(SurfaceL2(mesh, order=order, dual_mapping=True, complex=True, definedon=mesh.Boundaries("screen")))
print ("ndof_sphere = ", fes_sphere.ndof, "ndof_screen =", fes_screen.ndof)
ndof_sphere = 64530 ndof_screen = 125040
[ ]:
[5]:
with TaskManager():
C = HelmholtzCF(u*ds("sphere"), kappa)*v*ds
Id = BilinearForm(u*v*ds).Assemble()
[6]:
with TaskManager():
lhs = 0.5 * Id.mat + C.mat
source = exp(1j * kappa * x)
rhs = LinearForm(-source*v*ds).Assemble()
[7]:
gfu = GridFunction(fes_sphere)
pre = BilinearForm(u*v*ds, diagonal=True).Assemble().mat.Inverse()
with TaskManager():
gfu.vec[:] = solvers.GMRes(A=lhs, b=rhs.vec, pre=pre, maxsteps=40, tol=1e-8)
GMRes iteration 1, residual = 119.67056861548505
GMRes iteration 2, residual = 54.528905739702545
GMRes iteration 3, residual = 27.43583083323107
GMRes iteration 4, residual = 14.630009765759283
GMRes iteration 5, residual = 8.056881046688135
GMRes iteration 6, residual = 4.5092880982842045
GMRes iteration 7, residual = 2.5644035932336773
GMRes iteration 8, residual = 1.4899517632937291
GMRes iteration 9, residual = 0.902481528749911
GMRes iteration 10, residual = 0.5758201669293833
GMRes iteration 11, residual = 0.3560656008770565
GMRes iteration 12, residual = 0.23167701957578646
GMRes iteration 13, residual = 0.14055733272467313
GMRes iteration 14, residual = 0.08643232123193911
GMRes iteration 15, residual = 0.047855814472185336
GMRes iteration 16, residual = 0.025506372225821947
GMRes iteration 17, residual = 0.014746624560045246
GMRes iteration 18, residual = 0.008510791402233344
GMRes iteration 19, residual = 0.005510369124219757
GMRes iteration 20, residual = 0.00361870745439895
GMRes iteration 21, residual = 0.0024306137823114768
GMRes iteration 22, residual = 0.0015919306847560811
GMRes iteration 23, residual = 0.001008337332940999
GMRes iteration 24, residual = 0.0006321666784875867
GMRes iteration 25, residual = 0.0004009038155143265
GMRes iteration 26, residual = 0.00026043698947311495
GMRes iteration 27, residual = 0.00017182267762168868
GMRes iteration 28, residual = 0.00010868824022759562
GMRes iteration 29, residual = 6.329089401439324e-05
GMRes iteration 30, residual = 3.863213272799688e-05
GMRes iteration 31, residual = 2.1392678367079166e-05
GMRes iteration 32, residual = 1.229005253222999e-05
GMRes iteration 33, residual = 7.293316630860099e-06
GMRes iteration 34, residual = 4.796204871532231e-06
GMRes iteration 35, residual = 3.1141623903534824e-06
GMRes iteration 36, residual = 2.008425763402976e-06
GMRes iteration 37, residual = 1.2865040185392884e-06
GMRes iteration 38, residual = 7.999753400619245e-07
GMRes iteration 39, residual = 5.195825919563383e-07
GMRes iteration 40, residual = 3.362708290271326e-07
WARNING: GMRes did not converge to TOL
[8]:
Draw (gfu, order=5, min=-1, max=1);
prostprocessing on screen¶
[9]:
uscat = GridFunction(fes_screen)
with TaskManager(pajetrace=10**8):
intop = 1j*kappa*HelmholtzSL(u*ds("sphere"),kappa) + (-1)*HelmholtzDL(u*ds("sphere"), kappa)
uscat.Set(intop(gfu, mesh.Boundaries("screen")), definedon=mesh.Boundaries("screen"))
[10]:
print ("Scattered field")
Draw (uscat, mesh, min=-1,max=1, animate_complex=True, order=4);
Scattered field
[11]:
uin = mesh.BoundaryCF( {"screen": source }, default=0)
print ("Total field")
Draw (uin-uscat, mesh, min=-1,max=1, animate_complex=True, order=4);
Total field
Scattering from sphere with \(D = 50 \lambda\). About 5 min on Macbook Apple M4 Pro

[ ]: