LINK180 — 3D 2-node truss#

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.

MAPDL compatibility. Reproduces LINK180 at its default settings.


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]:

\[\varepsilon_\text{axial} = \frac{(\hat{d} \cdot u_J) - (\hat{d} \cdot u_I)}{L}.\]

Interpolating linearly between the two nodes gives the 1 × 6 strain-displacement row

\[B = \frac{1}{L} \begin{bmatrix} -l & -m & -n & l & m & n \end{bmatrix}.\]

Closed-form stiffness#

The internal strain energy is \(\tfrac{1}{2}\, u_e^\top K_e\, u_e\) with

\[K_e = \int_0^{L} B^\top E\, B\, A\,\mathrm{d}s = \frac{E A}{L} B^\top B.\]

Substituting \(B\) and collecting the 3 × 3 projector \(C_{jk} = \hat{d}_j \hat{d}_k\):

\[\begin{split}K_e = \frac{E A}{L} \begin{bmatrix} \phantom{-}C & -C \\ -C & \phantom{-}C \end{bmatrix}, \qquad C = \hat{d} \hat{d}^\top.\end{split}\]

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]:

\[\begin{split}M_e = \frac{\rho A L}{6} \begin{bmatrix} 2 I_3 & I_3 \\ I_3 & 2 I_3 \end{bmatrix},\end{split}\]

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 (LUMPM,ON in MAPDL, lumped=True in femorph-solver):

\[M_e^\text{lump} = \frac{\rho A L}{2} I_6.\]

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

real[0]

\(A\)

Cross-sectional area (mandatory).

real[1]

ADDMAS (added mass per unit length); parsed but not yet consumed by the mass kernel.

real[2]

ISTRN (initial strain); parsed but not yet consumed.

Validation#

Analytical natural frequency. A uniform bar clamped at one end has

\[f_n = \frac{(2n - 1)}{4 L}\,\sqrt{\frac{E}{\rho}}, \qquad n = 1, 2, 3, \dots\]

The mesh-refined LINK180 modal benchmark under tests.elements.link180 reproduces the first five frequencies to sub-percent.

API reference#

class femorph_solver.elements.link180.LINK180[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 MAPDL 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.

References#

[Cook2001] (1,2)

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

[ZTZ2013]

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