Laplace mixed boundary values#

keys: homogeneous Laplace bvp, Calderon projector, Dirichlet data, Neumann data

from ngsolve import *
from netgen.occ import *
from ngsolve.krylovspace import CG, GMRes
from ngsolve.webgui import Draw
from libbem import *
Loading ngbem library

We consider an interior boundary value problem with mixed boundary conditions like this:

\( \begin{array}{r rcl r} & \Delta u &=& 0 &\mathrm{in}\; \Omega\,,\\ \textnormal{Dirichlet condition} & \gamma_0 u &=& u_0 & \mathrm{on}\; \Gamma_0\,,\\ \textnormal{Neumann condition} & \gamma_1 u &=& u_1 & \mathrm{on}\; \Gamma_1\,. \end{array} \)

\(\quad\quad\quad\)

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"
botsphere.faces.name = "dirichlet"
shape = Fuse( [topsphere,botsphere] )

order = 3
mesh = Mesh(OCCGeometry(shape).GenerateMesh(maxh=0.25)).Curve(order)
Draw (mesh);

The following representation formula for the solution \(u\) holds:

\[ x \in \Omega: \quad u(x) = \displaystyle{ \int\limits_\Gamma} \displaystyle{\frac{1}{4\,\pi}\, \frac{1}{\| x-y\|} } \, \gamma_1 u (y)\, \mathrm{d}\sigma_y - \displaystyle{ \int\limits_\Gamma} \displaystyle{\frac{1}{4\,\pi}\, \frac{\langle n(y) , x-y\rangle }{\| x-y\|^3} } \, \gamma_0 u (y)\, \mathrm{d}\sigma_y\,. \]

The Calderon projector relates the Dirichlet and the Neumann traces of the solution \(u\), i.e.,

\[\begin{split} \left( \begin{array}{c} \gamma_0 u \\ \gamma_1 u \end{array}\right) = \left( \begin{array}{cc} V & \frac12 - K \\ \frac12 + K^\intercal & D \end{array} \right) \left( \begin{array}{c} \gamma_1 u \\ \gamma_0 u \end{array}\right)\,, \end{split}\]

and we use it to solve for the Dirchlet trace on \(\Gamma_1\) and the Neumann trace on \(\Gamma_0\).

Define the finite element spaces for the Dirichlet and the Neumann data:

# use L2 conform elements for the Neumann trace, label all dofs where Neumann data is given:
fesL2 = SurfaceL2(mesh, order=order-1, dirichlet="neumann")
u,v = fesL2.TnT()

# use H1 conform elements for the Dirichlet trace, label all dofs where Dirichlet data is given:
fesH1 = H1(mesh, order=order, dirichlet="dirichlet", definedon=mesh.Boundaries(".*"))
u1,v1 = fesH1.TnT()

print ("ndofL2 = ", fesL2.ndof, "ndofH1 = ", fesH1.ndof)
ndofL2 =  2736 ndofH1 =  2128

Compute and set the Dirichlet and the Neumann data:

uexa = CF(x)
ud = GridFunction(fesH1)
ud.Interpolate(uexa, definedon=mesh.Boundaries("dirichlet"))

n = specialcf.normal(3)
gradn_uexa = CF((uexa.Diff(x), uexa.Diff(y), uexa.Diff(z))) * n
un = GridFunction(fesL2)
un.Interpolate(gradn_uexa, definedon=mesh.Boundaries("neumann"))

Draw(ud, mesh, draw_vol=False);

Compute boundary integral operators \(\mathrm{V}, \mathrm{K}, \mathrm{W}\) and the mass matrix \(\mathrm{M}\):

eps = 1e-6
intorder = 2 * order + 6
with TaskManager():
    V = SingleLayerPotentialOperator(fesL2, intorder=intorder, eps=eps)
    K = DoubleLayerPotentialOperator(fesH1, fesL2, intorder=intorder, eps=eps)
    W = HypersingularOperator(fesH1, intorder=intorder, eps=eps)
    M = BilinearForm(u1.Trace() * v.Trace() * ds(bonus_intorder=3)).Assemble()

