"""Truss2 — 3D 2-node spar (truss).
Linear shape functions in the natural coordinate ``s ∈ [-1, 1]`` for
all three translations. Only axial strain contributes to internal
work::
B = (1/L) [-l, -m, -n, l, m, n] (1 × 6)
K_e = (E·A/L) · [[ C, -C],
[-C, C]] where C_ij = d_i·d_j (3 × 3)
with ``(l, m, n)`` the unit vector along the element axis.
Consistent mass::
M_e = ρ·A·L / 6 · [[2 I₃, I₃],
[ I₃, 2 I₃]]
Lumped mass: ``M_e = ρ·A·L/2 · I₆``.
Real constants
--------------
real[0]: AREA — cross-sectional area (mandatory)
real[1]: ADDMAS — added mass per unit length (optional; unused here)
real[2]: ISTRN — initial strain (optional; unused here — stress path only)
Material properties
-------------------
EX : Young's modulus (mandatory)
DENS : mass density (required for M_e)
References
----------
* Shape functions — linear Lagrange on a 1-D 2-node element
(``N₁ = ½(1 − s)``, ``N₂ = ½(1 + s)`` on ``s ∈ [−1, 1]``):
Zienkiewicz, O.C. and Taylor, R.L. *The Finite Element Method:
Its Basis and Fundamentals*, 7th ed., Butterworth-Heinemann,
2013, §2.3. Cook, Malkus, Plesha, Witt, *Concepts and
Applications of Finite Element Analysis*, 4th ed., Wiley, 2002,
§2.4.
* Axial-only strain (spar / truss element, no bending / shear
contribution, only longitudinal stiffness): Cook et al. §2.3.
* Consistent mass ``ρAL/6 · [[2, 1], [1, 2]]`` on the two axial
DOFs (and its 6×6 extension to 3-D for the three translational
pairs): Cook et al. Table 16.3-1; Zienkiewicz & Taylor §16.2.1.
* Row-sum lumped mass ``ρAL/2`` per translational DOF: Hinton,
Rock & Zienkiewicz, *Earthquake Engng. Struct. Dyn.* 4 (1976),
pp. 245–249.
"""
from __future__ import annotations
import numpy as np
from femorph_solver.elements._base import ElementBase
from femorph_solver.elements._registry import register
def _axial_frame(coords: np.ndarray) -> tuple[float, np.ndarray]:
"""Return element length and unit direction vector (I → J)."""
d = coords[1] - coords[0]
L = float(np.linalg.norm(d))
if L == 0.0:
raise ValueError("TRUSS2: coincident nodes, element length is zero")
return L, d / L
[docs]
@register
class Truss2(ElementBase):
name = "TRUSS2"
n_nodes = 2
n_dof_per_node = 3 # UX, UY, UZ
vtk_cell_type = 3 # VTK_LINE
[docs]
@staticmethod
def ke(
coords: np.ndarray,
material: dict[str, float],
real: np.ndarray,
) -> np.ndarray:
L, d = _axial_frame(np.asarray(coords, dtype=np.float64))
E = float(material["EX"])
A = float(np.asarray(real)[0])
C = np.outer(d, d) # 3 × 3
k = (E * A / L) * np.block([[C, -C], [-C, C]])
return np.ascontiguousarray(k)
[docs]
@staticmethod
def me(
coords: np.ndarray,
material: dict[str, float],
real: np.ndarray,
lumped: bool = False,
) -> np.ndarray:
L, _ = _axial_frame(np.asarray(coords, dtype=np.float64))
rho = float(material["DENS"])
A = float(np.asarray(real)[0])
m_total = rho * A * L
if lumped:
return (m_total / 2.0) * np.eye(6)
I3 = np.eye(3)
return (m_total / 6.0) * np.block([[2.0 * I3, I3], [I3, 2.0 * I3]])