from netgen.occ import *
import netgen.meshing as meshing
from ngsolve import *
from ngsolve.webgui import Draw
from ngsolve.bem import *
from ngsolve import Projector, Preconditioner
from ngsolve.krylovspace import CG, GMRes
keys: PEC Scattering with mixed conditions, Calderon projector, Dirichlet data, Neumann data
Maxwell Mixed Conditions#
Consider the mixed problem:
\( \begin{array}{l rcl r} & \mathrm{curl}\,\mathrm{curl} \boldsymbol E - \kappa^2 \boldsymbol E &=& \boldsymbol 0 &\mathrm{in}\; \Omega^c \,,\\ \textnormal{Dirichlet condition} & \gamma_R \boldsymbol E &=& \boldsymbol m & \mathrm{on}\; \Gamma_0\,,\\ \textnormal{Neumann condition} & \gamma_N \boldsymbol E &=& \boldsymbol j & \mathrm{on}\; \Gamma_1\,, \\ \textnormal{radiation condition} & | \mathrm{curl}\boldsymbol E \times \dfrac{x}{|x|} - i \omega\epsilon \boldsymbol E(x)| &=& \mathcal O(1/|x|^2), & |x|\to \infty\,. \end{array} \) |
\(\quad\quad\quad\) |
lower half sphere: Dirichlet boundary with given \(\gamma_D \, E = \boldsymbol m\)
upper half sphere: Neumann boundary with given \(\gamma_N \, E = \boldsymbol j\)
NGSolve solution
Define geometry and mesh:
topsphere = Sphere((0,0,0), 1) * Box((-1,-1,0),(1,1,1))
botsphere = Sphere((0,0,0), 1) - Box((-1,-1,0),(1,1,1))
topsphere.faces.name = "neumann" # 1/kappa n x curl E
botsphere.faces.name = "dirichlet" # nxE
shape = Fuse( [topsphere,botsphere] )
h = 0.5
order = 3
mesh = Mesh( OCCGeometry(shape).GenerateMesh(maxh=h, perfstepsend=meshing.MeshingStep.MESHSURFACE)).Curve(order)
#Draw (mesh);
Define a manufactured solution \(\boldsymbol E\):
# manufactured solution: kernel with source in point s and direction e
kappa = 1. / h
s = CF((0., 0., 0.))
e = (1.192, -0.189, 2.745)
e = CF(e) / sqrt(e[0] ** 2 + e[1] ** 2 + e[2] ** 2)
xms = CF((x, y , z)) - s
r = Norm(xms)
f = (-kappa ** 2 * r ** 2 - 1j * kappa * r + 1.) / r ** 3
g = (e * xms) * (-kappa ** 2 * r ** 2 - 3. * 1j * kappa * r + 3.) / r ** 5
E = exp(1j * kappa * r) * (f * e - g * xms)
curlE = CF((E[2].Diff(y) - E[1].Diff(z), E[0].Diff(z) - E[2].Diff(x), E[1].Diff(x) - E[0].Diff(y)))
print("kappa = ", kappa)
kappa = 2.0
Boundary Integral Equations
The Calderon projector relates the Dirichlet and the Neumann traces of the solution \(u\), i.e.,
and we use it to solve for the Dirchlet trace on \(\Gamma_1\) and the Neumann trace on \(\Gamma_0\).
Define the finite element product space with \(H^{-\frac12}(\mathrm{curl}_\Gamma, \Gamma)\) and \(H^{-\frac12}(\mathrm{div}_\Gamma, \Gamma)\) conforming finte elements:
order = 3
fesHDiv = HDivSurface(mesh, order=order, complex=True, dirichlet="neumann")
fesHCurl = HCurl(mesh, order=order, complex=True, dirichlet="dirichlet")
fes = fesHCurl * fesHDiv
(uHCurl,uHDiv),(vHCurl,vHDiv) = fes.TnT()
Get the Dirichlet trace \(\boldsymbol m\) and exact Neumann data \(\boldsymbol j\) as grid functions defined everywhere on \(\Gamma\):
n = specialcf.normal(3)
mjexa = GridFunction(fes)
mj = GridFunction(fes)
mjexa.components[0].Set( Cross( Cross(n, E), n) , definedon=mesh.Boundaries(".*"), dual=True) # Hcurl, gamma_R = (nxE)xn
mjexa.components[1].Set( 1/kappa*Cross(n,curlE) , definedon=mesh.Boundaries(".*"), dual=True) # Hdiv, gamma_N = 1/k nx curl E
Define Dirichlet boundary condition on the lower half sphere \(\Gamma_0\) and have a look at it:
mj.components[0].Set( mjexa.components[0], definedon=mesh.Boundaries("dirichlet"), dual=True) # given Dirichlet data
Draw (Norm(mj.components[0]), mesh, draw_vol=False, order=3);
Define Neumann boundary condition on the upper half sphere \(\Gamma_1\) and have a look at it:
mj.components[1].Set( mjexa.components[1], definedon=mesh.Boundaries("neumann"), dual=True) # given Neumann data
Draw (Norm(mj.components[1]), mesh, draw_vol=False, order=3);
Compute boundary integral operators \(\mathrm{V}\), \(\mathrm{K}\), \(\mathrm{D}\) and the mass matrix \(\mathrm M\)
# M, V, K and D
intorder = 2 * order + 6
with TaskManager():
M = BilinearForm( uHCurl.Trace() * vHDiv.Trace()* ds(bonus_intorder=3)).Assemble() # <Hcurl, Hdiv>
V1 = HelmholtzSL(uHDiv.Trace()*ds(bonus_intorder=6), kappa)*vHDiv.Trace()*ds(bonus_intorder=6)
V2 = HelmholtzSL(div(uHDiv.Trace())*ds(bonus_intorder=6), kappa)*div(vHDiv.Trace())*ds(bonus_intorder=6)
V = kappa * V1.mat - 1/kappa * V2.mat
K = HelmholtzDL(uHCurl.Operator("rotated_trace") * ds(bonus_intorder=6), kappa) * vHDiv.Trace() * ds(bonus_intorder=6)
D1 = HelmholtzSL(uHCurl.Operator("rotated_trace")*ds(bonus_intorder=6), kappa)*uHCurl.Operator("rotated_trace")*ds(bonus_intorder=6)
D2 = HelmholtzSL(uHCurl.Operator("surface_curl")*ds(bonus_intorder=6), kappa)*vHCurl.Operator("surface_curl")*ds(bonus_intorder=6)
D = kappa * D1.mat - 1/kappa * D2.mat
Generate right hand side, block matrix and solve the linear system of equations:
and we use it to solve for the Dirchlet trace on \(\Gamma_1\) and the Neumann trace on \(\Gamma_0\).
with TaskManager():
lhs = V + K.mat - K.mat.T - D
rhs = ((0.5 * M.mat - K.mat - V + 0.5 * M.mat.T + K.mat.T + D ) * mj.vec).Evaluate()
fesHDivN = HDivSurface(mesh, order=order, dirichlet="neumann") # Dirichlet nodes free dofs
fesHCurlD = HCurl(mesh, order=order, dirichlet="dirichlet") # Neumann nodes free dofs
bfpre = BilinearForm(uHDiv.Trace() * vHDiv.Trace() *ds + uHCurl.Trace() * vHCurl.Trace() *ds).Assemble()
pre = bfpre.mat.Inverse(freedofs=fes.FreeDofs())
with TaskManager():
sol = GMRes(A=lhs, b=rhs, pre=pre, maxsteps=3000, tol=1e-8)
GMRes iteration 1, residual = 59.03364496742984
GMRes iteration 2, residual = 34.60392533234917
GMRes iteration 3, residual = 31.396373978279804
GMRes iteration 4, residual = 20.366185558777005
GMRes iteration 5, residual = 20.15196935823866
GMRes iteration 6, residual = 14.598615604897471
GMRes iteration 7, residual = 14.577009704944333
GMRes iteration 8, residual = 11.496582346023116
GMRes iteration 9, residual = 11.37726266330068
GMRes iteration 10, residual = 10.090925774776247
GMRes iteration 11, residual = 9.833380337813571
GMRes iteration 12, residual = 9.538129904000545
GMRes iteration 13, residual = 9.135079961273414
GMRes iteration 14, residual = 9.012107836315268
GMRes iteration 15, residual = 8.472899899277525
GMRes iteration 16, residual = 8.457119728187186
GMRes iteration 17, residual = 8.055121423334915
GMRes iteration 18, residual = 8.023055567242485
GMRes iteration 19, residual = 7.768891547752717
GMRes iteration 20, residual = 7.655081621462767
GMRes iteration 21, residual = 7.542847300178079
GMRes iteration 22, residual = 7.449033166884139
GMRes iteration 23, residual = 7.271849030471795
GMRes iteration 24, residual = 7.182198409820015
GMRes iteration 25, residual = 6.821672562948277
GMRes iteration 26, residual = 6.744206212586918
GMRes iteration 27, residual = 6.420722968570941
GMRes iteration 28, residual = 6.330316839665963
GMRes iteration 29, residual = 6.001172307659864
GMRes iteration 30, residual = 5.924547862943652
GMRes iteration 31, residual = 5.541323092926721
GMRes iteration 32, residual = 5.447853762490792
GMRes iteration 33, residual = 5.1953795187148515
GMRes iteration 34, residual = 5.081793282654341
GMRes iteration 35, residual = 4.8984444268256055
GMRes iteration 36, residual = 4.6880249007915715
GMRes iteration 37, residual = 4.611425672979815
GMRes iteration 38, residual = 4.400542091799011
GMRes iteration 39, residual = 4.357349625973269
GMRes iteration 40, residual = 4.167943441124538
GMRes iteration 41, residual = 4.146382995697436
GMRes iteration 42, residual = 4.005138296189357
GMRes iteration 43, residual = 4.000918750009959
GMRes iteration 44, residual = 3.8284273716507364
GMRes iteration 45, residual = 3.825117087791785
GMRes iteration 46, residual = 3.67386394428674
GMRes iteration 47, residual = 3.666391287735034
GMRes iteration 48, residual = 3.555348254671442
GMRes iteration 49, residual = 3.5543781477852203
GMRes iteration 50, residual = 3.463350591093729
GMRes iteration 51, residual = 3.4629128735403905
GMRes iteration 52, residual = 3.3552202171641983
GMRes iteration 53, residual = 3.351088120800925
GMRes iteration 54, residual = 3.2655246955594173
GMRes iteration 55, residual = 3.2619147344257615
GMRes iteration 56, residual = 3.177389479339349
GMRes iteration 57, residual = 3.175257552473688
GMRes iteration 58, residual = 3.102332616174091
GMRes iteration 59, residual = 3.0936783828942205
GMRes iteration 60, residual = 3.007651161376715
GMRes iteration 61, residual = 2.9896254277038072
GMRes iteration 62, residual = 2.888804547019958
GMRes iteration 63, residual = 2.8493498169161655
GMRes iteration 64, residual = 2.7789663128847746
GMRes iteration 65, residual = 2.7352634343777824
GMRes iteration 66, residual = 2.7002545993833054
GMRes iteration 67, residual = 2.647593007614043
GMRes iteration 68, residual = 2.6025212938105886
GMRes iteration 69, residual = 2.5256565924278447
GMRes iteration 70, residual = 2.4986723389183494
GMRes iteration 71, residual = 2.416149088454926
GMRes iteration 72, residual = 2.4027864530456253
GMRes iteration 73, residual = 2.3282157411534015
GMRes iteration 74, residual = 2.312412475616787
GMRes iteration 75, residual = 2.244560900228166
GMRes iteration 76, residual = 2.212693580603435
GMRes iteration 77, residual = 2.162173954583258
GMRes iteration 78, residual = 2.142213867024194
GMRes iteration 79, residual = 2.111187139237609
GMRes iteration 80, residual = 2.1068688265698845
GMRes iteration 81, residual = 2.0782879287433516
GMRes iteration 82, residual = 2.075551853897131
GMRes iteration 83, residual = 2.0463473727269266
GMRes iteration 84, residual = 2.0460170927633903
GMRes iteration 85, residual = 2.0140232376548775
GMRes iteration 86, residual = 2.013558067187594
GMRes iteration 87, residual = 1.977217589435647
GMRes iteration 88, residual = 1.9745033700396262
GMRes iteration 89, residual = 1.9451098308799375
GMRes iteration 90, residual = 1.9371648572303077
GMRes iteration 91, residual = 1.9100209762897673
GMRes iteration 92, residual = 1.9017026522243947
GMRes iteration 93, residual = 1.8580183555392504
GMRes iteration 94, residual = 1.8502157787647755
GMRes iteration 95, residual = 1.8063241188949324
GMRes iteration 96, residual = 1.80292542705256
GMRes iteration 97, residual = 1.7471416623909128
GMRes iteration 98, residual = 1.745393577481051
GMRes iteration 99, residual = 1.7166751436859435
GMRes iteration 100, residual = 1.7166212064394875
GMRes iteration 101, residual = 1.6968121552547906
GMRes iteration 102, residual = 1.6962293731887668
GMRes iteration 103, residual = 1.6757166887018569
GMRes iteration 104, residual = 1.6738416621912557
GMRes iteration 105, residual = 1.6583127544319856
GMRes iteration 106, residual = 1.654197699420285
GMRes iteration 107, residual = 1.6372324521335537
GMRes iteration 108, residual = 1.6266120247763278
GMRes iteration 109, residual = 1.6105437449699622
GMRes iteration 110, residual = 1.5942698772641668
GMRes iteration 111, residual = 1.5785311741633379
GMRes iteration 112, residual = 1.56405199496855
GMRes iteration 113, residual = 1.5479891772935863
GMRes iteration 114, residual = 1.5390207870543937
GMRes iteration 115, residual = 1.5144816706840585
GMRes iteration 116, residual = 1.5107407172480738
GMRes iteration 117, residual = 1.4793156461456973
GMRes iteration 118, residual = 1.472800420857648
GMRes iteration 119, residual = 1.4388055126672712
GMRes iteration 120, residual = 1.4283616642701846
GMRes iteration 121, residual = 1.4111066797201333
GMRes iteration 122, residual = 1.40090991643511
GMRes iteration 123, residual = 1.3885989001072203
GMRes iteration 124, residual = 1.3749113539657292
GMRes iteration 125, residual = 1.3672312561078392
GMRes iteration 126, residual = 1.3471501089593896
GMRes iteration 127, residual = 1.3433574624577849
GMRes iteration 128, residual = 1.3296156115169608
GMRes iteration 129, residual = 1.3254103279716345
GMRes iteration 130, residual = 1.310521233851445
GMRes iteration 131, residual = 1.309917563501425
GMRes iteration 132, residual = 1.2926011522703456
GMRes iteration 133, residual = 1.2916524792972441
GMRes iteration 134, residual = 1.2747417938787813
GMRes iteration 135, residual = 1.2746601470141716
GMRes iteration 136, residual = 1.2571581691230698
GMRes iteration 137, residual = 1.2564535835799786
GMRes iteration 138, residual = 1.2428287376996159
GMRes iteration 139, residual = 1.2405203792210433
GMRes iteration 140, residual = 1.227046334488343
GMRes iteration 141, residual = 1.2241123500131088
GMRes iteration 142, residual = 1.2115552109916792
GMRes iteration 143, residual = 1.2093786012679215
GMRes iteration 144, residual = 1.1920746515650316
GMRes iteration 145, residual = 1.1910836679204386
GMRes iteration 146, residual = 1.172206698726013
GMRes iteration 147, residual = 1.1713008119549508
GMRes iteration 148, residual = 1.1565961063712848
GMRes iteration 149, residual = 1.155978184025668
GMRes iteration 150, residual = 1.1427740980852068
GMRes iteration 151, residual = 1.1414294175454702
GMRes iteration 152, residual = 1.1323326692705662
GMRes iteration 153, residual = 1.1304356067966137
GMRes iteration 154, residual = 1.1196310677448804
GMRes iteration 155, residual = 1.1179396121271008
GMRes iteration 156, residual = 1.1025920761189656
GMRes iteration 157, residual = 1.102486264846769
GMRes iteration 158, residual = 1.090870131160797
GMRes iteration 159, residual = 1.0900883092359195
GMRes iteration 160, residual = 1.0749844526288117
GMRes iteration 161, residual = 1.071014176220022
GMRes iteration 162, residual = 1.0475194415827893
GMRes iteration 163, residual = 1.037897517444736
GMRes iteration 164, residual = 1.017567628827744
GMRes iteration 165, residual = 1.0068907566691092
GMRes iteration 166, residual = 0.9950979262129316
GMRes iteration 167, residual = 0.9924199315057521
GMRes iteration 168, residual = 0.983489014379227
GMRes iteration 169, residual = 0.9819919406892204
GMRes iteration 170, residual = 0.9745463118876063
GMRes iteration 171, residual = 0.9719512156111725
GMRes iteration 172, residual = 0.9629058876664627
GMRes iteration 173, residual = 0.9600541966583012
GMRes iteration 174, residual = 0.9541730814953516
GMRes iteration 175, residual = 0.9514298665227457
GMRes iteration 176, residual = 0.9480023225456244
GMRes iteration 177, residual = 0.9451800839260921
GMRes iteration 178, residual = 0.9415202294642727
GMRes iteration 179, residual = 0.937673381821984
GMRes iteration 180, residual = 0.9327875587878602
GMRes iteration 181, residual = 0.9294460303538905
GMRes iteration 182, residual = 0.9274288615398905
GMRes iteration 183, residual = 0.9226917580139251
GMRes iteration 184, residual = 0.9173862413564655
GMRes iteration 185, residual = 0.9110234327439176
GMRes iteration 186, residual = 0.908789704381793
GMRes iteration 187, residual = 0.9041885304261958
GMRes iteration 188, residual = 0.8978483075851992
GMRes iteration 189, residual = 0.8950173258493018
GMRes iteration 190, residual = 0.8869815766777017
GMRes iteration 191, residual = 0.8833127045048346
GMRes iteration 192, residual = 0.8756273716152919
GMRes iteration 193, residual = 0.870497316141558
GMRes iteration 194, residual = 0.8658542821078472
GMRes iteration 195, residual = 0.8633025150314724
GMRes iteration 196, residual = 0.8627513312337832
GMRes iteration 197, residual = 0.8607973029667135
GMRes iteration 198, residual = 0.8607010759840342
GMRes iteration 199, residual = 0.858429553104729
GMRes iteration 200, residual = 0.857910753973286
GMRes iteration 201, residual = 0.8544386521273154
GMRes iteration 202, residual = 0.8531988722887688
GMRes iteration 203, residual = 0.848724432657956
GMRes iteration 204, residual = 0.8473885886310242
GMRes iteration 205, residual = 0.8439492025017704
GMRes iteration 206, residual = 0.841879206851401
GMRes iteration 207, residual = 0.8361292957787072
GMRes iteration 208, residual = 0.8351506287737985
GMRes iteration 209, residual = 0.8318466531552048
GMRes iteration 210, residual = 0.8317943492856363
GMRes iteration 211, residual = 0.8293881606366217
GMRes iteration 212, residual = 0.8292702243769355
GMRes iteration 213, residual = 0.8269652855066191
GMRes iteration 214, residual = 0.8268830221803236
GMRes iteration 215, residual = 0.822570166423503
GMRes iteration 216, residual = 0.8225695784176121
GMRes iteration 217, residual = 0.8171508657677836
GMRes iteration 218, residual = 0.8170469591560663
GMRes iteration 219, residual = 0.8108011672254435
GMRes iteration 220, residual = 0.8105343764900048
GMRes iteration 221, residual = 0.8041530589464218
GMRes iteration 222, residual = 0.8036943935539725
GMRes iteration 223, residual = 0.8008324207139653
GMRes iteration 224, residual = 0.8001776907092212
GMRes iteration 225, residual = 0.795936908396567
GMRes iteration 226, residual = 0.7958790703732961
GMRes iteration 227, residual = 0.7896114330251371
GMRes iteration 228, residual = 0.7887947646561662
GMRes iteration 229, residual = 0.7839336765807687
GMRes iteration 230, residual = 0.7812100744412964
GMRes iteration 231, residual = 0.7722432612857303
GMRes iteration 232, residual = 0.770167045581368
GMRes iteration 233, residual = 0.7622236129357802
GMRes iteration 234, residual = 0.7592067390582712
GMRes iteration 235, residual = 0.7554917543647065
GMRes iteration 236, residual = 0.7539670252512684
GMRes iteration 237, residual = 0.7524791584456303
GMRes iteration 238, residual = 0.7518008664462434
GMRes iteration 239, residual = 0.7493333872480269
GMRes iteration 240, residual = 0.7486403580985223
GMRes iteration 241, residual = 0.7424871985544138
GMRes iteration 242, residual = 0.7398852780999804
GMRes iteration 243, residual = 0.733620595752443
GMRes iteration 244, residual = 0.7280784854313642
GMRes iteration 245, residual = 0.7225040933333997
GMRes iteration 246, residual = 0.7141224465527257
GMRes iteration 247, residual = 0.7104748679626018
GMRes iteration 248, residual = 0.7026732614948676
GMRes iteration 249, residual = 0.7005148707068658
GMRes iteration 250, residual = 0.6924330586518151
GMRes iteration 251, residual = 0.691586619321591
GMRes iteration 252, residual = 0.6839042548112347
GMRes iteration 253, residual = 0.6822954628739302
GMRes iteration 254, residual = 0.6693712951188008
GMRes iteration 255, residual = 0.6690499902089141
GMRes iteration 256, residual = 0.6574358630497993
GMRes iteration 257, residual = 0.6574090513643146
GMRes iteration 258, residual = 0.6476650498877927
GMRes iteration 259, residual = 0.6466563196654331
GMRes iteration 260, residual = 0.6355996448327448
GMRes iteration 261, residual = 0.6351319815170348
GMRes iteration 262, residual = 0.6290790186652806
GMRes iteration 263, residual = 0.6282267084065131
GMRes iteration 264, residual = 0.6239377538489314
GMRes iteration 265, residual = 0.6235810661514085
GMRes iteration 266, residual = 0.6215335500967081
GMRes iteration 267, residual = 0.6213131248480004
GMRes iteration 268, residual = 0.6200812512081416
GMRes iteration 269, residual = 0.6196997072991296
GMRes iteration 270, residual = 0.6183119258207911
GMRes iteration 271, residual = 0.6183070053119729
GMRes iteration 272, residual = 0.6154632838174017
GMRes iteration 273, residual = 0.6154464397953593
GMRes iteration 274, residual = 0.6091912647910324
GMRes iteration 275, residual = 0.6088243955856872
GMRes iteration 276, residual = 0.596633032958513
GMRes iteration 277, residual = 0.5962218388572537
GMRes iteration 278, residual = 0.577649186587718
GMRes iteration 279, residual = 0.5771887088034808
GMRes iteration 280, residual = 0.5572877555813233
GMRes iteration 281, residual = 0.5567457276733612
GMRes iteration 282, residual = 0.5412298140699398
GMRes iteration 283, residual = 0.5387081146209673
GMRes iteration 284, residual = 0.5335593460950119
GMRes iteration 285, residual = 0.5311629531317149
GMRes iteration 286, residual = 0.5283905778155445
GMRes iteration 287, residual = 0.5234933111645169
GMRes iteration 288, residual = 0.5226124579893024
GMRes iteration 289, residual = 0.514801626697861
GMRes iteration 290, residual = 0.5120740644035021
GMRes iteration 291, residual = 0.5084381084235841
GMRes iteration 292, residual = 0.5044814999633506
GMRes iteration 293, residual = 0.5019635257123055
GMRes iteration 294, residual = 0.49789444385770915
GMRes iteration 295, residual = 0.49638788807960404
GMRes iteration 296, residual = 0.4935146363059558
GMRes iteration 297, residual = 0.4928732844869763
GMRes iteration 298, residual = 0.4876477710238178
GMRes iteration 299, residual = 0.4875483466727772
GMRes iteration 300, residual = 0.4824318186889675
GMRes iteration 301, residual = 0.4821560523236557
GMRes iteration 302, residual = 0.47694481539237615
GMRes iteration 303, residual = 0.4768373965429496
GMRes iteration 304, residual = 0.4710420530459894
GMRes iteration 305, residual = 0.47017356908505015
GMRes iteration 306, residual = 0.4660957126862427
GMRes iteration 307, residual = 0.46568744049338257
GMRes iteration 308, residual = 0.46108982513450264
GMRes iteration 309, residual = 0.45853386675068164
GMRes iteration 310, residual = 0.4543995950288249
GMRes iteration 311, residual = 0.4512875721758125
GMRes iteration 312, residual = 0.44889556767759015
GMRes iteration 313, residual = 0.4468955958729493
GMRes iteration 314, residual = 0.4434826343959677
GMRes iteration 315, residual = 0.43829506676778257
GMRes iteration 316, residual = 0.4371849647513269
GMRes iteration 317, residual = 0.4311434317067343
GMRes iteration 318, residual = 0.4288668368814687
GMRes iteration 319, residual = 0.4242721716436136
GMRes iteration 320, residual = 0.42263471098522126
GMRes iteration 321, residual = 0.4201330509762861
GMRes iteration 322, residual = 0.41861313921675725
GMRes iteration 323, residual = 0.41656766183624616
GMRes iteration 324, residual = 0.415239486607628
GMRes iteration 325, residual = 0.4120537090252377
GMRes iteration 326, residual = 0.41086332101924933
GMRes iteration 327, residual = 0.40873896551714123
GMRes iteration 328, residual = 0.408155937128319
GMRes iteration 329, residual = 0.4058954471362874
GMRes iteration 330, residual = 0.4056624771207148
GMRes iteration 331, residual = 0.4035914180673997
GMRes iteration 332, residual = 0.4034189005626761
GMRes iteration 333, residual = 0.4022874157676227
GMRes iteration 334, residual = 0.40182732044360225
GMRes iteration 335, residual = 0.4009943882363879
GMRes iteration 336, residual = 0.4006235637018285
GMRes iteration 337, residual = 0.3990932186027976
GMRes iteration 338, residual = 0.3987671017948297
GMRes iteration 339, residual = 0.395147617149015
GMRes iteration 340, residual = 0.39479759745893356
GMRes iteration 341, residual = 0.3927361124892041
GMRes iteration 342, residual = 0.3923995328509086
GMRes iteration 343, residual = 0.3904254030907426
GMRes iteration 344, residual = 0.39034743640067554
GMRes iteration 345, residual = 0.3884846920603609
GMRes iteration 346, residual = 0.3882968839178699
GMRes iteration 347, residual = 0.38657148137153335
GMRes iteration 348, residual = 0.3864996763186194
GMRes iteration 349, residual = 0.38316201607554434
GMRes iteration 350, residual = 0.3830812781512058
GMRes iteration 351, residual = 0.37901650588518193
GMRes iteration 352, residual = 0.3787899276697622
GMRes iteration 353, residual = 0.3767808756504366
GMRes iteration 354, residual = 0.3761297830174
GMRes iteration 355, residual = 0.37522246606311255
GMRes iteration 356, residual = 0.37494241586643806
GMRes iteration 357, residual = 0.37381088923530165
GMRes iteration 358, residual = 0.37263519000601464
GMRes iteration 359, residual = 0.36941644398485973
GMRes iteration 360, residual = 0.36633385703857674
GMRes iteration 361, residual = 0.36300187737602385
GMRes iteration 362, residual = 0.3597888738095794
GMRes iteration 363, residual = 0.35718251097899695
GMRes iteration 364, residual = 0.3538861870982611
GMRes iteration 365, residual = 0.35130924903971283
GMRes iteration 366, residual = 0.35017062219361544
GMRes iteration 367, residual = 0.34835687914180175
GMRes iteration 368, residual = 0.3473462161217514
GMRes iteration 369, residual = 0.34633789473019283
GMRes iteration 370, residual = 0.34562104188023157
GMRes iteration 371, residual = 0.34439392408024616
GMRes iteration 372, residual = 0.34421637524664067
GMRes iteration 373, residual = 0.3431015099247594
GMRes iteration 374, residual = 0.34274776615895747
GMRes iteration 375, residual = 0.3414427858510462
GMRes iteration 376, residual = 0.33807773211038405
GMRes iteration 377, residual = 0.33722211613624814
GMRes iteration 378, residual = 0.3273608498911106
GMRes iteration 379, residual = 0.32695576374471746
GMRes iteration 380, residual = 0.31297997102392366
GMRes iteration 381, residual = 0.31274532395637933
GMRes iteration 382, residual = 0.3028653123188283
GMRes iteration 383, residual = 0.30284415052314756
GMRes iteration 384, residual = 0.2988093684907523
GMRes iteration 385, residual = 0.2986662050374014
GMRes iteration 386, residual = 0.29509874168756506
GMRes iteration 387, residual = 0.29484073787962994
GMRes iteration 388, residual = 0.2929474985524429
GMRes iteration 389, residual = 0.29292563928129917
GMRes iteration 390, residual = 0.2911500712978511
GMRes iteration 391, residual = 0.29107088505130835
GMRes iteration 392, residual = 0.2894180936918138
GMRes iteration 393, residual = 0.28937215170000014
GMRes iteration 394, residual = 0.2885473201573094
GMRes iteration 395, residual = 0.2885204854476686
GMRes iteration 396, residual = 0.2880407821322877
GMRes iteration 397, residual = 0.28800264835922756
GMRes iteration 398, residual = 0.2870332893580639
GMRes iteration 399, residual = 0.2870172414239868
GMRes iteration 400, residual = 0.2858089796435719
GMRes iteration 401, residual = 0.2857748598007791
GMRes iteration 402, residual = 0.28456282134538746
GMRes iteration 403, residual = 0.28454736920043333
GMRes iteration 404, residual = 0.2829727172443052
GMRes iteration 405, residual = 0.2829315780197664
GMRes iteration 406, residual = 0.28170161690333917
GMRes iteration 407, residual = 0.28152310040166373
GMRes iteration 408, residual = 0.2804258258699756
GMRes iteration 409, residual = 0.2800895876205559
GMRes iteration 410, residual = 0.2789978707210909
GMRes iteration 411, residual = 0.2785547820485104
GMRes iteration 412, residual = 0.27684893283919493
GMRes iteration 413, residual = 0.27622770817992565
GMRes iteration 414, residual = 0.274994334906267
GMRes iteration 415, residual = 0.2745747037406607
GMRes iteration 416, residual = 0.2742116593075265
GMRes iteration 417, residual = 0.27344144943113996
GMRes iteration 418, residual = 0.27254835138103084
GMRes iteration 419, residual = 0.2722096668058126
GMRes iteration 420, residual = 0.271872629850898
GMRes iteration 421, residual = 0.2716720137014439
GMRes iteration 422, residual = 0.2714184807189578
GMRes iteration 423, residual = 0.27118315000713183
GMRes iteration 424, residual = 0.27082399379405786
GMRes iteration 425, residual = 0.2707786134518052
GMRes iteration 426, residual = 0.270682528614853
GMRes iteration 427, residual = 0.2700849775307015
GMRes iteration 428, residual = 0.26977637906406027
GMRes iteration 429, residual = 0.2688456957623627
GMRes iteration 430, residual = 0.26840489671731016
GMRes iteration 431, residual = 0.2673570963770047
GMRes iteration 432, residual = 0.2669729354190618
GMRes iteration 433, residual = 0.26416794204280064
GMRes iteration 434, residual = 0.26368296063266133
GMRes iteration 435, residual = 0.2609483976017693
GMRes iteration 436, residual = 0.26075030009601846
GMRes iteration 437, residual = 0.25973376711078267
GMRes iteration 438, residual = 0.2595715403859431
GMRes iteration 439, residual = 0.258538044427641
GMRes iteration 440, residual = 0.25835848127664823
GMRes iteration 441, residual = 0.2577624832494598
GMRes iteration 442, residual = 0.2575863922124834
GMRes iteration 443, residual = 0.2567547089735073
GMRes iteration 444, residual = 0.25636389399526205
GMRes iteration 445, residual = 0.2552030934907909
GMRes iteration 446, residual = 0.2551653758616645
GMRes iteration 447, residual = 0.2521337963622942
GMRes iteration 448, residual = 0.25184331420736455
GMRes iteration 449, residual = 0.24513508285220617
GMRes iteration 450, residual = 0.24476724878413844
GMRes iteration 451, residual = 0.23956034507698837
GMRes iteration 452, residual = 0.238877686666464
GMRes iteration 453, residual = 0.2371424305573599
GMRes iteration 454, residual = 0.23562628951722858
GMRes iteration 455, residual = 0.23371359798652025
GMRes iteration 456, residual = 0.23156386540048685
GMRes iteration 457, residual = 0.22960221065529302
GMRes iteration 458, residual = 0.22855414747657543
GMRes iteration 459, residual = 0.22713113490160117
GMRes iteration 460, residual = 0.22652917422894042
GMRes iteration 461, residual = 0.22595623683644428
GMRes iteration 462, residual = 0.2252441179313972
GMRes iteration 463, residual = 0.22517453425177528
GMRes iteration 464, residual = 0.2242533928892299
GMRes iteration 465, residual = 0.22425155051047754
GMRes iteration 466, residual = 0.22334111850882551
GMRes iteration 467, residual = 0.22334101349053137
GMRes iteration 468, residual = 0.22244913195007882
GMRes iteration 469, residual = 0.22214379648774502
GMRes iteration 470, residual = 0.22060528410541083
GMRes iteration 471, residual = 0.21858042886814996
GMRes iteration 472, residual = 0.21782577098740064
GMRes iteration 473, residual = 0.21524788615020737
GMRes iteration 474, residual = 0.21422344073415087
GMRes iteration 475, residual = 0.2113041532078068
GMRes iteration 476, residual = 0.2103423193920767
GMRes iteration 477, residual = 0.2077686210551783
GMRes iteration 478, residual = 0.20531783030749595
GMRes iteration 479, residual = 0.20384319345758098
GMRes iteration 480, residual = 0.20243949037992964
GMRes iteration 481, residual = 0.2021991414611724
GMRes iteration 482, residual = 0.20142363951557526
GMRes iteration 483, residual = 0.20130509623432877
GMRes iteration 484, residual = 0.20031049834734427
GMRes iteration 485, residual = 0.20028087418285237
GMRes iteration 486, residual = 0.19932347155705457
GMRes iteration 487, residual = 0.1992127763447957
GMRes iteration 488, residual = 0.1976625475229835
GMRes iteration 489, residual = 0.197659039655863
GMRes iteration 490, residual = 0.19663014102015344
GMRes iteration 491, residual = 0.1962192388980599
GMRes iteration 492, residual = 0.19475524981128856
GMRes iteration 493, residual = 0.19219747966967485
GMRes iteration 494, residual = 0.18701126187444952
GMRes iteration 495, residual = 0.18475562976654644
GMRes iteration 496, residual = 0.1799168293281288
GMRes iteration 497, residual = 0.178366287972596
GMRes iteration 498, residual = 0.1777455995369002
GMRes iteration 499, residual = 0.1771563547140096
GMRes iteration 500, residual = 0.17694724612222357
GMRes iteration 501, residual = 0.17557949494552866
GMRes iteration 502, residual = 0.17552644093163647
GMRes iteration 503, residual = 0.17387990755645652
GMRes iteration 504, residual = 0.1738570281150347
GMRes iteration 505, residual = 0.1728841203076171
GMRes iteration 506, residual = 0.1717737017949795
GMRes iteration 507, residual = 0.17127086848397557
GMRes iteration 508, residual = 0.16675848428494586
GMRes iteration 509, residual = 0.1667123281302717
GMRes iteration 510, residual = 0.16391148097240918
GMRes iteration 511, residual = 0.1639047255620484
GMRes iteration 512, residual = 0.16200305058173073
GMRes iteration 513, residual = 0.16142134914151407
GMRes iteration 514, residual = 0.16029731667752964
GMRes iteration 515, residual = 0.15826028092554498
GMRes iteration 516, residual = 0.15659514289346635
GMRes iteration 517, residual = 0.15284341564319515
GMRes iteration 518, residual = 0.15044422024799464
GMRes iteration 519, residual = 0.14853014046526497
GMRes iteration 520, residual = 0.14674087057347837
GMRes iteration 521, residual = 0.14658678336924647
GMRes iteration 522, residual = 0.14475940441148083
GMRes iteration 523, residual = 0.14465155096869828
GMRes iteration 524, residual = 0.14359145867483095
GMRes iteration 525, residual = 0.1432770641849256
GMRes iteration 526, residual = 0.1419494001160668
GMRes iteration 527, residual = 0.13718472599848602
GMRes iteration 528, residual = 0.13596545696586557
GMRes iteration 529, residual = 0.13122067529471035
GMRes iteration 530, residual = 0.13033785100588935
GMRes iteration 531, residual = 0.1265531927137526
GMRes iteration 532, residual = 0.1263915340278582
GMRes iteration 533, residual = 0.1246328893596629
GMRes iteration 534, residual = 0.1245230033494897
GMRes iteration 535, residual = 0.1233493087429893
GMRes iteration 536, residual = 0.1226559543715694
GMRes iteration 537, residual = 0.12103680352717203
GMRes iteration 538, residual = 0.11869227135649586
GMRes iteration 539, residual = 0.11596953915159117
GMRes iteration 540, residual = 0.11567315736160179
GMRes iteration 541, residual = 0.11400146481913485
GMRes iteration 542, residual = 0.11400034389243373
GMRes iteration 543, residual = 0.11361270017540888
GMRes iteration 544, residual = 0.11359539399201117
GMRes iteration 545, residual = 0.11346997295009105
GMRes iteration 546, residual = 0.11338209083265086
GMRes iteration 547, residual = 0.11312138448119198
GMRes iteration 548, residual = 0.1128960653696426
GMRes iteration 549, residual = 0.11184513097150939
GMRes iteration 550, residual = 0.11024415075266786
GMRes iteration 551, residual = 0.10809707876966007
GMRes iteration 552, residual = 0.10535899366895636
GMRes iteration 553, residual = 0.10364408754363796
GMRes iteration 554, residual = 0.10232834534193558
GMRes iteration 555, residual = 0.10212440890175328
GMRes iteration 556, residual = 0.10177301505966381
GMRes iteration 557, residual = 0.10173826474442337
GMRes iteration 558, residual = 0.10164521066532199
GMRes iteration 559, residual = 0.10156115776497825
GMRes iteration 560, residual = 0.10118462317848687
GMRes iteration 561, residual = 0.10117468471220227
GMRes iteration 562, residual = 0.10054443241540444
GMRes iteration 563, residual = 0.10052120673299672
GMRes iteration 564, residual = 0.09435738130350899
GMRes iteration 565, residual = 0.09428057874992228
GMRes iteration 566, residual = 0.0908816742660973
GMRes iteration 567, residual = 0.09083803956651468
GMRes iteration 568, residual = 0.09036409218529182
GMRes iteration 569, residual = 0.0896588530697943
GMRes iteration 570, residual = 0.08961506690585114
GMRes iteration 571, residual = 0.0882954190312161
GMRes iteration 572, residual = 0.0882825090807802
GMRes iteration 573, residual = 0.08767885172620998
GMRes iteration 574, residual = 0.08743594666803904
GMRes iteration 575, residual = 0.08697703205297248
GMRes iteration 576, residual = 0.0848472948113492
GMRes iteration 577, residual = 0.08426904326965434
GMRes iteration 578, residual = 0.08278226700082882
GMRes iteration 579, residual = 0.08261421138289517
GMRes iteration 580, residual = 0.08236538238963234
GMRes iteration 581, residual = 0.08180192544946917
GMRes iteration 582, residual = 0.08177728243638462
GMRes iteration 583, residual = 0.08060735883508682
GMRes iteration 584, residual = 0.08041542592386629
GMRes iteration 585, residual = 0.079855531746912
GMRes iteration 586, residual = 0.07983180595715747
GMRes iteration 587, residual = 0.0787894792400113
GMRes iteration 588, residual = 0.07810238365717222
GMRes iteration 589, residual = 0.07420756538443582
GMRes iteration 590, residual = 0.07351087662535517
GMRes iteration 591, residual = 0.07120050743276317
GMRes iteration 592, residual = 0.0711415041472412
GMRes iteration 593, residual = 0.0709557362010502
GMRes iteration 594, residual = 0.07093683511767299
GMRes iteration 595, residual = 0.07072138336927913
GMRes iteration 596, residual = 0.0707209679583454
GMRes iteration 597, residual = 0.07005843761163243
GMRes iteration 598, residual = 0.0698853110454026
GMRes iteration 599, residual = 0.06720480026839883
GMRes iteration 600, residual = 0.06472305246926092
GMRes iteration 601, residual = 0.06443263946036697
GMRes iteration 602, residual = 0.06424841748057986
GMRes iteration 603, residual = 0.06422368397429185
GMRes iteration 604, residual = 0.06413055763221762
GMRes iteration 605, residual = 0.0637460665458121
GMRes iteration 606, residual = 0.06275076287168357
GMRes iteration 607, residual = 0.05876588743822288
GMRes iteration 608, residual = 0.05832631129043257
GMRes iteration 609, residual = 0.05830059025711854
GMRes iteration 610, residual = 0.058265096788742164
GMRes iteration 611, residual = 0.052387821901912525
GMRes iteration 612, residual = 0.04917371631316275
GMRes iteration 613, residual = 0.04917334933154909
GMRes iteration 614, residual = 0.048940755293880145
GMRes iteration 615, residual = 0.044343516681711104
GMRes iteration 616, residual = 0.040416672495725
GMRes iteration 617, residual = 0.04037511892591613
GMRes iteration 618, residual = 0.040367165718826606
GMRes iteration 619, residual = 0.039807155801954416
GMRes iteration 620, residual = 0.03178672305344636
GMRes iteration 621, residual = 0.031775891213650806
GMRes iteration 622, residual = 0.02747499735148958
GMRes iteration 623, residual = 0.027469500417076076
GMRes iteration 624, residual = 0.02101182098649043
GMRes iteration 625, residual = 0.02093813004365113
GMRes iteration 626, residual = 0.017380402728779134
GMRes iteration 627, residual = 0.017211756556916428
GMRes iteration 628, residual = 0.01333789527707615
GMRes iteration 629, residual = 0.013096144604082054
GMRes iteration 630, residual = 0.011106675468781923
GMRes iteration 631, residual = 0.010968567434000286
GMRes iteration 632, residual = 0.009090277904527535
GMRes iteration 633, residual = 0.0090230604122296
GMRes iteration 634, residual = 0.007852612333539146
GMRes iteration 635, residual = 0.007684597268645358
GMRes iteration 636, residual = 0.006306822200556493
GMRes iteration 637, residual = 0.006290971531691319
GMRes iteration 638, residual = 0.005064739907370749
GMRes iteration 639, residual = 0.005063866503691141
GMRes iteration 640, residual = 0.004141723505134663
GMRes iteration 641, residual = 0.003976134137587752
GMRes iteration 642, residual = 0.0035225374184612154
GMRes iteration 643, residual = 0.0035172420592543326
GMRes iteration 644, residual = 0.002403847011383781
GMRes iteration 645, residual = 0.0023877234454769227
GMRes iteration 646, residual = 0.0018990507659398328
GMRes iteration 647, residual = 0.0018925425943039993
GMRes iteration 648, residual = 0.0015281194696384896
GMRes iteration 649, residual = 0.001522027207872133
GMRes iteration 650, residual = 0.0012060141902295461
GMRes iteration 651, residual = 0.0011851549883525982
GMRes iteration 652, residual = 0.0009542555251787795
GMRes iteration 653, residual = 0.0009482489835218494
GMRes iteration 654, residual = 0.0007685433238801558
GMRes iteration 655, residual = 0.0007674665278888415
GMRes iteration 656, residual = 0.0005928504266472813
GMRes iteration 657, residual = 0.0005924114772339488
GMRes iteration 658, residual = 0.0004824162369247501
GMRes iteration 659, residual = 0.0004813816953866833
GMRes iteration 660, residual = 0.00039083142404151665
GMRes iteration 661, residual = 0.00038960531095737905
GMRes iteration 662, residual = 0.00031157068382865684
GMRes iteration 663, residual = 0.00030802397206851654
GMRes iteration 664, residual = 0.0002500800241902727
GMRes iteration 665, residual = 0.0002492440882894102
GMRes iteration 666, residual = 0.00019493051112650476
GMRes iteration 667, residual = 0.0001803390518294062
GMRes iteration 668, residual = 0.00014733696531070022
GMRes iteration 669, residual = 0.0001431490633147749
GMRes iteration 670, residual = 0.00011407976299645825
GMRes iteration 671, residual = 0.00010736480161074876
GMRes iteration 672, residual = 8.573240870918364e-05
GMRes iteration 673, residual = 8.543709732597381e-05
GMRes iteration 674, residual = 6.691185368277701e-05
GMRes iteration 675, residual = 6.66696303307648e-05
GMRes iteration 676, residual = 5.3992672927675206e-05
GMRes iteration 677, residual = 4.8700198309469355e-05
GMRes iteration 678, residual = 4.0969944952722224e-05
GMRes iteration 679, residual = 3.8470723553031335e-05
GMRes iteration 680, residual = 3.2194351665232915e-05
GMRes iteration 681, residual = 3.055219448800798e-05
GMRes iteration 682, residual = 2.4327440044773192e-05
GMRes iteration 683, residual = 2.35343041060687e-05
GMRes iteration 684, residual = 1.9141575333058215e-05
GMRes iteration 685, residual = 1.8321669714437562e-05
GMRes iteration 686, residual = 1.6032066038626422e-05
GMRes iteration 687, residual = 1.5567685629711072e-05
GMRes iteration 688, residual = 1.2735880858823126e-05
GMRes iteration 689, residual = 1.1724234405675214e-05
GMRes iteration 690, residual = 9.496909104890058e-06
GMRes iteration 691, residual = 8.748053974823705e-06
GMRes iteration 692, residual = 7.106907053740631e-06
GMRes iteration 693, residual = 7.0249855431573875e-06
GMRes iteration 694, residual = 5.445586837674106e-06
GMRes iteration 695, residual = 5.423018644850816e-06
GMRes iteration 696, residual = 4.206833371592613e-06
GMRes iteration 697, residual = 4.185774923734363e-06
GMRes iteration 698, residual = 3.4069654758828417e-06
GMRes iteration 699, residual = 3.3716154571963974e-06
GMRes iteration 700, residual = 2.624988973020391e-06
GMRes iteration 701, residual = 2.5804850260526183e-06
GMRes iteration 702, residual = 2.0965730381792133e-06
GMRes iteration 703, residual = 1.9859599830974252e-06
GMRes iteration 704, residual = 1.6708215812718427e-06
GMRes iteration 705, residual = 1.5749818236399247e-06
GMRes iteration 706, residual = 1.2013053607195803e-06
GMRes iteration 707, residual = 1.0966901288815747e-06
GMRes iteration 708, residual = 7.830399843234347e-07
GMRes iteration 709, residual = 7.35778952891837e-07
GMRes iteration 710, residual = 5.551026648539047e-07
GMRes iteration 711, residual = 5.507165726605375e-07
GMRes iteration 712, residual = 4.139785289953667e-07
GMRes iteration 713, residual = 4.0800549139797e-07
GMRes iteration 714, residual = 3.3635678375358203e-07
GMRes iteration 715, residual = 3.3123773708543883e-07
GMRes iteration 716, residual = 2.494435390891933e-07
GMRes iteration 717, residual = 2.412983607068223e-07
GMRes iteration 718, residual = 1.9388407544859712e-07
GMRes iteration 719, residual = 1.8913405399422878e-07
GMRes iteration 720, residual = 1.372797092412445e-07
GMRes iteration 721, residual = 1.3542061153849683e-07
GMRes iteration 722, residual = 1.0358198705245888e-07
GMRes iteration 723, residual = 9.984490773934667e-08
GMRes iteration 724, residual = 7.781027848566133e-08
GMRes iteration 725, residual = 7.751730638011497e-08
GMRes iteration 726, residual = 5.990578541031411e-08
GMRes iteration 727, residual = 5.833528603159932e-08
GMRes iteration 728, residual = 4.578923297599001e-08
GMRes iteration 729, residual = 4.5025971340150395e-08
GMRes iteration 730, residual = 3.2744630578177445e-08
GMRes iteration 731, residual = 3.2732274417705146e-08
GMRes iteration 732, residual = 2.3446301883279962e-08
GMRes iteration 733, residual = 2.2461195675241048e-08
GMRes iteration 734, residual = 1.669307164057589e-08
GMRes iteration 735, residual = 1.574175712298691e-08
GMRes iteration 736, residual = 1.2200285986606366e-08
GMRes iteration 737, residual = 1.2060310752976231e-08
GMRes iteration 738, residual = 8.67426687466792e-09
Compute the Error
Have a look at Neumann data on \(\Gamma_0\) and compute the \(L_2\)-error:
# compare with exact Dirichlet trace
gf = GridFunction(fes)
gf.vec[:] = sol
gfdiv = gf.components[1]
print ("Neumann error =", sqrt(Integrate(Norm(gfdiv + mj.components[1] - mjexa.components[1])**2, mesh.Boundaries(".*"), BND)))
Draw (Norm(gfdiv), mesh.Boundaries("dirichlet"), order=3);
Neumann error = 0.016204016503892308
Have a look at Dirichlet data on \(\Gamma_1\) and compute the \(L_2\)-error:
# compare with the exact Neuman trace
gfcurl = gf.components[0]
print ("Dirichlet error =", sqrt(Integrate(Norm(gfcurl +mj.components[0] - mjexa.components[0])**2, mesh.Boundaries(".*"), BND)))
Draw (Norm(gfcurl), mesh.Boundaries("neumann"), draw_vol=False, order=3);
Dirichlet error = 0.02859337100450199
References (theory results):
Boundary Element Methods for Maxwell Transmission Problems in Lipschitz Domains
Note that the Calderon-projector in the referenced paper is slightly different from ours as only trace space \(H^{-\frac12}(\mathrm{div}_\Gamma,\Gamma)\) is used. Theoretical results stay the same.