Insert all given data in the Dirichlet-to-Neumann map and compute the right hand side vector:

with TaskManager():
    fd = ((0.5 * M.mat + K.mat) * ud.vec - V.mat * un.vec).Evaluate()
    fn = ((0.5 * M.mat.T - K.mat.T) * un.vec - W.mat * ud.vec).Evaluate()

Solve for the missing trace data:

with TaskManager():
    pred = BilinearForm(u.Trace() * v.Trace() * ds(bonus_intorder=3), check_unused=False).Assemble()
    pren = BilinearForm(u1.Trace() * v1.Trace() * ds(bonus_intorder=3), check_unused=False).Assemble()

lhs = BlockMatrix([[V.mat, - K.mat], [K.mat.T, W.mat]])
rhs = BlockVector([fd, fn])
pre = BlockMatrix ([[pred.mat.Inverse(freedofs=fesL2.FreeDofs()), None], [None, pren.mat.Inverse(freedofs=fesH1.FreeDofs())]])

sol = GMRes(A=lhs, b=rhs, pre=pre, tol=1e-8, maxsteps=500)
GMRes iteration 1, residual = 109.56977428387343     
GMRes iteration 2, residual = 53.578088986467456     
GMRes iteration 3, residual = 20.02971529993236     
GMRes iteration 4, residual = 7.7263872104400155     
GMRes iteration 5, residual = 5.904571084454537     
GMRes iteration 6, residual = 5.2252761934918395     
GMRes iteration 7, residual = 4.695885280589776     
GMRes iteration 8, residual = 3.845742237212549     
GMRes iteration 9, residual = 2.934693441298961     
GMRes iteration 10, residual = 2.147365296182972     
GMRes iteration 11, residual = 1.7984822439562025     
GMRes iteration 12, residual = 1.6096599099508464     
GMRes iteration 13, residual = 1.4221319492240632     
GMRes iteration 14, residual = 1.2833854974484056     
GMRes iteration 15, residual = 1.1815881666909398     
GMRes iteration 16, residual = 1.0835527239355542     
GMRes iteration 17, residual = 0.9660125018594619     
GMRes iteration 18, residual = 0.8510309673852604     
GMRes iteration 19, residual = 0.7248868740095112     
GMRes iteration 20, residual = 0.5415679170954206     
GMRes iteration 21, residual = 0.3602651898065477     
GMRes iteration 22, residual = 0.24574930210050622     
GMRes iteration 23, residual = 0.16671279818692858     
GMRes iteration 24, residual = 0.11627420287939133     
GMRes iteration 25, residual = 0.08564097546140019     
GMRes iteration 26, residual = 0.0758561831109907     
GMRes iteration 27, residual = 0.07280030611206426     
GMRes iteration 28, residual = 0.07111902237833044     
GMRes iteration 29, residual = 0.06974609122355316     
GMRes iteration 30, residual = 0.06677310051883685     
GMRes iteration 31, residual = 0.05875038157473041     
GMRes iteration 32, residual = 0.04548871576559034     
GMRes iteration 33, residual = 0.03231955789774897     
GMRes iteration 34, residual = 0.023750620126478925     
GMRes iteration 35, residual = 0.020036457915274163     
GMRes iteration 36, residual = 0.018765556446985618     
GMRes iteration 37, residual = 0.01818819127464733     
GMRes iteration 38, residual = 0.017664650298430284     
GMRes iteration 39, residual = 0.01732086232269319     
GMRes iteration 40, residual = 0.01707632550640583     
GMRes iteration 41, residual = 0.01656318703744669     
GMRes iteration 42, residual = 0.016004452831295424     
GMRes iteration 43, residual = 0.015337970331199826     
GMRes iteration 44, residual = 0.013876278672407328     
GMRes iteration 45, residual = 0.012477739890169466     
GMRes iteration 46, residual = 0.011011806260897724     
GMRes iteration 47, residual = 0.00925126073303309     
GMRes iteration 48, residual = 0.008024978559541034     
GMRes iteration 49, residual = 0.007398338859470811     
GMRes iteration 50, residual = 0.00709056706095909     
GMRes iteration 51, residual = 0.006951365559012586     
GMRes iteration 52, residual = 0.006886446539741309     
GMRes iteration 53, residual = 0.006835147742570182     
GMRes iteration 54, residual = 0.0067461082459643286     
GMRes iteration 55, residual = 0.006528485457066561     
GMRes iteration 56, residual = 0.006065521233534668     
GMRes iteration 57, residual = 0.005512292482027019     
GMRes iteration 58, residual = 0.0051718733154409785     
GMRes iteration 59, residual = 0.0049687421511334495     
GMRes iteration 60, residual = 0.004768523235849397     
GMRes iteration 61, residual = 0.004482629863739388     
GMRes iteration 62, residual = 0.004153366193813053     
GMRes iteration 63, residual = 0.0039036545920718836     
GMRes iteration 64, residual = 0.00373164607519697     
GMRes iteration 65, residual = 0.0035920096499237864     
GMRes iteration 66, residual = 0.003433144943178115     
GMRes iteration 67, residual = 0.0032991632011996753     
GMRes iteration 68, residual = 0.003203710531547274     
GMRes iteration 69, residual = 0.003132116792552964     
GMRes iteration 70, residual = 0.003056409131120033     
GMRes iteration 71, residual = 0.002968547780125638     
GMRes iteration 72, residual = 0.002913940232309002     
GMRes iteration 73, residual = 0.002882554120250164     
GMRes iteration 74, residual = 0.00283902399392659     
GMRes iteration 75, residual = 0.002773799962717094     
GMRes iteration 76, residual = 0.002724327798536615     
GMRes iteration 77, residual = 0.002674270630770131     
GMRes iteration 78, residual = 0.0025934001721240234     
GMRes iteration 79, residual = 0.002502067561479456     
GMRes iteration 80, residual = 0.002419853792526647     
GMRes iteration 81, residual = 0.0023256128688907774     
GMRes iteration 82, residual = 0.0021279539514550113     
GMRes iteration 83, residual = 0.0018496195277642106     
GMRes iteration 84, residual = 0.00164904433874347     
GMRes iteration 85, residual = 0.0014651241182137545     
GMRes iteration 86, residual = 0.0013769601603737971     
GMRes iteration 87, residual = 0.0013527776251833132     
GMRes iteration 88, residual = 0.0013342752661683944     
GMRes iteration 89, residual = 0.001320713512216962     
GMRes iteration 90, residual = 0.001310138799080009     
GMRes iteration 91, residual = 0.0012964795481426358     
GMRes iteration 92, residual = 0.0012800685821756999     
GMRes iteration 93, residual = 0.0012555716118835116     
GMRes iteration 94, residual = 0.0012050893547945944     
GMRes iteration 95, residual = 0.001129063177477714     
GMRes iteration 96, residual = 0.001013674388941282     
GMRes iteration 97, residual = 0.0008735506743615239     
GMRes iteration 98, residual = 0.0007328742896250309     
GMRes iteration 99, residual = 0.0006221024427603799     
GMRes iteration 100, residual = 0.00057304127249464     
GMRes iteration 101, residual = 0.0005531219034365028     
GMRes iteration 102, residual = 0.0005448435654116934     
GMRes iteration 103, residual = 0.0005413806361590752     
GMRes iteration 104, residual = 0.0005384367090834324     
GMRes iteration 105, residual = 0.0005335287513642109     
GMRes iteration 106, residual = 0.0005207662185522647     
GMRes iteration 107, residual = 0.0004931003788896067     
GMRes iteration 108, residual = 0.0004312925508196811     
GMRes iteration 109, residual = 0.00035859605989118395     
GMRes iteration 110, residual = 0.0003124861255359668     
GMRes iteration 111, residual = 0.00029677393117933136     
GMRes iteration 112, residual = 0.0002905735382176985     
GMRes iteration 113, residual = 0.00028677160524709746     
GMRes iteration 114, residual = 0.000281258424847006     
GMRes iteration 115, residual = 0.0002699616697054488     
GMRes iteration 116, residual = 0.0002554968801013063     
GMRes iteration 117, residual = 0.00024122834699909213     
GMRes iteration 118, residual = 0.00023166813931738207     
GMRes iteration 119, residual = 0.0002257610019030631     
GMRes iteration 120, residual = 0.00022107896617206094     
GMRes iteration 121, residual = 0.00021803817056023986     
GMRes iteration 122, residual = 0.00021636170364115287     
GMRes iteration 123, residual = 0.00021517045893728503     
GMRes iteration 124, residual = 0.00021287723413407683     
GMRes iteration 125, residual = 0.0002071670865210616     
GMRes iteration 126, residual = 0.00019270421406906112     
GMRes iteration 127, residual = 0.00017473174289475373     
GMRes iteration 128, residual = 0.00015861044030472958     
GMRes iteration 129, residual = 0.00015080016323220017     
GMRes iteration 130, residual = 0.0001462218453872712     
GMRes iteration 131, residual = 0.00014338115451992622     
GMRes iteration 132, residual = 0.00014178309224767747     
GMRes iteration 133, residual = 0.00013940001917242865     
GMRes iteration 134, residual = 0.00013321179606264177     
GMRes iteration 135, residual = 0.00012059801866744118     
GMRes iteration 136, residual = 0.00010686892638408754     
GMRes iteration 137, residual = 0.00010017975479853343     
GMRes iteration 138, residual = 9.730692558440012e-05     
GMRes iteration 139, residual = 9.507084509285511e-05     
GMRes iteration 140, residual = 9.360492974615095e-05     
GMRes iteration 141, residual = 9.267813122211062e-05     
GMRes iteration 142, residual = 9.205029970256798e-05     
GMRes iteration 143, residual = 9.148690261953726e-05     
GMRes iteration 144, residual = 9.083882540847361e-05     
GMRes iteration 145, residual = 8.992648426783839e-05     
GMRes iteration 146, residual = 8.77275165131721e-05     
GMRes iteration 147, residual = 8.194963463980959e-05     
GMRes iteration 148, residual = 7.259970096266416e-05     
GMRes iteration 149, residual = 6.380993287204227e-05     
GMRes iteration 150, residual = 5.6212194916275296e-05     
GMRes iteration 151, residual = 5.1297401074797614e-05     
GMRes iteration 152, residual = 4.910219476177334e-05     
GMRes iteration 153, residual = 4.845710202115158e-05     
GMRes iteration 154, residual = 4.823089732998881e-05     
GMRes iteration 155, residual = 4.805289160581213e-05     
GMRes iteration 156, residual = 4.7824873526917997e-05     
GMRes iteration 157, residual = 4.735477623671503e-05     
GMRes iteration 158, residual = 4.611097727366553e-05     
GMRes iteration 159, residual = 4.300893531163364e-05     
GMRes iteration 160, residual = 4.0562658923454674e-05     
GMRes iteration 161, residual = 3.962909895924205e-05     
GMRes iteration 162, residual = 3.9165646783051706e-05     
GMRes iteration 163, residual = 3.888716090733904e-05     
GMRes iteration 164, residual = 3.856838712954146e-05     
GMRes iteration 165, residual = 3.749562493306458e-05     
GMRes iteration 166, residual = 3.450460875622774e-05     
GMRes iteration 167, residual = 3.1483739020182964e-05     
GMRes iteration 168, residual = 2.977485011115027e-05     
GMRes iteration 169, residual = 2.8781280544594152e-05     
GMRes iteration 170, residual = 2.841118392169541e-05     
GMRes iteration 171, residual = 2.809561594134719e-05     
GMRes iteration 172, residual = 2.790574817114556e-05     
GMRes iteration 173, residual = 2.7817196373760218e-05     
GMRes iteration 174, residual = 2.7724767859045117e-05     
GMRes iteration 175, residual = 2.7522790432006366e-05     
GMRes iteration 176, residual = 2.690507430624193e-05     
GMRes iteration 177, residual = 2.5395473029836956e-05     
GMRes iteration 178, residual = 2.3820225726383374e-05     
GMRes iteration 179, residual = 2.3151593387874277e-05     
GMRes iteration 180, residual = 2.2886477447236217e-05     
GMRes iteration 181, residual = 2.2766872794685194e-05     
GMRes iteration 182, residual = 2.2628020710578327e-05     
GMRes iteration 183, residual = 2.2392803008878936e-05     
GMRes iteration 184, residual = 2.194027901400437e-05     
GMRes iteration 185, residual = 2.0867344793669864e-05     
GMRes iteration 186, residual = 1.8593058046349692e-05     
GMRes iteration 187, residual = 1.5431318245309762e-05     
GMRes iteration 188, residual = 1.3608477365707445e-05     
GMRes iteration 189, residual = 1.3030763926547275e-05     
GMRes iteration 190, residual = 1.2892132798946974e-05     
GMRes iteration 191, residual = 1.2833211134322818e-05     
GMRes iteration 192, residual = 1.2768960092100174e-05     
GMRes iteration 193, residual = 1.2660414808311673e-05     
GMRes iteration 194, residual = 1.244308687653562e-05     
GMRes iteration 195, residual = 1.2020533987830851e-05     
GMRes iteration 196, residual = 1.1273119252995644e-05     
GMRes iteration 197, residual = 1.0287353963177704e-05     
GMRes iteration 198, residual = 9.48768008962981e-06     
GMRes iteration 199, residual = 9.121838307200509e-06     
GMRes iteration 200, residual = 8.963639095769533e-06     
GMRes iteration 201, residual = 8.876033846632993e-06     
GMRes iteration 202, residual = 8.81269913488846e-06     
GMRes iteration 203, residual = 8.712832627813084e-06     
GMRes iteration 204, residual = 8.396911224890718e-06     
GMRes iteration 205, residual = 7.213177838100287e-06     
GMRes iteration 206, residual = 5.5862934465884875e-06     
GMRes iteration 207, residual = 4.749618915342769e-06     
GMRes iteration 208, residual = 4.398703197417364e-06     
GMRes iteration 209, residual = 4.301314556803851e-06     
GMRes iteration 210, residual = 4.258754613133445e-06     
GMRes iteration 211, residual = 4.238099623095168e-06     
GMRes iteration 212, residual = 4.226490936566946e-06     
GMRes iteration 213, residual = 4.203000645552451e-06     
GMRes iteration 214, residual = 4.12313990857982e-06     
GMRes iteration 215, residual = 3.952928327603729e-06     
GMRes iteration 216, residual = 3.636817402007184e-06     
GMRes iteration 217, residual = 3.2137155667482297e-06     
GMRes iteration 218, residual = 2.825810347939307e-06     
GMRes iteration 219, residual = 2.5961710691969028e-06     
GMRes iteration 220, residual = 2.5006113612396026e-06     
GMRes iteration 221, residual = 2.4487086839501255e-06     
GMRes iteration 222, residual = 2.42485872410723e-06     
GMRes iteration 223, residual = 2.404700243207949e-06     
GMRes iteration 224, residual = 2.361413238520266e-06     
GMRes iteration 225, residual = 2.2085520954291453e-06     
GMRes iteration 226, residual = 1.9261843019901497e-06     
GMRes iteration 227, residual = 1.7505255439150946e-06     
GMRes iteration 228, residual = 1.6863026989666392e-06     
GMRes iteration 229, residual = 1.6595696301780483e-06     
GMRes iteration 230, residual = 1.6488654182748682e-06     
GMRes iteration 231, residual = 1.642480214397099e-06     
GMRes iteration 232, residual = 1.6342475245706617e-06     
GMRes iteration 233, residual = 1.6189638817212505e-06     
GMRes iteration 234, residual = 1.5907892695972014e-06     
GMRes iteration 235, residual = 1.50570722326555e-06     
GMRes iteration 236, residual = 1.321737940582364e-06     
GMRes iteration 237, residual = 1.2154174107189966e-06     
GMRes iteration 238, residual = 1.1885513196108266e-06     
GMRes iteration 239, residual = 1.1790368383613325e-06     
GMRes iteration 240, residual = 1.1725314245873247e-06     
GMRes iteration 241, residual = 1.1618976498956049e-06     
GMRes iteration 242, residual = 1.1259289083842638e-06     
GMRes iteration 243, residual = 1.0071739829132192e-06     
GMRes iteration 244, residual = 8.851221473432089e-07     
GMRes iteration 245, residual = 8.044450385452296e-07     
GMRes iteration 246, residual = 7.454426644853448e-07     
GMRes iteration 247, residual = 7.291711541165116e-07     
GMRes iteration 248, residual = 7.249716113210452e-07     
GMRes iteration 249, residual = 7.2263963508252e-07     
GMRes iteration 250, residual = 7.201730344167685e-07     
GMRes iteration 251, residual = 7.116041196993898e-07     
GMRes iteration 252, residual = 6.751986887123247e-07     
GMRes iteration 253, residual = 6.103247393605541e-07     
GMRes iteration 254, residual = 5.64740838126283e-07     
GMRes iteration 255, residual = 5.513634891186011e-07     
GMRes iteration 256, residual = 5.465598495011697e-07     
GMRes iteration 257, residual = 5.432479363238409e-07     
GMRes iteration 258, residual = 5.40517226165194e-07     
GMRes iteration 259, residual = 5.369780136024469e-07     
GMRes iteration 260, residual = 5.281700100575385e-07     
GMRes iteration 261, residual = 5.010582431872546e-07     
GMRes iteration 262, residual = 4.39161835906652e-07     
GMRes iteration 263, residual = 3.739983477676239e-07     
GMRes iteration 264, residual = 3.358182324443866e-07     
GMRes iteration 265, residual = 3.192524544414752e-07     
GMRes iteration 266, residual = 3.1313457867431974e-07     
GMRes iteration 267, residual = 3.1097162079984716e-07     
GMRes iteration 268, residual = 3.0924454645012883e-07     
GMRes iteration 269, residual = 3.0521397984019766e-07     
GMRes iteration 270, residual = 2.913232356909367e-07     
GMRes iteration 271, residual = 2.698203696954536e-07     
GMRes iteration 272, residual = 2.5903239922270553e-07     
GMRes iteration 273, residual = 2.5530488655287477e-07     
GMRes iteration 274, residual = 2.538321027827013e-07     
GMRes iteration 275, residual = 2.5285643022733265e-07     
GMRes iteration 276, residual = 2.5107350672819025e-07     
GMRes iteration 277, residual = 2.45729976531329e-07     
GMRes iteration 278, residual = 2.26070537524473e-07     
GMRes iteration 279, residual = 1.863259496771153e-07     
GMRes iteration 280, residual = 1.6189432390388383e-07     
GMRes iteration 281, residual = 1.54042473033192e-07     
GMRes iteration 282, residual = 1.5130987762260686e-07     
GMRes iteration 283, residual = 1.4981113279408446e-07     
GMRes iteration 284, residual = 1.4873696971754102e-07     
GMRes iteration 285, residual = 1.4671880088356592e-07     
GMRes iteration 286, residual = 1.4008973674925218e-07     
GMRes iteration 287, residual = 1.27702493702094e-07     
GMRes iteration 288, residual = 1.0702141400568658e-07     
GMRes iteration 289, residual = 9.175872738783187e-08     
GMRes iteration 290, residual = 8.647687877271999e-08     
GMRes iteration 291, residual = 8.493890810817643e-08     
GMRes iteration 292, residual = 8.463404925805552e-08     
GMRes iteration 293, residual = 8.44032644576603e-08     
GMRes iteration 294, residual = 8.368825549356162e-08     
GMRes iteration 295, residual = 8.134248245932006e-08     
GMRes iteration 296, residual = 7.499681896402192e-08     
GMRes iteration 297, residual = 6.571969807542046e-08     
GMRes iteration 298, residual = 5.961345334184864e-08     
GMRes iteration 299, residual = 5.739359730886974e-08     
GMRes iteration 300, residual = 5.677529668027571e-08     
GMRes iteration 301, residual = 5.6373895404074145e-08     
GMRes iteration 302, residual = 5.589758037690308e-08     
GMRes iteration 303, residual = 5.481047304579451e-08     
GMRes iteration 304, residual = 5.0623741985508214e-08     
GMRes iteration 305, residual = 4.0700325365337355e-08     
GMRes iteration 306, residual = 3.486679695681286e-08     
GMRes iteration 307, residual = 3.329060903310628e-08     
GMRes iteration 308, residual = 3.291509287866332e-08     
GMRes iteration 309, residual = 3.273642803780152e-08     
GMRes iteration 310, residual = 3.258679443456137e-08     
GMRes iteration 311, residual = 3.228919307518225e-08     
GMRes iteration 312, residual = 3.154634444249595e-08     
GMRes iteration 313, residual = 2.984244021170185e-08     
GMRes iteration 314, residual = 2.6288360580884104e-08     
GMRes iteration 315, residual = 2.3372534917595815e-08     
GMRes iteration 316, residual = 2.248993123815727e-08     
GMRes iteration 317, residual = 2.2265598407014498e-08     
GMRes iteration 318, residual = 2.2161941491037488e-08     
GMRes iteration 319, residual = 2.2098288495216932e-08     
GMRes iteration 320, residual = 2.1980859778677595e-08     
GMRes iteration 321, residual = 2.1434256977448923e-08     
GMRes iteration 322, residual = 1.9612451532793384e-08     
GMRes iteration 323, residual = 1.5901327207379285e-08     
GMRes iteration 324, residual = 1.3451447960219688e-08     
GMRes iteration 325, residual = 1.2935341449011167e-08     
GMRes iteration 326, residual = 1.2831404731980754e-08     
GMRes iteration 327, residual = 1.2785627749787023e-08     
GMRes iteration 328, residual = 1.2713858638499269e-08     
GMRes iteration 329, residual = 1.2341982790999762e-08     
GMRes iteration 330, residual = 1.0867695326295311e-08     
GMRes iteration 331, residual = 9.512413411196047e-09     

Have a look at the Neumann data on \(\Gamma_0\) and compute the error:

gfu = GridFunction(fesL2)
gfu.vec[:] = sol[0]
print ("L2 error of surface gradients =", sqrt(Integrate((gradn_uexa - gfu -un)**2, mesh.Boundaries(".*"), BND)))
Draw(gfu)
L2 error of surface gradients = 0.0006228098891096159
BaseWebGuiScene

Have a look at the Dirichlet data on \(\Gamma_1\) and compute the error:

gfu = GridFunction(fesH1)
gfu.vec[:] = sol[1]
print ("L2 error Dirichlet =", sqrt(Integrate((uexa - gfu -ud)**2, mesh.Boundaries(".*"), BND)))
Draw(gfu)
L2 error Dirichlet = 9.200427517814712e-06
BaseWebGuiScene