.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "gallery/elements/solid185/example_hex8_shape_function_surfaces.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code. .. rst-class:: sphx-glr-example-title .. _sphx_glr_gallery_elements_solid185_example_hex8_shape_function_surfaces.py: HEX8 shape-function surfaces on the reference cube ================================================== Renders all eight trilinear Lagrange shape functions :math:`N_i(\xi, \eta, \zeta)` for the 8-node hexahedron on the reference cube :math:`\hat\Omega = [-1, +1]^3`. The basis functions are .. math:: N_i(\xi, \eta, \zeta) = \tfrac{1}{8} (1 + \xi_i \xi)(1 + \eta_i \eta)(1 + \zeta_i \zeta), \qquad i = 1, \ldots, 8, with :math:`(\xi_i, \eta_i, \zeta_i) \in \{-1, +1\}^3` the nodal coordinates. Each :math:`N_i` is a tensor product of 1-D linear "hat" functions, equal to 1 at corner :math:`i` and 0 at every other corner — the Kronecker-delta interpolation property. The plot below shows :math:`N_i` evaluated on the :math:`\zeta = -1` cube face (:math:`i = 1, \ldots, 4`) and the :math:`\zeta = +1` face (:math:`i = 5, \ldots, 8`). References ---------- * Hughes, T. J. R. (2000) *The Finite Element Method — Linear Static and Dynamic Finite Element Analysis*, Dover, §3.1 — trilinear shape functions. * Cook, R. D., Malkus, D. S., Plesha, M. E., Witt, R. J. (2002) *Concepts and Applications of Finite Element Analysis*, 4th ed., Wiley, §6.2. Implementation: :class:`femorph_solver.elements.hex8.Hex8`. .. GENERATED FROM PYTHON SOURCE LINES 38-44 .. code-block:: Python from __future__ import annotations import matplotlib.pyplot as plt import numpy as np .. GENERATED FROM PYTHON SOURCE LINES 45-46 Corner natural-coordinate map (matches VTK_HEXAHEDRON ordering) .. GENERATED FROM PYTHON SOURCE LINES 46-67 .. code-block:: Python corners = np.array( [ [-1, -1, -1], [+1, -1, -1], [+1, +1, -1], [-1, +1, -1], [-1, -1, +1], [+1, -1, +1], [+1, +1, +1], [-1, +1, +1], ], dtype=float, ) def n_i(xi, eta, zeta, i): """Trilinear Lagrange shape function for corner i.""" xi_i, eta_i, zeta_i = corners[i] return 0.125 * (1 + xi_i * xi) * (1 + eta_i * eta) * (1 + zeta_i * zeta) .. GENERATED FROM PYTHON SOURCE LINES 68-69 Sample on a 41×41 grid on each ζ = ±1 face .. GENERATED FROM PYTHON SOURCE LINES 69-74 .. code-block:: Python n = 41 xi = eta = np.linspace(-1.0, 1.0, n) XI, ETA = np.meshgrid(xi, eta, indexing="ij") .. GENERATED FROM PYTHON SOURCE LINES 75-76 4×2 surface plot — one column per ζ-face, 4 rows for the 4 corners .. GENERATED FROM PYTHON SOURCE LINES 76-108 .. code-block:: Python fig, axes = plt.subplots( 4, 2, figsize=(10.0, 14.0), subplot_kw={"projection": "3d"}, constrained_layout=True, ) fig.suptitle("HEX8 shape functions $N_i(\\xi, \\eta, \\zeta)$ on cube faces", fontsize=13) face_zeta = (-1.0, +1.0) node_indices_per_face = ((0, 1, 2, 3), (4, 5, 6, 7)) for col, (zeta, idxs) in enumerate(zip(face_zeta, node_indices_per_face, strict=True)): for row, idx in enumerate(idxs): ax = axes[row, col] N = n_i(XI, ETA, zeta, idx) ax.plot_surface(XI, ETA, N, cmap="plasma", linewidth=0, antialiased=True) x_i, y_i, _ = corners[idx] ax.set_title( f"$N_{{{idx + 1}}}$ at corner $({int(x_i):+d}, {int(y_i):+d}, {int(zeta):+d})$", fontsize=10, ) ax.set_xlabel(r"$\xi$") ax.set_ylabel(r"$\eta$") ax.set_zlabel(r"$N$") ax.set_xticks([-1, 0, 1]) ax.set_yticks([-1, 0, 1]) ax.set_zticks([0, 0.5, 1]) ax.view_init(elev=22, azim=-55) fig.show() .. image-sg:: /gallery/elements/solid185/images/sphx_glr_example_hex8_shape_function_surfaces_001.png :alt: HEX8 shape functions $N_i(\xi, \eta, \zeta)$ on cube faces, $N_{1}$ at corner $(-1, -1, -1)$, $N_{5}$ at corner $(-1, -1, +1)$, $N_{2}$ at corner $(+1, -1, -1)$, $N_{6}$ at corner $(+1, -1, +1)$, $N_{3}$ at corner $(+1, +1, -1)$, $N_{7}$ at corner $(+1, +1, +1)$, $N_{4}$ at corner $(-1, +1, -1)$, $N_{8}$ at corner $(-1, +1, +1)$ :srcset: /gallery/elements/solid185/images/sphx_glr_example_hex8_shape_function_surfaces_001.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 109-114 Verify partition of unity at every interior sample point :math:`\sum_{i=1}^8 N_i(\xi, \eta, \zeta) = 1` for any :math:`(\xi, \eta, \zeta) \in \hat\Omega`. Sample on a 21×21×21 lattice and assert. .. GENERATED FROM PYTHON SOURCE LINES 114-124 .. code-block:: Python n3 = 21 xi3 = eta3 = zeta3 = np.linspace(-1.0, 1.0, n3) XI3, ETA3, ZETA3 = np.meshgrid(xi3, eta3, zeta3, indexing="ij") total = sum(n_i(XI3, ETA3, ZETA3, idx) for idx in range(8)) max_dev = float(np.max(np.abs(total - 1.0))) assert max_dev < 1.0e-14, f"partition-of-unity violated: max |Σ N_i - 1| = {max_dev:.3e}" print(f"max |Σ N_i - 1| over a 21^3 grid = {max_dev:.3e}") print("OK — eight trilinear HEX8 shape functions form a partition of unity.") .. rst-class:: sphx-glr-script-out .. code-block:: none max |Σ N_i - 1| over a 21^3 grid = 4.441e-16 OK — eight trilinear HEX8 shape functions form a partition of unity. .. GENERATED FROM PYTHON SOURCE LINES 125-131 Verify Kronecker-delta interpolation at every corner :math:`N_i(\boldsymbol{\xi}_j) = \delta_{ij}`: each shape function equals 1 at its own corner and 0 at all other corners — required so the nodal-displacement DOFs are literally the displacement values. .. GENERATED FROM PYTHON SOURCE LINES 131-138 .. code-block:: Python print() for i in range(8): row = [n_i(*corners[j], i) for j in range(8)] expected = [1.0 if j == i else 0.0 for j in range(8)] np.testing.assert_allclose(row, expected, atol=1.0e-15) print("OK — Kronecker-delta interpolation verified at every corner.") .. rst-class:: sphx-glr-script-out .. code-block:: none OK — Kronecker-delta interpolation verified at every corner. .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 1.943 seconds) .. _sphx_glr_download_gallery_elements_solid185_example_hex8_shape_function_surfaces.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: example_hex8_shape_function_surfaces.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: example_hex8_shape_function_surfaces.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: example_hex8_shape_function_surfaces.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_