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=10
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 = 17010 ndof_screen = 31110
[ ]:
[5]:
with TaskManager():
# V = HelmholtzSingleLayerPotentialOperator(fes_sphere, fes_sphere, kappa=kappa, intorder=10)
# K = HelmholtzDoubleLayerPotentialOperator(fes_sphere, fes_sphere, kappa=kappa, intorder=10)
# C = HelmholtzCombinedFieldOperator(fes_sphere, fes_sphere, kappa=kappa, intorder=10)
C = HelmholtzCF(u*ds("sphere"), kappa)*v*ds
u,v = fes_sphere.TnT()
Id = BilinearForm(u*v*ds).Assemble()
[6]:
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 = 60.52440579389443
GMRes iteration 2, residual = 24.05582183523761
GMRes iteration 3, residual = 10.866381612939417
GMRes iteration 4, residual = 5.207724837026282
GMRes iteration 5, residual = 2.5606451105593298
GMRes iteration 6, residual = 1.2928262518532614
GMRes iteration 7, residual = 0.6567675233068868
GMRes iteration 8, residual = 0.35921392193432033
GMRes iteration 9, residual = 0.20001557577759402
GMRes iteration 10, residual = 0.10890790559448302
GMRes iteration 11, residual = 0.05889051992112618
GMRes iteration 12, residual = 0.03274406127301688
GMRes iteration 13, residual = 0.013561377443551825
GMRes iteration 14, residual = 0.007286640763756341
GMRes iteration 15, residual = 0.0036319265850675553
GMRes iteration 16, residual = 0.002036492293296241
GMRes iteration 17, residual = 0.0012139471515801256
GMRes iteration 18, residual = 0.0007228533895894179
GMRes iteration 19, residual = 0.0004301332817236131
GMRes iteration 20, residual = 0.00025834108327070607
GMRes iteration 21, residual = 0.00014417783739587156
GMRes iteration 22, residual = 7.770652618062814e-05
GMRes iteration 23, residual = 4.468014818822641e-05
GMRes iteration 24, residual = 2.211722906069835e-05
GMRes iteration 25, residual = 1.1247609798545048e-05
GMRes iteration 26, residual = 5.988321048133865e-06
GMRes iteration 27, residual = 3.2611582435840496e-06
GMRes iteration 28, residual = 1.7721145247162111e-06
GMRes iteration 29, residual = 1.0683563619743223e-06
GMRes iteration 30, residual = 5.955571812726123e-07
GMRes iteration 31, residual = 3.451428639678431e-07
GMRes iteration 32, residual = 1.9976493608372992e-07
GMRes iteration 33, residual = 1.119930207900579e-07
GMRes iteration 34, residual = 6.802286055499613e-08
GMRes iteration 35, residual = 4.01633537521196e-08
GMRes iteration 36, residual = 1.8963293310897004e-08
GMRes iteration 37, residual = 1.0014341697476182e-08
GMRes iteration 38, residual = 5.596424976143244e-09
[8]:
Draw (gfu, order=5, min=-1, max=1);
prostprocessing on screen¶
[9]:
uscat = GridFunction(fes_screen)
with TaskManager(pajetrace=10**8):
# uscat.Set(1j*kappa*V.GetPotential(gfu)-K.GetPotential(gfu) , definedon=mesh.Boundaries("screen"))
# uscat.Set(HelmholtzCF(u*ds("sphere"), kappa)(gfu) , definedon=mesh.Boundaries("screen"))
uscat.Set(1j*kappa*HelmholtzSL(u*ds("sphere"),kappa)(gfu) - HelmholtzDL(u*ds("sphere"), kappa)(gfu), \
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

[ ]: