HEX20 — 20-node serendipity hexahedron ====================================== The quadratic-hex serendipity kernel. Eight corner nodes plus twelve mid-edge nodes — 20 total, 60 DOFs per element. The serendipity basis omits the body / face / mid-element nodes of the full Lagrange family while preserving exact reproduction of quadratic polynomials (Ergatoudis, Irons & Zienkiewicz 1968). Dispatches collapsed-corner input to :doc:`wedge15_pyr13` automatically. * **Spec:** ``ELEMENTS.HEX20`` Shape functions --------------- Two families on the reference cube :math:`\hat\Omega = [-1, +1]^{3}`: **Corner nodes** :math:`(\xi_i, \eta_i, \zeta_i) \in \{-1, +1\}^{3}`: .. math:: N^{c}_{i}(\xi, \eta, \zeta) = \tfrac{1}{8} (1 + \xi_i \xi)(1 + \eta_i \eta)(1 + \zeta_i \zeta) (\xi_i \xi + \eta_i \eta + \zeta_i \zeta - 2). **Mid-edge nodes** with one of :math:`\xi_j, \eta_j, \zeta_j` zero (the axis along which the edge runs). For :math:`\xi_j = 0` mid-edges: .. math:: N^{m}_{j}(\xi, \eta, \zeta) = \tfrac{1}{4} (1 - \xi^{2})(1 + \eta_j \eta)(1 + \zeta_j \zeta), with analogous expressions for :math:`\eta_j = 0` and :math:`\zeta_j = 0` mid-edges (the missing-variable axis carries the :math:`(1 - \cdot^{2})` factor). The basis satisfies Kronecker-delta interpolation at all 20 nodes and partition-of-unity on :math:`\hat\Omega`. See :ref:`sphx_glr_gallery_elements_solid186_example_hex20_shape_function_slices.py` for slice plots of the eight functions that don't vanish on the :math:`\zeta = -1` face. Integration ----------- * **Stiffness:** 2 × 2 × 2 reduced Gauss-Legendre (8 points) by default — Bedrosian (1992) and the standard practitioner choice for serendipity hex. Bending-mode response stays exact; spurious hourglass modes are absent because the 8 points still couple all 60 stiffness rows. * **Mass:** 14-point **Irons rule** (Irons 1971). Integrates polynomials of degree ≤ 5 exactly on the cube — a non-tensor- product cubature constructed by hand for the cube's symmetry group with the minimum point count for this degree. Preferred over the 27-point :math:`3 \times 3 \times 3` rule because it preserves the consistent-mass sparsity pattern without overshooting the integrand polynomial degree. The full 27-point :math:`3 \times 3 \times 3` rule is also available via ``ELEMENTS.HEX20(integration="full")``; useful only for cross-checking the reduced rule on a sensitive bending mode. Verification cross-references ----------------------------- * :ref:`sphx_glr_gallery_elements_solid186_example_solid186.py` — single-hex uniaxial via mid-edge consistent loads. * :ref:`sphx_glr_gallery_elements_solid186_example_hex20_reference_geometry.py` — corner + mid-edge node layout, stiffness and mass Gauss points. Implementation: :mod:`femorph_solver.elements.hex20`; collapsed-corner forms in :mod:`femorph_solver.elements._wedge15_pyr13`. References ---------- * Ergatoudis, J. G., Irons, B. M. and Zienkiewicz, O. C. (1968) "Curved isoparametric quadrilateral elements for finite element analysis," *International Journal of Solids and Structures* 4 (1), 31–42. * Hughes, T. J. R. (2000) *The Finite Element Method — Linear Static and Dynamic Finite Element Analysis*, Dover, §3.7 (serendipity family). * Zienkiewicz, O. C. and Taylor, R. L. (2013) *The Finite Element Method: Its Basis and Fundamentals*, 7th ed., Butterworth-Heinemann, §8.7.3. * Irons, B. M. (1971) "Quadrature rules for brick based finite elements," *International Journal for Numerical Methods in Engineering* 3 (2), 293–294 (14-point degree-5 rule for the consistent mass integral). * Bedrosian, G. (1992) "Shape functions and integration formulas for three-dimensional finite element analysis," *International Journal for Numerical Methods in Engineering* 35 (1), 95–108 (reduced 2 × 2 × 2 stiffness). * Cook, R. D., Malkus, D. S., Plesha, M. E., Witt, R. J. (2002) *Concepts and Applications of Finite Element Analysis*, 4th ed., Wiley, §6.5 (quadratic hex).