TRUSS2 — 3-D 2-node axial bar#
Kinematics. Two-node straight bar carrying axial load only; three translational DOFs per node; 6 DOFs per element. No bending, no torsion.
Stiffness. Closed-form \(K_e = \frac{EA}{L} \begin{bmatrix} C & -C \\ -C & C \end{bmatrix}\) with the projector \(C = \hat{d} \hat{d}^\top\).
Mass. Consistent or lumped, both closed-form.
Truss kinematics#
A straight bar between nodes \(I\) and \(J\) with unit direction vector \(\hat{d} = (l, m, n)\) and length \(L = \lVert x_J - x_I\rVert\) stores only axial strain [Cook2001] [ZTZ2013]:
Interpolating linearly between the two nodes gives the 1 × 6 strain-displacement row
Closed-form stiffness#
The internal strain energy is \(\tfrac{1}{2}\, u_e^\top K_e\, u_e\) with
Substituting \(B\) and collecting the 3 × 3 projector \(C_{jk} = \hat{d}_j \hat{d}_k\):
The rank of \(K_e\) is 1 (one axial stretching mode); the five remaining eigenvalues are exact zeros (two transverse translations of each end, plus a rigid-body axial translation). That’s expected for a bar — it carries no bending or shear, and any transverse force passes straight through.
Consistent mass#
Uniform-density integration of \(N^\top N\) along \(s\) gives [Cook2001]:
so each node gets translational inertia \(\frac{\rho A L}{3}\) from its own shape-function pair plus \(\frac{\rho A L}{6}\) from the coupling with the other node.
Lumped mass (lumped=True on the assembler):
Numpy walk-through#
import numpy as np
E, A = 2.1e11, 1e-4 # Pa, m^2
X = np.array([[0.0, 0.0, 0.0],
[1.0, 0.0, 0.0]])
d = X[1] - X[0]
L = float(np.linalg.norm(d))
dhat = d / L
# Projector onto the axial direction.
Cproj = np.outer(dhat, dhat) # 3x3
K = (E * A / L) * np.block([[ Cproj, -Cproj],
[-Cproj, Cproj]])
# K has rank 1; eigvals [0, 0, 0, 0, 0, 2 E A / L].
rho = 7850.0
M_cons = (rho * A * L / 6) * np.block([[2 * np.eye(3), np.eye(3)],
[ np.eye(3), 2 * np.eye(3)]])
M_lump = (rho * A * L / 2) * np.eye(6)
Real constants#
Slot |
Symbol |
Meaning |
|---|---|---|
|
\(A\) |
Cross-sectional area (mandatory). |
|
— |
|
|
— |
|
Validation#
Analytical natural frequency. A uniform bar clamped at one end has
The mesh-refined TRUSS2 modal benchmark under
tests.elements.link180 reproduces the first five frequencies
to sub-percent.
API reference#
- class femorph_solver.elements.truss2.Truss2[source]#
Bases:
ElementBase- static ke(coords: ndarray, material: dict[str, float], real: ndarray) ndarray[source]#
Return the element stiffness matrix in global DOF ordering.
- Parameters:
coords ((n_nodes, 3) float64)
material (mapping with property names as keys (
"EX","PRXY","DENS", …))real (1-D array of real constants)
- static me(coords: ndarray, material: dict[str, float], real: ndarray, lumped: bool = False) ndarray[source]#
Return the element mass matrix in global DOF ordering.
- static strain_batch(coords: ndarray, u_e: ndarray, material: dict[str, float] | None = None) ndarray | None[source]#
Per-element 3D Voigt strain at the two end nodes.
The truss carries only an axial strain
eps_a = (u2 - u1) . e / L- the perpendicular cross-section stretches under Poisson contraction when the bar is in uniaxial stress. Returning the full 3D continuum strain that corresponds to the pure-uniaxial stress statesigma = E*eps_a * (e^e)letscompute_nodal_stress()recover the right stress via the standard isotropicC @ epscontraction without having to special-case 1D kernels.For a uniaxial stress state along unit axis
ewith straineps_a(computed from the displacement field), the 3D strain tensor is:eps = eps_a * ((1 + nu) * (e^e) - nu * I)
Voigt order is
[exx, eyy, ezz, gxy, gyz, gxz](engineering shears), matchingMaterialTable.eval_c_isotropic(). Returns the per-element-node strain (same value at both nodes since the truss is constant-strain).
- static stress_update_static_nonlinear(coords: ndarray, material: dict[str, float], real: ndarray, u_e_global: ndarray, state, plastic)[source]#
Phase-2 plasticity hook — uniaxial return-map on a TRUSS2.
- Returns:
K_t_global ((6, 6)) – Consistent algorithmic tangent in the global frame (block
e_hat ⊗ e_hatscaled byD_t·A/L).F_int_global ((6,)) – Internal-force vector
±σ·A·e_hatat each node.new_state – Trial state after this Newton step’s return-map.
- static thermal_load(coords: ndarray, real: ndarray, material: dict[str, float], *, delta_t: float) ndarray[source]#
Equivalent nodal force vector for a uniform thermal pre-strain.
For an axial bar with thermal-expansion coefficient
α = material["ALPX"]and uniform temperature changeΔT, the thermal pre-strain isε_th = α·ΔT. The work-equivalent nodal force vector follows fromf_th = ∫ Bᵀ E ε_th dV:f_local_axial = E · A · α · ΔT · [-1, +1] (along axis)
i.e. an outward force
E·A·α·ΔTat each end pulling the constrained bar apart on heating (or pushing it together on cooling). The return value is the 6-component global-frame force vector at[ux_I, uy_I, uz_I, ux_J, uy_J, uz_J], already rotated by the same direction-cosine matrixke()uses, so the caller can scatter it directly intoF.Returns
np.zeros(6)when the material has noALPXset — a missing thermal-expansion coefficient just means thermally inert, not an error.References
Cook, R.D., Malkus, D.S., Plesha, M.E., Witt, R.J., Concepts and Applications of Finite Element Analysis, 4th ed., Wiley, 2002, §3.6.4 (thermal pre-strain), §6.10 (initial-strain equivalent nodal forces).
References#
Cook, R. D., Malkus, D. S., Plesha, M. E., and Witt, R. J. (2001). Concepts and Applications of Finite Element Analysis, 4th ed. Wiley. Ch. 2 (truss elements). https://www.wiley.com/en-us/Concepts+and+Applications+of+Finite+Element+Analysis%2C+4th+Edition-p-9780471356059
Zienkiewicz, O. C., Taylor, R. L., and Zhu, J. Z. (2013). The Finite Element Method: Its Basis and Fundamentals, 7th ed. Butterworth-Heinemann. Ch. 1 (discretisation of a bar). https://doi.org/10.1016/C2009-0-24909-9