Element kernels#
femorph-solver ships a structural element library. Each kernel is
independently implemented from the open literature; every
non-trivial algorithm (shape function, integration rule,
formulation variant) carries a References block in its source
module. The pre-processing pages under
Element library are the
user-facing how-to; this section is the engineering reference —
topology, shape functions, integration rules, and cites for every
element.
HEX8#
8-node trilinear hexahedron. 24 DOFs / element.
Shape functions. Trilinear Lagrange on the reference cube, \(N_i = \tfrac{1}{8}(1+\xi_i\xi)(1+\eta_i\eta)(1+\zeta_i\zeta)\) — Zienkiewicz & Taylor §6; Cook §6.2.
Integration. 2×2×2 Gauss-Legendre (product rule). K uses Hughes’ volumetric B-bar (default); Simo–Rifai 9-parameter enhanced-assumed-strain variant available via the
_SOLID185_FORMULATION = "enhanced_strain"flag.Mass. Consistent at 2×2×2; row-sum (HRZ) lumping optional.
MAPDL alias.
SOLID185. KEYOPT(2)=0 / 1 / 3 supported.
Full cites in femorph_solver.elements.solid185.
HEX20#
20-node quadratic serendipity hexahedron. 60 DOFs / element.
Dispatches collapsed-corner input to the WEDGE15 / PYR13
kernels automatically.
Shape functions. Corner-plus-mid-edge serendipity basis from Ergatoudis, Irons, Zienkiewicz (IJSS 1968); Zienkiewicz & Taylor §8.7.3.
Stiffness integration. 2×2×2 uniform reduced by default (Bedrosian 1992) — matches MAPDL KEYOPT(2)=0. Full 3×3×3 rule (27 points) available via
_SOLID186_FORMULATION = "full".Mass integration. 14-point Irons rule (Irons, IJNME 1971) — integrates polynomials through degree 5 exactly and matches the MAPDL KEYOPT(2)=0 mass kernel bit-exactly.
MAPDL alias.
SOLID186(plusSOLID186W/SOLID186Pfor collapsed-corner forms).
Full cites in femorph_solver.elements.solid186 and
femorph_solver.elements.solid186_degenerate.
TET10#
10-node quadratic tetrahedron. 30 DOFs / element.
Shape functions. Lagrange in volume coordinates,
Nᵢ = Lᵢ(2Lᵢ − 1)on corners,N_{ij} = 4LᵢLⱼon mid-edges — Zienkiewicz & Taylor §8.8.2; Cook Table 6.5-1.Integration. 4-point Keast Gauss rule (Keast, CMAME 1986, degree-2 exact, weights
1/24) — used for both K and M.MAPDL alias.
SOLID187.
Full cites in femorph_solver.elements.solid187.
WEDGE15 / PYR13#
15-node quadratic wedge and 13-node apex-singular pyramid — degenerate serendipity-hex forms.
WEDGE15 shape functions. Serendipity on the triangular base × quadratic through thickness (Bedrosian 1992; ZT §8.8.3).
WEDGE15 integration. 9-point (triangle × line) for K; 21-point for consistent mass. Triangle rules from Dunavant 1985.
PYR13 shape functions. Apex-singular Peterson / Bedrosian 1992 polynomial basis with a removable
1/(1−ζ)factor. Alternatives cited: Zgainski-Coulomb-Marechal (IEEE Mag. 1996), Wachspress, A Rational Finite Element Basis, 1975.PYR13 integration. 2×2×2 Duffy collapsed-hex rule (Duffy, SIAM J. Numer. Anal. 1982) — matches MAPDL’s SOLID186-pyramid mass kernel bit-exactly on
ptr.cdb.MAPDL alias.
SOLID186W/SOLID186P(auto-dispatched from theHEX20entry).
BEAM2#
3-D 2-node Euler–Bernoulli beam. 6 DOFs / node, 12 total.
Kinematics. Slender-beam limit of Timoshenko-beam theory (Φ = 12EI/(κAGL²) → 0); shear neglected. Timoshenko-beam form is roadmap.
Shape functions. Hermite cubic for transverse displacement + slope; linear for axial and torsion. ZT §2.5.1; Cook §2.4–§2.6.
Mass. Consistent from Hermite cubics (Cook Table 16.3-1); row-sum lumped available.
MAPDL alias.
BEAM188.
Full cites in femorph_solver.elements.beam188.
QUAD4_SHELL#
4-node flat Mindlin–Reissner shell. 6 DOFs / node, 24 total.
Kinematics. First-order shear-deformation theory (Mindlin JAM 1951; Reissner JAM 1945) with independent bending rotation.
Integration. 2×2 Gauss for membrane + bending; 1×1 selective-reduced Gauss on transverse shear (Malkus & Hughes CMAME 1978) to suppress shear locking on thin shells.
Drilling DOF stabilisation.
α·G·t·Aper-node penalty (Allman 1984; Hughes-Brezzi CMAME 1989) keeps the local 24×24 non-singular whenθ_zis free.MAPDL alias.
SHELL181.
Full cites in femorph_solver.elements.shell181. MITC4
variant (Bathe-Dvorkin IJNME 1985) is roadmap.
QUAD4_PLANE#
2-D 4-node bilinear quadrilateral, plane stress / plane strain.
Shape functions. Bilinear Q4 Lagrange (ZT §6.3.2).
Integration. 2×2 Gauss-Legendre.
Constitutive. Plane-stress / plane-strain selected by
_PLANE_MODEmaterial flag; Cook §3.5–§3.6.MAPDL alias.
PLANE182.
TRUSS2#
2-node 3-D axial bar. 3 translational DOFs / node.
Analytic stiffness + consistent mass (Cook §2.3 / Table 16.3-1).
MAPDL alias: LINK180.
SPRING#
2-node 3-D longitudinal spring. 3 translational DOFs / node.
Direction-cosine-rotated 6×6 stiffness block. No mass
contribution. MAPDL alias: COMBIN14.
POINT_MASS#
Single-node lumped translational mass. Diagonal 3×3 contribution
to global M. MAPDL alias: MASS21 (KEYOPT(3)=2 variant).
Full cites in femorph_solver.elements.mass21.
See also
Element library — per-element KEYOPT tables, material-property slots, and runnable demos.