Laplace with Mixed Conditions

Laplace with Mixed Conditions#

keys: homogeneous Laplace bvp with mixed conditions, 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 ngbem 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 =  2113

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.56977428387172     
GMRes iteration 2, residual = 53.57808898646729     
GMRes iteration 3, residual = 20.029715299932537     
GMRes iteration 4, residual = 7.726387210439835     
GMRes iteration 5, residual = 5.904571084454461     
GMRes iteration 6, residual = 5.22527619349181     
GMRes iteration 7, residual = 4.695885280589728     
GMRes iteration 8, residual = 3.8457422372125007     
GMRes iteration 9, residual = 2.934693441298911     
GMRes iteration 10, residual = 2.1473652961829495     
GMRes iteration 11, residual = 1.7984822439562096     
GMRes iteration 12, residual = 1.6096599099508475     
GMRes iteration 13, residual = 1.4221319492240607     
GMRes iteration 14, residual = 1.2833854974484025     
GMRes iteration 15, residual = 1.1815881666909396     
GMRes iteration 16, residual = 1.0835527239355596     
GMRes iteration 17, residual = 0.9660125018594671     
GMRes iteration 18, residual = 0.8510309673852675     
GMRes iteration 19, residual = 0.7248868740095177     
GMRes iteration 20, residual = 0.5415679170954301     
GMRes iteration 21, residual = 0.36026518980655226     
GMRes iteration 22, residual = 0.2457493021005106     
GMRes iteration 23, residual = 0.16671279818692938     
GMRes iteration 24, residual = 0.11627420287939019     
GMRes iteration 25, residual = 0.08564097546139941     
GMRes iteration 26, residual = 0.07585618311099034     
GMRes iteration 27, residual = 0.0728003061120642     
GMRes iteration 28, residual = 0.07111902237833037     
GMRes iteration 29, residual = 0.06974609122355305     
GMRes iteration 30, residual = 0.06677310051883639     
GMRes iteration 31, residual = 0.0587503815747299     
GMRes iteration 32, residual = 0.045488715765589625     
GMRes iteration 33, residual = 0.03231955789774842     
GMRes iteration 34, residual = 0.023750620126477835     
GMRes iteration 35, residual = 0.020036457915273056     
GMRes iteration 36, residual = 0.018765556446984483     
GMRes iteration 37, residual = 0.018188191274646184     
GMRes iteration 38, residual = 0.017664650298429164     
GMRes iteration 39, residual = 0.017320862322692077     
GMRes iteration 40, residual = 0.01707632550640471     
GMRes iteration 41, residual = 0.016563187037445515     
GMRes iteration 42, residual = 0.016004452831294202     
GMRes iteration 43, residual = 0.01533797033119851     
GMRes iteration 44, residual = 0.013876278672406026     
GMRes iteration 45, residual = 0.012477739890168009     
GMRes iteration 46, residual = 0.011011806260896324     
GMRes iteration 47, residual = 0.00925126073303191     
GMRes iteration 48, residual = 0.00802497855954001     
GMRes iteration 49, residual = 0.007398338859469842     
GMRes iteration 50, residual = 0.0070905670609581915     
GMRes iteration 51, residual = 0.006951365559011698     
GMRes iteration 52, residual = 0.006886446539740405     
GMRes iteration 53, residual = 0.00683514774256927     
GMRes iteration 54, residual = 0.006746108245963398     
GMRes iteration 55, residual = 0.00652848545706547     
GMRes iteration 56, residual = 0.006065521233533223     
GMRes iteration 57, residual = 0.005512292482025432     
GMRes iteration 58, residual = 0.005171873315439872     
GMRes iteration 59, residual = 0.0049687421511336     
GMRes iteration 60, residual = 0.004768523235851341     
GMRes iteration 61, residual = 0.004482629863743435     
GMRes iteration 62, residual = 0.004153366193817344     
GMRes iteration 63, residual = 0.0039036545920747125     
GMRes iteration 64, residual = 0.003731646075198437     
GMRes iteration 65, residual = 0.0035920096499243333     
GMRes iteration 66, residual = 0.003433144943178042     
GMRes iteration 67, residual = 0.0032991632011993917     
GMRes iteration 68, residual = 0.0032037105315469465     
GMRes iteration 69, residual = 0.0031321167925525374     
GMRes iteration 70, residual = 0.0030564091311192903     
GMRes iteration 71, residual = 0.002968547780124301     
GMRes iteration 72, residual = 0.002913940232307292     
GMRes iteration 73, residual = 0.002882554120248412     
GMRes iteration 74, residual = 0.0028390239939250737     
GMRes iteration 75, residual = 0.002773799962716115     
GMRes iteration 76, residual = 0.002724327798535985     
GMRes iteration 77, residual = 0.0026742706307696135     
GMRes iteration 78, residual = 0.0025934001721231548     
GMRes iteration 79, residual = 0.002502067561477762     
GMRes iteration 80, residual = 0.0024198537925243     
GMRes iteration 81, residual = 0.002325612868888091     
GMRes iteration 82, residual = 0.0021279539514522322     
GMRes iteration 83, residual = 0.001849619527761744     
GMRes iteration 84, residual = 0.0016490443387415716     
GMRes iteration 85, residual = 0.001465124118212395     
GMRes iteration 86, residual = 0.0013769601603725874     
GMRes iteration 87, residual = 0.0013527776251821297     
GMRes iteration 88, residual = 0.0013342752661672228     
GMRes iteration 89, residual = 0.0013207135122157902     
GMRes iteration 90, residual = 0.0013101387990788324     
GMRes iteration 91, residual = 0.001296479548141446     
GMRes iteration 92, residual = 0.0012800685821744847     
GMRes iteration 93, residual = 0.001255571611882244     
GMRes iteration 94, residual = 0.0012050893547932283     
GMRes iteration 95, residual = 0.0011290631774763665     
GMRes iteration 96, residual = 0.0010136743889400134     
GMRes iteration 97, residual = 0.0008735506743604157     
GMRes iteration 98, residual = 0.0007328742896240542     
GMRes iteration 99, residual = 0.0006221024427594619     
GMRes iteration 100, residual = 0.0005730412724937148     
GMRes iteration 101, residual = 0.0005531219034355668     
GMRes iteration 102, residual = 0.0005448435654107516     
GMRes iteration 103, residual = 0.0005413806361581299     
GMRes iteration 104, residual = 0.0005384367090824826     
GMRes iteration 105, residual = 0.0005335287513632513     
GMRes iteration 106, residual = 0.0005207662185512949     
GMRes iteration 107, residual = 0.0004931003788887135     
GMRes iteration 108, residual = 0.00043129255081938604     
GMRes iteration 109, residual = 0.00035859605989232995     
GMRes iteration 110, residual = 0.0003124861255384106     
GMRes iteration 111, residual = 0.00029677393118218893     
GMRes iteration 112, residual = 0.00029057353822058846     
GMRes iteration 113, residual = 0.00028677160524983024     
GMRes iteration 114, residual = 0.0002812584248492545     
GMRes iteration 115, residual = 0.0002699616697065845     
GMRes iteration 116, residual = 0.00025549688010131254     
GMRes iteration 117, residual = 0.0002412283469985304     
GMRes iteration 118, residual = 0.00023166813931683693     
GMRes iteration 119, residual = 0.00022576100190271173     
GMRes iteration 120, residual = 0.00022107896617193335     
GMRes iteration 121, residual = 0.00021803817056025317     
GMRes iteration 122, residual = 0.00021636170364122204     
GMRes iteration 123, residual = 0.00021517045893736745     
GMRes iteration 124, residual = 0.0002128772341341456     
GMRes iteration 125, residual = 0.0002071670865210676     
GMRes iteration 126, residual = 0.00019270421406891852     
GMRes iteration 127, residual = 0.00017473174289449076     
GMRes iteration 128, residual = 0.0001586104403044319     
GMRes iteration 129, residual = 0.0001508001632319099     
GMRes iteration 130, residual = 0.00014622184538699625     
GMRes iteration 131, residual = 0.0001433811545196677     
GMRes iteration 132, residual = 0.00014178309224743447     
GMRes iteration 133, residual = 0.000139400019172219     
GMRes iteration 134, residual = 0.00013321179606252717     
GMRes iteration 135, residual = 0.00012059801866750472     
GMRes iteration 136, residual = 0.00010686892638430285     
GMRes iteration 137, residual = 0.00010017975479879376     
GMRes iteration 138, residual = 9.730692558465834e-05     
GMRes iteration 139, residual = 9.507084509309706e-05     
GMRes iteration 140, residual = 9.36049297463807e-05     
GMRes iteration 141, residual = 9.267813122233455e-05     
GMRes iteration 142, residual = 9.205029970279062e-05     
GMRes iteration 143, residual = 9.148690261976094e-05     
GMRes iteration 144, residual = 9.08388254086989e-05     
GMRes iteration 145, residual = 8.992648426806402e-05     
GMRes iteration 146, residual = 8.772751651339395e-05     
GMRes iteration 147, residual = 8.194963464002108e-05     
GMRes iteration 148, residual = 7.259970096288794e-05     
GMRes iteration 149, residual = 6.380993287231848e-05     
GMRes iteration 150, residual = 5.6212194916631565e-05     
GMRes iteration 151, residual = 5.129740107521732e-05     
GMRes iteration 152, residual = 4.910219476221496e-05     
GMRes iteration 153, residual = 4.84571020215911e-05     
GMRes iteration 154, residual = 4.8230897330418064e-05     
GMRes iteration 155, residual = 4.805289160622076e-05     
GMRes iteration 156, residual = 4.782487352728753e-05     
GMRes iteration 157, residual = 4.735477623699916e-05     
GMRes iteration 158, residual = 4.6110977273809055e-05     
GMRes iteration 159, residual = 4.300893531206663e-05     
GMRes iteration 160, residual = 4.056265892475008e-05     
GMRes iteration 161, residual = 3.962909896097159e-05     
GMRes iteration 162, residual = 3.916564678495691e-05     
GMRes iteration 163, residual = 3.8887160909279616e-05     
GMRes iteration 164, residual = 3.8568387131419865e-05     
GMRes iteration 165, residual = 3.749562493456989e-05     
GMRes iteration 166, residual = 3.450460875685047e-05     
GMRes iteration 167, residual = 3.1483739020321166e-05     
GMRes iteration 168, residual = 2.9774850111195844e-05     
GMRes iteration 169, residual = 2.8781280544642717e-05     
GMRes iteration 170, residual = 2.8411183921751597e-05     
GMRes iteration 171, residual = 2.8095615941411662e-05     
GMRes iteration 172, residual = 2.790574817121746e-05     
GMRes iteration 173, residual = 2.7817196373833784e-05     
GMRes iteration 174, residual = 2.772476785911531e-05     
GMRes iteration 175, residual = 2.7522790432061758e-05     
GMRes iteration 176, residual = 2.6905074306252422e-05     
GMRes iteration 177, residual = 2.5395473029801825e-05     
GMRes iteration 178, residual = 2.3820225726407684e-05     
GMRes iteration 179, residual = 2.3151593387954854e-05     
GMRes iteration 180, residual = 2.28864774473411e-05     
GMRes iteration 181, residual = 2.2766872794797768e-05     
GMRes iteration 182, residual = 2.2628020710693544e-05     
GMRes iteration 183, residual = 2.239280300899226e-05     
GMRes iteration 184, residual = 2.1940279014106077e-05     
GMRes iteration 185, residual = 2.0867344793737874e-05     
GMRes iteration 186, residual = 1.8593058046358688e-05     
GMRes iteration 187, residual = 1.543131824528233e-05     
GMRes iteration 188, residual = 1.3608477365684173e-05     
GMRes iteration 189, residual = 1.303076392652798e-05     
GMRes iteration 190, residual = 1.2892132798928075e-05     
GMRes iteration 191, residual = 1.2833211134303372e-05     
GMRes iteration 192, residual = 1.2768960092079417e-05     
GMRes iteration 193, residual = 1.2660414808289072e-05     
GMRes iteration 194, residual = 1.2443086876511675e-05     
GMRes iteration 195, residual = 1.20205339878075e-05     
GMRes iteration 196, residual = 1.1273119252975977e-05     
GMRes iteration 197, residual = 1.0287353963170665e-05     
GMRes iteration 198, residual = 9.487680089648913e-06     
GMRes iteration 199, residual = 9.121838307237268e-06     
GMRes iteration 200, residual = 8.963639095814258e-06     
GMRes iteration 201, residual = 8.876033846681242e-06     
GMRes iteration 202, residual = 8.812699134937597e-06     
GMRes iteration 203, residual = 8.712832627860089e-06     
GMRes iteration 204, residual = 8.39691122492633e-06     
GMRes iteration 205, residual = 7.2131778381012095e-06     
GMRes iteration 206, residual = 5.5862934465709565e-06     
GMRes iteration 207, residual = 4.749618915328175e-06     
GMRes iteration 208, residual = 4.398703197406332e-06     
GMRes iteration 209, residual = 4.301314556793863e-06     
GMRes iteration 210, residual = 4.25875461312359e-06     
GMRes iteration 211, residual = 4.238099623085154e-06     
GMRes iteration 212, residual = 4.226490936556634e-06     
GMRes iteration 213, residual = 4.203000645541212e-06     
GMRes iteration 214, residual = 4.123139908566528e-06     
GMRes iteration 215, residual = 3.952928327591466e-06     
GMRes iteration 216, residual = 3.636817402010788e-06     
GMRes iteration 217, residual = 3.2137155667921047e-06     
GMRes iteration 218, residual = 2.8258103480416604e-06     
GMRes iteration 219, residual = 2.596171069343188e-06     
GMRes iteration 220, residual = 2.500611361403598e-06     
GMRes iteration 221, residual = 2.448708684120466e-06     
GMRes iteration 222, residual = 2.4248587242763213e-06     
GMRes iteration 223, residual = 2.4047002433696094e-06     
GMRes iteration 224, residual = 2.361413238657766e-06     
GMRes iteration 225, residual = 2.2085520954889878e-06     
GMRes iteration 226, residual = 1.926184301980817e-06     
GMRes iteration 227, residual = 1.75052554391299e-06     
GMRes iteration 228, residual = 1.6863026989764508e-06     
GMRes iteration 229, residual = 1.6595696301937407e-06     
GMRes iteration 230, residual = 1.6488654182924861e-06     
GMRes iteration 231, residual = 1.642480214415018e-06     
GMRes iteration 232, residual = 1.6342475245878175e-06     
GMRes iteration 233, residual = 1.6189638817357479e-06     
GMRes iteration 234, residual = 1.5907892696060879e-06     
GMRes iteration 235, residual = 1.5057072232610156e-06     
GMRes iteration 236, residual = 1.3217379405755436e-06     
GMRes iteration 237, residual = 1.2154174107277621e-06     
GMRes iteration 238, residual = 1.1885513196251603e-06     
GMRes iteration 239, residual = 1.1790368383772774e-06     
GMRes iteration 240, residual = 1.1725314246035392e-06     
GMRes iteration 241, residual = 1.1618976499110567e-06     
GMRes iteration 242, residual = 1.1259289083955135e-06     
GMRes iteration 243, residual = 1.0071739829128448e-06     
GMRes iteration 244, residual = 8.85122147336326e-07     
GMRes iteration 245, residual = 8.0444503853742e-07     
GMRes iteration 246, residual = 7.454426644784821e-07     
GMRes iteration 247, residual = 7.291711541100854e-07     
GMRes iteration 248, residual = 7.249716113147211e-07     
GMRes iteration 249, residual = 7.226396350762251e-07     
GMRes iteration 250, residual = 7.2017303441046e-07     
GMRes iteration 251, residual = 7.116041196929695e-07     
GMRes iteration 252, residual = 6.751986887056901e-07     
GMRes iteration 253, residual = 6.103247393548231e-07     
GMRes iteration 254, residual = 5.6474083812216e-07     
GMRes iteration 255, residual = 5.513634891150593e-07     
GMRes iteration 256, residual = 5.465598494977997e-07     
GMRes iteration 257, residual = 5.432479363205316e-07     
GMRes iteration 258, residual = 5.405172261618714e-07     
GMRes iteration 259, residual = 5.369780135990282e-07     
GMRes iteration 260, residual = 5.281700100538163e-07     
GMRes iteration 261, residual = 5.010582431830561e-07     
GMRes iteration 262, residual = 4.3916183590447767e-07     
GMRes iteration 263, residual = 3.7399834777228846e-07     
GMRes iteration 264, residual = 3.3581823245522816e-07     
GMRes iteration 265, residual = 3.1925245445517386e-07     
GMRes iteration 266, residual = 3.1313457868871474e-07     
GMRes iteration 267, residual = 3.109716208139693e-07     
GMRes iteration 268, residual = 3.0924454646329627e-07     
GMRes iteration 269, residual = 3.0521397985025993e-07     
GMRes iteration 270, residual = 2.9132323569237196e-07     
GMRes iteration 271, residual = 2.6982036969349137e-07     
GMRes iteration 272, residual = 2.590323992237061e-07     
GMRes iteration 273, residual = 2.5530488655551137e-07     
GMRes iteration 274, residual = 2.5383210278592465e-07     
GMRes iteration 275, residual = 2.528564302307273e-07     
GMRes iteration 276, residual = 2.5107350673152433e-07     
GMRes iteration 277, residual = 2.4572997653407723e-07     
GMRes iteration 278, residual = 2.2607053752537137e-07     
GMRes iteration 279, residual = 1.8632594967792607e-07     
GMRes iteration 280, residual = 1.6189432390714922e-07     
GMRes iteration 281, residual = 1.5404247303759677e-07     
GMRes iteration 282, residual = 1.5130987762735458e-07     
GMRes iteration 283, residual = 1.498111327989046e-07     
GMRes iteration 284, residual = 1.4873696972225248e-07     
GMRes iteration 285, residual = 1.4671880088782277e-07     
GMRes iteration 286, residual = 1.4008973675192804e-07     
GMRes iteration 287, residual = 1.2770249370260618e-07     
GMRes iteration 288, residual = 1.0702141400520125e-07     
GMRes iteration 289, residual = 9.175872738873346e-08     
GMRes iteration 290, residual = 8.647687877449822e-08     
GMRes iteration 291, residual = 8.49389081102102e-08     
GMRes iteration 292, residual = 8.463404926009913e-08     
GMRes iteration 293, residual = 8.440326445962579e-08     
GMRes iteration 294, residual = 8.368825549515328e-08     
GMRes iteration 295, residual = 8.13424824598137e-08     
GMRes iteration 296, residual = 7.499681896350918e-08     
GMRes iteration 297, residual = 6.571969807877304e-08     
GMRes iteration 298, residual = 5.961345335116281e-08     
GMRes iteration 299, residual = 5.7393597320880815e-08     
GMRes iteration 300, residual = 5.677529669287566e-08     
GMRes iteration 301, residual = 5.6373895416719447e-08     
GMRes iteration 302, residual = 5.589758038918744e-08     
GMRes iteration 303, residual = 5.481047305669175e-08     
GMRes iteration 304, residual = 5.06237419913053e-08     
GMRes iteration 305, residual = 4.070032536520796e-08     
GMRes iteration 306, residual = 3.4866796957271934e-08     
GMRes iteration 307, residual = 3.3290609034231475e-08     
GMRes iteration 308, residual = 3.2915092879955606e-08     
GMRes iteration 309, residual = 3.2736428039139494e-08     
GMRes iteration 310, residual = 3.258679443589558e-08     
GMRes iteration 311, residual = 3.2289193076434484e-08     
GMRes iteration 312, residual = 3.154634444348647e-08     
GMRes iteration 313, residual = 2.98424402121489e-08     
GMRes iteration 314, residual = 2.628836058068018e-08     
GMRes iteration 315, residual = 2.3372534917344726e-08     
GMRes iteration 316, residual = 2.248993123797571e-08     
GMRes iteration 317, residual = 2.2265598406854945e-08     
GMRes iteration 318, residual = 2.2161941490886226e-08     
GMRes iteration 319, residual = 2.20982884950682e-08     
GMRes iteration 320, residual = 2.1980859778528047e-08     
GMRes iteration 321, residual = 2.1434256977288073e-08     
GMRes iteration 322, residual = 1.9612451532598593e-08     
GMRes iteration 323, residual = 1.590132720715626e-08     
GMRes iteration 324, residual = 1.3451447960005195e-08     
GMRes iteration 325, residual = 1.2935341448800856e-08     
GMRes iteration 326, residual = 1.2831404731771259e-08     
GMRes iteration 327, residual = 1.2785627749577692e-08     
GMRes iteration 328, residual = 1.2713858638289756e-08     
GMRes iteration 329, residual = 1.2341982790788906e-08     
GMRes iteration 330, residual = 1.0867695326085594e-08     
GMRes iteration 331, residual = 9.512413410997898e-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.0006228098891097808
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.20042751824788e-06
BaseWebGuiScene