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.63000976575928
GMRes iteration 5, residual = 8.056881046688135
GMRes iteration 6, residual = 4.509288098284204
GMRes iteration 7, residual = 2.564403593233677
GMRes iteration 8, residual = 1.4899517632937291
GMRes iteration 9, residual = 0.9024815287499113
GMRes iteration 10, residual = 0.5758201669293832
GMRes iteration 11, residual = 0.3560656008770565
GMRes iteration 12, residual = 0.23167701957578665
GMRes iteration 13, residual = 0.1405573327246735
GMRes iteration 14, residual = 0.08643232123193934
GMRes iteration 15, residual = 0.04785581447218547
GMRes iteration 16, residual = 0.02550637222580812
GMRes iteration 17, residual = 0.01474662456004314
GMRes iteration 18, residual = 0.008510791402229797
GMRes iteration 19, residual = 0.0055103691242221545
GMRes iteration 20, residual = 0.003618707454400927
GMRes iteration 21, residual = 0.0024306137823109503
GMRes iteration 22, residual = 0.0015919306847563355
GMRes iteration 23, residual = 0.0010083373329430496
GMRes iteration 24, residual = 0.0006321666784882902
GMRes iteration 25, residual = 0.0004009038155148561
GMRes iteration 26, residual = 0.00026043698947360853
GMRes iteration 27, residual = 0.00017182267762203059
GMRes iteration 28, residual = 0.00010868824022773714
GMRes iteration 29, residual = 6.329089401456432e-05
GMRes iteration 30, residual = 3.863213272819033e-05
GMRes iteration 31, residual = 2.1392678367102154e-05
GMRes iteration 32, residual = 1.2290052532212022e-05
GMRes iteration 33, residual = 7.293316630848352e-06
GMRes iteration 34, residual = 4.796204871492728e-06
GMRes iteration 35, residual = 3.1141623903424515e-06
GMRes iteration 36, residual = 2.0084257702550194e-06
GMRes iteration 37, residual = 1.2865040236412238e-06
GMRes iteration 38, residual = 7.999753490173828e-07
GMRes iteration 39, residual = 5.195825979236753e-07
GMRes iteration 40, residual = 3.3627082886387124e-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

Alternative text

[ ]: