POINT_MASS — concentrated nodal mass#

Kinematics. Single-node element; three translational DOFs (\(u_x, u_y, u_z\)) in the current implementation.

Stiffness. Zero — a point mass contributes no stiffness.

Mass. Diagonal \(M_e = \text{diag}(m_x, m_y, m_z)\); lumped and consistent formulations coincide for a one-node element.

Variants. Translational-only is implemented today; the rotary- inertia variant (6 DOFs per node, with three angular-inertia components) is roadmap and will land alongside the shell / beam rotary-coupling path.


Discrete inertia as an element#

A point mass is not an element in the variational-FEM sense — there’s no shape function, no quadrature, no internal work. What it contributes is a discrete inertial term in the equation of motion [Cook2001]:

\[M_e\, \ddot{u}_e\ +\ (\text{stiffness contribution from other elements})\, =\ F_e(t),\]

with

\[M_e = \text{diag}\bigl(m_x,\ m_y,\ m_z\bigr),\qquad K_e = 0_{3\times 3}.\]

Each translational mass is passed as its own positive scalar:

  • real[0]\(m_x\) (mandatory).

  • real[1]\(m_y\) (defaults to \(m_x\) when omitted).

  • real[2]\(m_z\) (defaults to \(m_x\) when omitted).

That allows direction-dependent inertia (uncommon in structural mechanics but occasionally useful when a physical device with distinct axial / transverse inertia is being idealised as a lumped point).

Numpy walk-through#

import numpy as np

mass = 2.5     # kg, isotropic
M = mass * np.eye(3)             # 3 x 3 diagonal
K = np.zeros((3, 3))             # no stiffness

# Anisotropic variant:
M_aniso = np.diag([2.5, 2.5, 0.8])

In the global assembly the \(M_e\) block lands on the diagonal of the global \(M\) matrix at the DOFs of the attached node.

Validation#

Single-DOF oscillator. A POINT_MASS attached to a single SPRING element gives a one-DOF system with analytic natural frequency

\[f = \frac{1}{2\pi}\,\sqrt{\frac{K}{m}}.\]

tests.elements.mass21 checks this against the modal solver to machine precision.

Lumped-mass equivalence. For any element \(M_e^\text{lump} = \text{diag}(\sum_j M_{e,ij})\) is the row-sum lumping — for POINT_MASS the lumped and consistent forms are identical by construction.

API reference#

class femorph_solver.elements.point_mass.PointMass[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.

References#

[Cook2001]

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. 11 (consistent vs lumped mass; discrete inertia). https://www.wiley.com/en-us/Concepts+and+Applications+of+Finite+Element+Analysis%2C+4th+Edition-p-9780471356059