Element Boundaries#

Visualize CoefficientFunctions evaluated on element boundaries (facets). Each element face is rendered independently, so inter-element jumps become visible — useful for DG/HDG methods, L2 spaces, and FacetFESpace functions.

  • 2D: FacetFunctionData + FacetCFRenderer draw colored thick lines on element edges.

  • 3D: FacetCFRenderer3D renders all tet faces with slight shrink, using per-face Bernstein evaluation.

2D: Continuous Function#

Even for a smooth CF, element boundary rendering shows the value along every edge.

[1]:
from ngsolve import *
from ngsolve_webgpu import *
from webgpu.jupyter import Draw

mesh = Mesh(unit_square.GenerateMesh(maxh=0.15))
mdata = MeshData(mesh)

facet_data = FacetFunctionData(mdata, sin(5 * x) * y, order=3)
colormap = Colormap(minval=-1, maxval=1)
renderer = FacetCFRenderer(facet_data, colormap=colormap)
wireframe = MeshWireframe2d(mdata)

Draw([renderer, wireframe, Colorbar(colormap)])
[1]:

2D: L2 Function — Inter-Element Jumps#

For discontinuous (L2) functions, each edge is rendered twice (once per adjacent element), naturally revealing jumps.

[2]:
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
mdata = MeshData(mesh)

fes = L2(mesh, order=0)
gf = GridFunction(fes)
gf.vec.FV().NumPy()[:] = [i/len(gf.vec) for i in range(len(gf.vec))]

facet_data = FacetFunctionData(mdata, gf, order=2)
colormap = Colormap()
renderer = FacetCFRenderer(facet_data, colormap=colormap)
wireframe = MeshWireframe2d(mdata)

Draw([renderer, wireframe, Colorbar(colormap)])
[2]:

3D: Continuous Function on a Tet Mesh#

FacetCFRenderer3D renders all faces of all 3D elements with a slight shrink so adjacent element faces separate, revealing the function on each face. Clipping cuts through individual faces (not whole elements).

[3]:
from ngsolve_webgpu.facet_cf import FacetCFRenderer3D
from webgpu.clipping import Clipping

mesh = Mesh(unit_cube.GenerateMesh(maxh=0.3))
mdata = MeshData(mesh)

colormap = Colormap()
clipping = Clipping()
clipping.mode = clipping.Mode.PLANE
clipping.center = [0.5, 0.5, 0.5]
clipping.normal = [1, 0, 0]

renderer = FacetCFRenderer3D(mdata, x + 2 * y, order=2,
                              colormap=colormap, clipping=clipping)

Draw([renderer, Colorbar(colormap)])
[3]:

3D: FacetFESpace — HDG Example#

FacetFESpace functions live on element interfaces and cannot be evaluated inside elements. Element boundary rendering handles this naturally via MapToAllElements with element_boundary=BND.

[4]:
mesh = Mesh(unit_cube.GenerateMesh(maxh=0.3))
mdata = MeshData(mesh)

order = 2
V = L2(mesh, order=order)
F = FacetFESpace(mesh, order=order, dirichlet=".*")
fes = V * F

u, uhat = fes.TrialFunction()
v, vhat = fes.TestFunction()
n = specialcf.normal(3)
h = specialcf.mesh_size
alpha = 4 * order ** 2 / h

a = BilinearForm(fes, condense=True)
a += grad(u) * grad(v) * dx
a += (-n * grad(u) * (v - vhat) - n * grad(v) * (u - uhat)
      + alpha * (u - uhat) * (v - vhat)) * dx(element_boundary=True)

f = LinearForm(fes)
f += 1 * v * dx

gf = GridFunction(fes)
with TaskManager():
    a.Assemble()
    f.Assemble()
    solvers.BVP(bf=a, lf=f, gf=gf)

# Render the facet component (uhat)
colormap = Colormap()
clipping = Clipping()
clipping.mode = clipping.Mode.PLANE
clipping.center = [0.5, 0.5, 0.5]
clipping.normal = [1, 0, 0]

renderer = FacetCFRenderer3D(mdata, gf.components[1], order=order,
                              colormap=colormap, clipping=clipping)

Draw([renderer, Colorbar(colormap)])
[4]:

[ ]: