{ "cells": [ { "cell_type": "markdown", "id": "a1", "metadata": {}, "source": [ "# Element Boundaries\n", "\n", "Visualize CoefficientFunctions evaluated on element boundaries (facets).\n", "Each element face is rendered independently, so inter-element jumps become visible —\n", "useful for DG/HDG methods, L2 spaces, and FacetFESpace functions.\n", "\n", "- **2D**: `FacetFunctionData` + `FacetCFRenderer` draw colored thick lines on element edges.\n", "- **3D**: `FacetCFRenderer3D` renders all tet faces with slight shrink, using per-face Bernstein evaluation." ] }, { "cell_type": "markdown", "id": "a2", "metadata": {}, "source": [ "## 2D: Continuous Function\n", "\n", "Even for a smooth CF, element boundary rendering shows the value along every edge." ] }, { "cell_type": "code", "execution_count": null, "id": "a3", "metadata": {}, "outputs": [], "source": [ "from ngsolve import *\n", "from ngsolve_webgpu import *\n", "from webgpu.jupyter import Draw\n", "\n", "mesh = Mesh(unit_square.GenerateMesh(maxh=0.15))\n", "mdata = MeshData(mesh)\n", "\n", "facet_data = FacetFunctionData(mdata, sin(5 * x) * y, order=3)\n", "colormap = Colormap(minval=-1, maxval=1)\n", "renderer = FacetCFRenderer(facet_data, colormap=colormap)\n", "wireframe = MeshWireframe2d(mdata)\n", "\n", "Draw([renderer, wireframe, Colorbar(colormap)])" ] }, { "cell_type": "markdown", "id": "a4", "metadata": {}, "source": [ "## 2D: L2 Function — Inter-Element Jumps\n", "\n", "For discontinuous (L2) functions, each edge is rendered twice (once per adjacent element),\n", "naturally revealing jumps." ] }, { "cell_type": "code", "execution_count": null, "id": "a5", "metadata": {}, "outputs": [], "source": [ "mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))\n", "mdata = MeshData(mesh)\n", "\n", "fes = L2(mesh, order=0)\n", "gf = GridFunction(fes)\n", "gf.vec.FV().NumPy()[:] = [i/len(gf.vec) for i in range(len(gf.vec))]\n", "\n", "facet_data = FacetFunctionData(mdata, gf, order=2)\n", "colormap = Colormap()\n", "renderer = FacetCFRenderer(facet_data, colormap=colormap)\n", "wireframe = MeshWireframe2d(mdata)\n", "\n", "Draw([renderer, wireframe, Colorbar(colormap)])" ] }, { "cell_type": "markdown", "id": "a6", "metadata": {}, "source": [ "## 3D: Continuous Function on a Tet Mesh\n", "\n", "`FacetCFRenderer3D` renders all faces of all 3D elements with a slight shrink\n", "so adjacent element faces separate, revealing the function on each face.\n", "Clipping cuts through individual faces (not whole elements)." ] }, { "cell_type": "code", "execution_count": null, "id": "a7", "metadata": {}, "outputs": [], "source": [ "from ngsolve_webgpu.facet_cf import FacetCFRenderer3D\n", "from webgpu.clipping import Clipping\n", "\n", "mesh = Mesh(unit_cube.GenerateMesh(maxh=0.3))\n", "mdata = MeshData(mesh)\n", "\n", "colormap = Colormap()\n", "clipping = Clipping()\n", "clipping.mode = clipping.Mode.PLANE\n", "clipping.center = [0.5, 0.5, 0.5]\n", "clipping.normal = [1, 0, 0]\n", "\n", "renderer = FacetCFRenderer3D(mdata, x + 2 * y, order=2,\n", " colormap=colormap, clipping=clipping)\n", "\n", "Draw([renderer, Colorbar(colormap)])" ] }, { "cell_type": "markdown", "id": "a8", "metadata": {}, "source": [ "## 3D: FacetFESpace — HDG Example\n", "\n", "FacetFESpace functions live on element interfaces and cannot be evaluated inside elements.\n", "Element boundary rendering handles this naturally via `MapToAllElements` with `element_boundary=BND`." ] }, { "cell_type": "code", "execution_count": null, "id": "a9", "metadata": {}, "outputs": [], "source": [ "mesh = Mesh(unit_cube.GenerateMesh(maxh=0.3))\n", "mdata = MeshData(mesh)\n", "\n", "order = 2\n", "V = L2(mesh, order=order)\n", "F = FacetFESpace(mesh, order=order, dirichlet=\".*\")\n", "fes = V * F\n", "\n", "u, uhat = fes.TrialFunction()\n", "v, vhat = fes.TestFunction()\n", "n = specialcf.normal(3)\n", "h = specialcf.mesh_size\n", "alpha = 4 * order ** 2 / h\n", "\n", "a = BilinearForm(fes, condense=True)\n", "a += grad(u) * grad(v) * dx\n", "a += (-n * grad(u) * (v - vhat) - n * grad(v) * (u - uhat)\n", " + alpha * (u - uhat) * (v - vhat)) * dx(element_boundary=True)\n", "\n", "f = LinearForm(fes)\n", "f += 1 * v * dx\n", "\n", "gf = GridFunction(fes)\n", "with TaskManager():\n", " a.Assemble()\n", " f.Assemble()\n", " solvers.BVP(bf=a, lf=f, gf=gf)\n", "\n", "# Render the facet component (uhat)\n", "colormap = Colormap()\n", "clipping = Clipping()\n", "clipping.mode = clipping.Mode.PLANE\n", "clipping.center = [0.5, 0.5, 0.5]\n", "clipping.normal = [1, 0, 0]\n", "\n", "renderer = FacetCFRenderer3D(mdata, gf.components[1], order=order,\n", " colormap=colormap, clipping=clipping)\n", "\n", "Draw([renderer, Colorbar(colormap)])" ] }, { "cell_type": "code", "execution_count": null, "id": "93442706-d076-492e-8c08-7ca04754b57a", "metadata": {}, "outputs": [], "source": [] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.14.5" } }, "nbformat": 4, "nbformat_minor": 5 }