Skip to content

Masses — g.masses

Solver-agnostic nodal-mass definitions, records, and resolver. Masses are declared on geometry as concentrated values or densities (linear, areal, volumetric) and accumulated to nodes after meshing by g.mesh.queries.get_fem_data.

Two-stage pipeline

Stage 1 — declare before meshing. The four factories on g.masses (point, line, surface, volume) store MassDef dataclasses on geometric targets.

Stage 2 — resolve after meshing. MassResolver converts each def to per-node contributions, then the composite accumulates contributions across overlapping defs so each node ends up with at most one MassRecord. Records land on fem.nodes.masses as a MassSet (an iterable of MassRecord with total_mass() / by_node() helpers — see FEM Broker).

Each record carries a length-6 vector (mx, my, mz, Ixx, Iyy, Izz); the OpenSees bridge slices it to the model's ndf when emitting ops.mass(...) commands (rotational components are dropped for ndf < 4).

No patterns

Unlike loads, masses are not grouped under named patterns. Mass is intrinsic to the model — there is one nodal mass per node regardless of which load pattern is active. Multiple mass definitions targeting overlapping nodes simply accumulate.

Lumped vs. consistent reduction

Each factory accepts reduction="lumped" (default) or reduction="consistent":

  • lumped — element total mass split equally among the corner nodes. Diagonal mass matrix; cheap and stable for explicit dynamics.
  • consistent — proper consistent mass matrix integration. Today this is only meaningful for line elements (returns the ρ_l L / 6 · [[2, 1], [1, 2]] consistent matrix). Surface and volume paths fall through to lumped because tri3 / quad4 / tet4 / hex8 with constant density have the same diagonal-summed per-node share. The separate paths exist so higher-order elements (tri6, quad8, tet10, hex20) can be wired in without changing the public API.

Avoiding double-counting

apeGmsh emits explicit ops.mass(...) commands. If your OpenSees material or section also carries a non-zero rho, those contributions add to whatever this composite emits. Either:

  • keep rho=0 on the material and let g.masses carry all inertia, or
  • skip the matching volume / surface call and let the material handle it.

Pair g.masses.volume(..., density=ρ) with g.loads.gravity(..., density=ρ) for matching gravity body weight.

Target identification

Targets follow the same flexible scheme as LoadsComposite. Pass pg= / label= / tag= to bypass auto-resolution.

Worked example

from apeGmsh import apeGmsh

with apeGmsh(model_name="frame") as g:
    # ... geometry + Parts ...

    # Lumped mass at the top of a tower
    g.masses.point("Antenna", mass=350.0)

    # Cladding mass spread along an exterior edge (kg/m)
    g.masses.line("PerimeterEdge", linear_density=85.0)

    # Slab self-mass via shell areal density (ρ·t = kg/m²)
    g.masses.surface("Slab", areal_density=2400.0 * 0.20)

    # Steel column self-mass via material density (kg/m³)
    g.masses.volume("Columns", density=7850.0)

    g.mesh.generation.generate(dim=3)
    fem = g.mesh.queries.get_fem_data(dim=3)

    for m in fem.nodes.masses:
        ops.mass(m.node_id, *m.mass[:3])    # ndm=3 only

    print("Total mass:", fem.nodes.masses.total_mass())

Composite

apeGmsh.core.MassesComposite.MassesComposite

MassesComposite(parent: '_ApeGmshSession')

Solver-agnostic nodal-mass composite — declare on geometry, accumulate per-node mass after meshing.

Two-stage pipeline
  1. Declare (pre-mesh): the factory methods on this composite (:meth:point, :meth:line, :meth:surface, :meth:volume) store :class:~apeGmsh.solvers.Masses.MassDef dataclasses describing intent on geometric targets — concentrated lumps, line densities, areal densities, or material density on volumes.
  2. Resolve (post-mesh): :meth:resolve (called automatically by :meth:Mesh.queries.get_fem_data) walks the def list, hands each one to :class:~apeGmsh.solvers.Masses.MassResolver, and accumulates contributions across overlapping targets so each node ends up with at most one :class:~apeGmsh.solvers.Masses.MassRecord.

Resolved records land on fem.nodes.masses as a :class:~apeGmsh.mesh._record_set.MassSet. Each record carries a length-6 mass vector (mx, my, mz, Ixx, Iyy, Izz); downstream solver bridges slice it to the model's ndf (the rotational components are dropped for ndf < 4).

No patterns

Unlike :class:LoadsComposite, masses are not grouped by pattern. Mass is intrinsic to the model — there is one nodal mass per node regardless of which load pattern is active.

Reduction modes

Each factory accepts reduction="lumped" (default) or reduction="consistent":

  • lumped — each element's total mass is split equally among its corner nodes. Diagonal mass matrix; cheap and stable for explicit dynamics.
  • consistent — line elements use the proper ρ_l L / 6 · [[2,1],[1,2]] consistent matrix; surface and volume paths currently fall through to lumped because, for tri3 / quad4 / tet4 / hex8 with constant density, the consistent diagonal sum equals the lumped per-node share. The separate paths are kept so higher-order types (tri6, quad8, tet10, hex20) can be wired in without changing the public API.
Avoiding double-counting

apeGmsh always emits explicit ops.mass(node, mx, my, mz, …) commands. If your OpenSees material or section also carries a non-zero rho, those contributions add to whatever this composite emits. Either:

  • keep rho=0 on the material and let this composite carry all inertia, or
  • skip the matching :meth:volume / :meth:surface call and let the material handle it.
Target identification

Targets follow the same flexible scheme as :class:LoadsComposite:

  • a list of (dim, tag) tuples
  • a part label (g.parts.instances[label])
  • a physical group name (g.physical)
  • a mesh-selection name (g.mesh_selection)
  • a Tier 1 internal label

Pass pg= / label= / tag= to bypass auto-resolution and pin a specific source.

Examples

Declare three sources of mass and resolve to nodes::

with apeGmsh(model_name="frame") as g:
    # ... geometry + Parts ...

    # Lumped mass at the top of a tower
    g.masses.point("Antenna", mass=350.0)

    # Cladding mass spread along an exterior edge
    g.masses.line("Cladding", linear_density=120.0)

    # Self-mass of all column volumes
    g.masses.volume("Columns", density=7850.0)

    g.mesh.generation.generate(dim=3)
    fem = g.mesh.queries.get_fem_data(dim=3)

    for m in fem.nodes.masses:
        ops.mass(m.node_id, *m.mass[:3])    # ndm=3 only
    print("Total mass:", fem.nodes.masses.total_mass())
Source code in src/apeGmsh/core/MassesComposite.py
def __init__(self, parent: "_ApeGmshSession") -> None:
    self._parent = parent
    self.mass_defs: list[MassDef] = []
    self.mass_records: list[MassRecord] = []

point

point(target=None, *, pg=None, label=None, tag=None, mass: float, rotational: tuple | None = None, dofs: list[int] | None = None, reduction: str = 'lumped', name: str | None = None) -> PointMassDef

Concentrated mass at every node of target.

The same scalar mass (and optional rotational inertia triple) is applied to every targeted node. Useful for equipment, lumped fixtures, or any localised inertial contribution that doesn't come from a material density.

Resolution emits one :class:~apeGmsh.solvers.Masses.MassRecord per targeted node, accumulated with any other mass contributions on the same node.

Parameters

target : str or list of (dim, tag), optional Target node(s). See class docstring for the lookup order. pg, label, tag : Explicit-source overrides. mass : float Translational mass (per node) in model mass units. rotational : (Ixx, Iyy, Izz), optional Rotational inertia triple. Required for ndf >= 4 models that carry rotational DOFs. Default None stores zero rotational inertia. dofs : list[int], optional Translational DOF mask (subset of {1, 2, 3}). The scalar mass is applied only to the listed DOFs; the others receive zero translational mass. Default None applies mass to all three (1=ux, 2=uy, 3=uz). Rotational inertia (DOFs 4-6) is independent — set it via rotational=. reduction : "lumped" or "consistent", default "lumped" Has no effect for point masses (already lumped) — the argument is accepted for API symmetry with the distributed factories. name : str, optional Friendly name.

Returns

PointMassDef

Raises

KeyError If target doesn't resolve.

Examples

Equipment mass on a slab corner::

g.masses.point(
    "Equipment", mass=2500.0,
    rotational=(800.0, 800.0, 1200.0),
)
Source code in src/apeGmsh/core/MassesComposite.py
def point(
    self,
    target=None,
    *,
    pg=None, label=None, tag=None,
    mass: float,
    rotational: tuple | None = None,
    dofs: list[int] | None = None,
    reduction: str = "lumped",
    name: str | None = None,
) -> PointMassDef:
    """Concentrated mass at every node of *target*.

    The same scalar ``mass`` (and optional rotational inertia
    triple) is applied to every targeted node. Useful for
    equipment, lumped fixtures, or any localised inertial
    contribution that doesn't come from a material density.

    Resolution emits one
    :class:`~apeGmsh.solvers.Masses.MassRecord` per targeted
    node, accumulated with any other mass contributions on the
    same node.

    Parameters
    ----------
    target : str or list of (dim, tag), optional
        Target node(s). See class docstring for the lookup
        order.
    pg, label, tag :
        Explicit-source overrides.
    mass : float
        Translational mass (per node) in model mass units.
    rotational : (Ixx, Iyy, Izz), optional
        Rotational inertia triple. Required for ``ndf >= 4``
        models that carry rotational DOFs. Default ``None``
        stores zero rotational inertia.
    dofs : list[int], optional
        Translational DOF mask (subset of ``{1, 2, 3}``). The
        scalar ``mass`` is applied only to the listed DOFs; the
        others receive zero translational mass. Default ``None``
        applies mass to all three (1=ux, 2=uy, 3=uz). Rotational
        inertia (DOFs 4-6) is independent — set it via
        ``rotational=``.
    reduction : ``"lumped"`` or ``"consistent"``, default
        ``"lumped"``
        Has no effect for point masses (already lumped) — the
        argument is accepted for API symmetry with the
        distributed factories.
    name : str, optional
        Friendly name.

    Returns
    -------
    PointMassDef

    Raises
    ------
    KeyError
        If ``target`` doesn't resolve.

    Examples
    --------
    Equipment mass on a slab corner::

        g.masses.point(
            "Equipment", mass=2500.0,
            rotational=(800.0, 800.0, 1200.0),
        )
    """
    t, src = self._coalesce_target(target, pg=pg, label=label, tag=tag)
    return self._add_def(PointMassDef(
        target=t, target_source=src, name=name, reduction=reduction,
        dofs=_validate_translational_dofs(dofs),
        mass=mass, rotational=_validate_rotational(rotational),
    ))

line

line(target=None, *, pg=None, label=None, tag=None, linear_density: float, rotational: tuple | None = None, dofs: list[int] | None = None, reduction: str = 'lumped', name: str | None = None) -> LineMassDef

Distributed line mass on the curve(s) of target.

linear_density is in mass per unit length (e.g. kg/m). Each curve element contributes linear_density × edge_length of mass, distributed to its end nodes:

  • "lumped": split equally between the two end nodes (½ × ρ_l × L each).
  • "consistent": shape-function-integrated 2-node line mass matrix ρ_l L / 6 · [[2, 1], [1, 2]]. The diagonal-summed nodal share equals ρ_l L / 2 — same as lumped at the diagonal level — so the visible difference today is only in the coupling term emitted for solvers that consume it.

Use cases: cladding mass per unit length on shells, ducts or piping along a beam, façade panels along an edge.

Parameters

target : str or list of (dim, tag) Curve(s) carrying the line mass. pg, label, tag : Explicit-source overrides. linear_density : float Mass per unit length. rotational : (Ixx, Iyy, Izz), optional Fixed rotational inertia attached to every node receiving mass from this def. apeGmsh does not derive rotational inertia from linear_density — the user supplies it. Default None (no rotational mass). dofs : list[int], optional Translational DOF mask (subset of {1, 2, 3}). Default None applies to all three. reduction : "lumped" or "consistent", default "lumped" Lumping scheme. name : str, optional Friendly name.

Returns

LineMassDef

Raises

KeyError If target doesn't resolve.

Examples

Curtain-wall cladding along a building edge::

g.masses.line("PerimeterEdge", linear_density=85.0)
Source code in src/apeGmsh/core/MassesComposite.py
def line(
    self,
    target=None,
    *,
    pg=None, label=None, tag=None,
    linear_density: float,
    rotational: tuple | None = None,
    dofs: list[int] | None = None,
    reduction: str = "lumped",
    name: str | None = None,
) -> LineMassDef:
    """Distributed line mass on the curve(s) of *target*.

    ``linear_density`` is in mass per unit length (e.g. kg/m).
    Each curve element contributes ``linear_density × edge_length``
    of mass, distributed to its end nodes:

    * ``"lumped"``: split equally between the two end nodes
      (``½ × ρ_l × L`` each).
    * ``"consistent"``: shape-function-integrated 2-node line
      mass matrix ``ρ_l L / 6 · [[2, 1], [1, 2]]``. The
      diagonal-summed nodal share equals ``ρ_l L / 2`` — same
      as lumped at the diagonal level — so the visible
      difference today is only in the coupling term emitted
      for solvers that consume it.

    Use cases: cladding mass per unit length on shells, ducts
    or piping along a beam, façade panels along an edge.

    Parameters
    ----------
    target : str or list of (dim, tag)
        Curve(s) carrying the line mass.
    pg, label, tag :
        Explicit-source overrides.
    linear_density : float
        Mass per unit length.
    rotational : (Ixx, Iyy, Izz), optional
        Fixed rotational inertia attached to every node receiving
        mass from this def. apeGmsh does not derive rotational
        inertia from ``linear_density`` — the user supplies it.
        Default ``None`` (no rotational mass).
    dofs : list[int], optional
        Translational DOF mask (subset of ``{1, 2, 3}``). Default
        ``None`` applies to all three.
    reduction : ``"lumped"`` or ``"consistent"``, default
        ``"lumped"``
        Lumping scheme.
    name : str, optional
        Friendly name.

    Returns
    -------
    LineMassDef

    Raises
    ------
    KeyError
        If ``target`` doesn't resolve.

    Examples
    --------
    Curtain-wall cladding along a building edge::

        g.masses.line("PerimeterEdge", linear_density=85.0)
    """
    t, src = self._coalesce_target(target, pg=pg, label=label, tag=tag)
    return self._add_def(LineMassDef(
        target=t, target_source=src, name=name, reduction=reduction,
        dofs=_validate_translational_dofs(dofs),
        linear_density=linear_density,
        rotational=_validate_rotational(rotational),
    ))

surface

surface(target=None, *, pg=None, label=None, tag=None, areal_density: float, rotational: tuple | None = None, derive_rotational: bool = False, dofs: list[int] | None = None, reduction: str = 'lumped', name: str | None = None) -> SurfaceMassDef

Distributed surface mass on the face(s) of target.

areal_density is in mass per unit area (e.g. kg/m² or ρ × t for a shell of thickness t). Each face element contributes areal_density × face_area of mass, distributed to its corner nodes:

  • "lumped": equal split among corner nodes (ρ_a × A / n_corners each).
  • "consistent": for tri3 / quad4 with constant areal density the consistent diagonal sum equals lumped, so the resolver currently falls through to the lumped implementation. The separate path is reserved for future higher-order shells.

Use cases: slab self-mass when modelled as a 2-D shell, cladding panels, water mass on a deck, pavement.

Parameters

target : str or list of (dim, tag) Surface(s) carrying the areal mass. pg, label, tag : Explicit-source overrides. areal_density : float Mass per unit area. For a shell of constant thickness t and material density ρ, pass ρ * t. rotational : (Ixx, Iyy, Izz), optional Fixed rotational inertia attached to every node receiving mass from this def. apeGmsh does not derive rotational inertia from areal_density — the user supplies it. Default None (no rotational mass). dofs : list[int], optional Translational DOF mask (subset of {1, 2, 3}). Default None applies to all three. reduction : "lumped" or "consistent", default "lumped" Lumping scheme. name : str, optional Friendly name.

Returns

SurfaceMassDef

Examples

Concrete slab self-mass via shell areal density::

g.masses.surface(
    "Slab",
    areal_density=2400.0 * 0.20,   # ρ·t (kg/m²)
)
Source code in src/apeGmsh/core/MassesComposite.py
def surface(
    self,
    target=None,
    *,
    pg=None, label=None, tag=None,
    areal_density: float,
    rotational: tuple | None = None,
    derive_rotational: bool = False,
    dofs: list[int] | None = None,
    reduction: str = "lumped",
    name: str | None = None,
) -> SurfaceMassDef:
    """Distributed surface mass on the face(s) of *target*.

    ``areal_density`` is in mass per unit area (e.g. kg/m² or
    ``ρ × t`` for a shell of thickness ``t``). Each face element
    contributes ``areal_density × face_area`` of mass,
    distributed to its corner nodes:

    * ``"lumped"``: equal split among corner nodes
      (``ρ_a × A / n_corners`` each).
    * ``"consistent"``: for tri3 / quad4 with constant areal
      density the consistent diagonal sum equals lumped, so
      the resolver currently falls through to the lumped
      implementation. The separate path is reserved for
      future higher-order shells.

    Use cases: slab self-mass when modelled as a 2-D shell,
    cladding panels, water mass on a deck, pavement.

    Parameters
    ----------
    target : str or list of (dim, tag)
        Surface(s) carrying the areal mass.
    pg, label, tag :
        Explicit-source overrides.
    areal_density : float
        Mass per unit area. For a shell of constant thickness
        ``t`` and material density ``ρ``, pass ``ρ * t``.
    rotational : (Ixx, Iyy, Izz), optional
        Fixed rotational inertia attached to every node receiving
        mass from this def. apeGmsh does not derive rotational
        inertia from ``areal_density`` — the user supplies it.
        Default ``None`` (no rotational mass).
    dofs : list[int], optional
        Translational DOF mask (subset of ``{1, 2, 3}``). Default
        ``None`` applies to all three.
    reduction : ``"lumped"`` or ``"consistent"``, default
        ``"lumped"``
        Lumping scheme.
    name : str, optional
        Friendly name.

    Returns
    -------
    SurfaceMassDef

    Examples
    --------
    Concrete slab self-mass via shell areal density::

        g.masses.surface(
            "Slab",
            areal_density=2400.0 * 0.20,   # ρ·t (kg/m²)
        )
    """
    t, src = self._coalesce_target(target, pg=pg, label=label, tag=tag)
    return self._add_def(SurfaceMassDef(
        target=t, target_source=src, name=name, reduction=reduction,
        dofs=_validate_translational_dofs(dofs),
        areal_density=areal_density,
        rotational=_validate_rotational(rotational),
        derive_rotational=_validate_derive_rotational(
            derive_rotational, rotational=rotational,
            reduction=reduction,
        ),
    ))

volume

volume(target=None, *, pg=None, label=None, tag=None, density: float, rotational: tuple | None = None, derive_rotational: bool = False, dofs: list[int] | None = None, reduction: str = 'lumped', name: str | None = None) -> VolumeMassDef

Distributed mass from material density over the volume(s) of target.

density is in mass per unit volume (e.g. kg/m³). Each volume element contributes density × element_volume of mass, distributed to its corner nodes:

  • "lumped": equal split among the element's corner nodes (ρ × V / n_corners each).
  • "consistent": for tet4 / hex8 with constant density the consistent diagonal sum equals lumped, so the resolver currently falls through to the lumped implementation. The separate path is reserved for future higher-order solid types.

Pair this with a matching :meth:loads.gravity call for gravitational body weight. Don't double-count — see the class docstring on avoiding double-counting when the OpenSees material also carries rho.

Parameters

target : str or list of (dim, tag) Volume(s) carrying the mass. pg, label, tag : Explicit-source overrides. density : float Material density (mass per unit volume). rotational : (Ixx, Iyy, Izz), optional Fixed rotational inertia attached to every node receiving mass from this def. apeGmsh does not derive rotational inertia from density × volume — the user supplies it (deriving it properly would need a moment-of-inertia integration over the element). Default None. dofs : list[int], optional Translational DOF mask (subset of {1, 2, 3}). Default None applies to all three (1=ux, 2=uy, 3=uz). reduction : "lumped" or "consistent", default "lumped" Lumping scheme. name : str, optional Friendly name.

Returns

VolumeMassDef

See Also

loads.gravity : Apply matching gravitational body weight.

Examples

Steel column self-mass::

g.masses.volume("Columns", density=7850.0)
Source code in src/apeGmsh/core/MassesComposite.py
def volume(
    self,
    target=None,
    *,
    pg=None, label=None, tag=None,
    density: float,
    rotational: tuple | None = None,
    derive_rotational: bool = False,
    dofs: list[int] | None = None,
    reduction: str = "lumped",
    name: str | None = None,
) -> VolumeMassDef:
    """Distributed mass from material density over the volume(s)
    of *target*.

    ``density`` is in mass per unit volume (e.g. kg/m³). Each
    volume element contributes ``density × element_volume`` of
    mass, distributed to its corner nodes:

    * ``"lumped"``: equal split among the element's corner
      nodes (``ρ × V / n_corners`` each).
    * ``"consistent"``: for tet4 / hex8 with constant density
      the consistent diagonal sum equals lumped, so the
      resolver currently falls through to the lumped
      implementation. The separate path is reserved for
      future higher-order solid types.

    Pair this with a matching :meth:`loads.gravity` call for
    gravitational body weight. **Don't double-count** — see the
    class docstring on avoiding double-counting when the
    OpenSees material also carries ``rho``.

    Parameters
    ----------
    target : str or list of (dim, tag)
        Volume(s) carrying the mass.
    pg, label, tag :
        Explicit-source overrides.
    density : float
        Material density (mass per unit volume).
    rotational : (Ixx, Iyy, Izz), optional
        Fixed rotational inertia attached to every node receiving
        mass from this def. apeGmsh does not derive rotational
        inertia from ``density × volume`` — the user supplies it
        (deriving it properly would need a moment-of-inertia
        integration over the element). Default ``None``.
    dofs : list[int], optional
        Translational DOF mask (subset of ``{1, 2, 3}``). Default
        ``None`` applies to all three (1=ux, 2=uy, 3=uz).
    reduction : ``"lumped"`` or ``"consistent"``, default
        ``"lumped"``
        Lumping scheme.
    name : str, optional
        Friendly name.

    Returns
    -------
    VolumeMassDef

    See Also
    --------
    loads.gravity : Apply matching gravitational body weight.

    Examples
    --------
    Steel column self-mass::

        g.masses.volume("Columns", density=7850.0)
    """
    t, src = self._coalesce_target(target, pg=pg, label=label, tag=tag)
    return self._add_def(VolumeMassDef(
        target=t, target_source=src, name=name, reduction=reduction,
        dofs=_validate_translational_dofs(dofs),
        density=density,
        rotational=_validate_rotational(rotational),
        derive_rotational=_validate_derive_rotational(
            derive_rotational, rotational=rotational,
            reduction=reduction,
        ),
    ))

validate_pre_mesh

validate_pre_mesh() -> None

Validate every registered mass's target can be resolved.

Called by :meth:Mesh.generate before meshing so typos fail fast. Raw (dim, tag) lists are skipped.

Source code in src/apeGmsh/core/MassesComposite.py
def validate_pre_mesh(self) -> None:
    """Validate every registered mass's target can be resolved.

    Called by :meth:`Mesh.generate` before meshing so typos fail
    fast.  Raw ``(dim, tag)`` lists are skipped.
    """
    for defn in self.mass_defs:
        target = defn.target
        if not isinstance(target, str):
            continue
        self._resolve_target(target, source=defn.target_source)

resolve

resolve(node_tags, node_coords, elem_tags=None, connectivity=None, *, node_map=None, face_map=None, ndf: int = 6) -> MassSet

Resolve all stored MassDefs into a :class:MassSet.

Multiple definitions targeting overlapping nodes are accumulated — each node ends up with at most one :class:MassRecord whose vector is the sum of contributions.

Source code in src/apeGmsh/core/MassesComposite.py
def resolve(
    self,
    node_tags,
    node_coords,
    elem_tags=None,
    connectivity=None,
    *,
    node_map=None,
    face_map=None,
    ndf: int = 6,
) -> MassSet:
    """Resolve all stored MassDefs into a :class:`MassSet`.

    Multiple definitions targeting overlapping nodes are
    accumulated — each node ends up with at most one
    :class:`MassRecord` whose vector is the sum of contributions.
    """
    resolver = MassResolver(
        node_tags, node_coords, elem_tags, connectivity, ndf=ndf,
    )
    all_nodes = set(int(t) for t in node_tags)

    # Per-node accumulator across all defs
    accum: dict[int, np.ndarray] = {}

    for defn in self.mass_defs:
        cfg = _DISPATCH[type(defn)]
        method_name = cfg.get(defn.reduction)
        if method_name is None:
            raise ValueError(
                f"{type(defn).__name__} does not support "
                f"reduction={defn.reduction!r}"
            )
        method = getattr(self, method_name)
        raw_records = method(resolver, defn, node_map, all_nodes)
        for r in raw_records:
            vec = np.asarray(r.mass, dtype=float)
            if r.node_id in accum:
                accum[r.node_id] += vec
            else:
                accum[r.node_id] = vec.copy()

    # Build final flattened MassRecord list (one per node)
    records: list[MassRecord] = [
        MassRecord(
            node_id=int(nid),
            mass=tuple(float(v) for v in vec),
        )
        for nid, vec in sorted(accum.items())
    ]
    self.mass_records = records
    return MassSet(records)

Definitions

apeGmsh._kernel.defs.masses.MassDef dataclass

MassDef(kind: str, target: object, name: str | None = None, reduction: str = 'lumped', target_source: str = 'auto', dofs: list[int] | None = None)

Base class for all mass definitions.

dofs masks which translational DOFs (subset of {1, 2, 3}) receive the resolved mass. None (default) applies mass to all three translational components. Rotational DOFs (4, 5, 6) are controlled separately via rotational=(Ixx, Iyy, Izz) on the individual def subclasses — dofs deliberately does not cover them because density × volume (and the line / area analogues) cannot derive rotational inertia without element-shape awareness.

apeGmsh._kernel.defs.masses.PointMassDef dataclass

PointMassDef(kind: str, target: object, name: str | None = None, reduction: str = 'lumped', target_source: str = 'auto', dofs: list[int] | None = None, mass: float = 0.0, rotational: tuple[float, float, float] | None = None)

Bases: MassDef

Concentrated mass at a node (or set of nodes).

The same scalar mass (and optional rotational inertia) is applied to every node in the target. Useful for representing equipment, point fixtures, or any localised lumped mass.

apeGmsh._kernel.defs.masses.LineMassDef dataclass

LineMassDef(kind: str, target: object, name: str | None = None, reduction: str = 'lumped', target_source: str = 'auto', dofs: list[int] | None = None, linear_density: float = 0.0, rotational: tuple[float, float, float] | None = None)

Bases: MassDef

Distributed line mass — linear density along curve(s).

linear_density is in mass per unit length (e.g. kg/m). The total mass per curve is linear_density × curve_length and is distributed to the curve's mesh nodes.

rotational (optional) attaches a fixed rotational inertia triple (Ixx, Iyy, Izz) to every node receiving mass from this def. apeGmsh does not derive rotational inertia from linear_density — the user supplies the value. Useful when an ndf=6 model needs rotational inertia along an edge (e.g. cable with cross-sectional moment of inertia).

apeGmsh._kernel.defs.masses.SurfaceMassDef dataclass

SurfaceMassDef(kind: str, target: object, name: str | None = None, reduction: str = 'lumped', target_source: str = 'auto', dofs: list[int] | None = None, areal_density: float = 0.0, rotational: tuple[float, float, float] | None = None, derive_rotational: bool = False)

Bases: MassDef

Distributed surface mass — areal density on face(s).

areal_density is in mass per unit area (e.g. kg/m²). The total mass per face is areal_density × face_area distributed to the face's mesh nodes.

rotational (optional) attaches a fixed rotational inertia triple (Ixx, Iyy, Izz) to every node receiving mass from this def. apeGmsh does not derive rotational inertia from areal_density — the user supplies the value.

apeGmsh._kernel.defs.masses.VolumeMassDef dataclass

VolumeMassDef(kind: str, target: object, name: str | None = None, reduction: str = 'lumped', target_source: str = 'auto', dofs: list[int] | None = None, density: float = 0.0, rotational: tuple[float, float, float] | None = None, derive_rotational: bool = False)

Bases: MassDef

Distributed volume mass — material density on volume(s).

density is in mass per unit volume (e.g. kg/m³). The total mass per element is density × element_volume distributed to the element's nodes.

rotational (optional) attaches a fixed rotational inertia triple (Ixx, Iyy, Izz) to every node receiving mass from this def. apeGmsh does not derive rotational inertia from density × volume — that would need a moment-of-inertia integration over the element, which is element-shape aware. Until that's implemented, this kwarg lets the user attach a pre-computed rotational inertia uniformly.

derive_rotational (optional, default False) instead computes per-node rotational inertia from the element shape functions: I_xx^(I) = ∫ρ N_I (y²+z²) dV − M_I (y_I²+z_I²) (and cyclic), the about-node parallel-axis form whose assembled sum reproduces the continuum rigid-rotation kinetic energy exactly. Requires reduction='consistent' and is mutually exclusive with a fixed rotational= tuple. Per-node emitted values may be negative (they are parallel-axis corrections, not standalone inertias); the physical guarantee is on the assembled total.

Note

The user is responsible for setting the OpenSees material's rho=0 (or equivalent in other solvers) to avoid double counting. This composite always emits explicit nodal mass via ops.mass(...) commands.

Resolved record

apeGmsh._kernel.records._masses.MassRecord dataclass

MassRecord(node_id: int = 0, mass: tuple = (0.0, 0.0, 0.0, 0.0, 0.0, 0.0), name: str | None = None)

Resolved per-node mass entry.

Always length 6: (mx, my, mz, Ixx, Iyy, Izz). The OpenSees bridge slices to ndf when emitting commands (the rotational components are dropped for ndf<4 models).

Multiple :class:MassDef may contribute to the same node — the composite accumulates them so each node gets at most one :class:MassRecord in the final :class:MassSet.

Resolver

apeGmsh._kernel.resolvers._mass_resolver.MassResolver

MassResolver(node_tags: ndarray, node_coords: ndarray, elem_tags: ndarray | None = None, connectivity: ndarray | None = None, ndf: int = 6)

Convert :class:MassDef instances to :class:MassRecord lists.

Pure mesh math — receives raw arrays, returns record lists. The composite handles target -> mesh-entity resolution before calling these methods.

Source code in src/apeGmsh/_kernel/resolvers/_mass_resolver.py
def __init__(
    self,
    node_tags: ndarray,
    node_coords: ndarray,
    elem_tags: ndarray | None = None,
    connectivity: ndarray | None = None,
    ndf: int = 6,
) -> None:
    self.node_tags = np.asarray(node_tags, dtype=np.int64)
    self.node_coords = np.asarray(node_coords, dtype=np.float64).reshape(-1, 3)
    self.elem_tags = (
        np.asarray(elem_tags, dtype=np.int64) if elem_tags is not None else None
    )
    self.connectivity = (
        np.asarray(connectivity, dtype=np.int64) if connectivity is not None else None
    )
    self.ndf = int(ndf)
    self._node_to_idx = {int(t): i for i, t in enumerate(self.node_tags)}

element_volume

element_volume(conn_row: ndarray) -> float

Volume of one solid element.

  • n == 4 (tet4): exact analytic scalar triple product.
  • n == 8 (hex8): exact 6-tetrahedron decomposition.
  • any other catalog type (wedge6, tet10, hex20, hex27): isoparametric V = ∫_{Ω_ref} |J(ξ)| dξ via the shared shape-function Jacobian + reference quadrature — exact for affine elements, the standard high-accuracy approximation for curved higher-order ones.
  • unknown element type: bounding-box last resort.
Source code in src/apeGmsh/_kernel/resolvers/_mass_resolver.py
def element_volume(self, conn_row: ndarray) -> float:
    """Volume of one solid element.

    * ``n == 4`` (tet4): exact analytic scalar triple product.
    * ``n == 8`` (hex8): exact 6-tetrahedron decomposition.
    * any other catalog type (wedge6, tet10, hex20, hex27):
      isoparametric ``V = ∫_{Ω_ref} |J(ξ)| dξ`` via the shared
      shape-function Jacobian + reference quadrature — exact for
      affine elements, the standard high-accuracy approximation
      for curved higher-order ones.
    * unknown element type: bounding-box last resort.
    """
    n = len(conn_row)
    pts = np.array([self.coords_of(int(nid)) for nid in conn_row])
    if n == 4:
        v = np.abs(np.dot(pts[1] - pts[0], np.cross(pts[2] - pts[0], pts[3] - pts[0]))) / 6.0
        return float(v)
    if n == 8:
        tets = [
            (0, 1, 2, 5), (0, 2, 3, 7), (0, 5, 2, 6),
            (0, 5, 6, 7), (0, 7, 2, 6), (0, 4, 5, 7),
        ]
        tot = 0.0
        for a, b, c, d in tets:
            tot += np.abs(np.dot(
                pts[b] - pts[a],
                np.cross(pts[c] - pts[a], pts[d] - pts[a]),
            )) / 6.0
        return float(tot)
    code = volume_code(n)
    if code is not None:
        catalog = get_shape_functions(code)
        if catalog is not None:
            _, dN_fn, geom_kind, _ = catalog
            qp, qw = reference_quadrature(code)
            detJ = compute_jacobian_dets(
                qp, pts[None, :, :], dN_fn, geom_kind,
            )[0]
            return float(np.sum(qw * detJ))
    # Unknown element type (pyramid, wedge15, …) — bbox last resort.
    mn, mx = pts.min(axis=0), pts.max(axis=0)
    return float(np.prod(mx - mn))

resolve_point_lumped

resolve_point_lumped(defn: PointMassDef, node_set: set[int]) -> list[MassRecord]

Apply the same point mass to every node in node_set.

Honors defn.dofs (translational mask) and defn.rotational (rotational inertia tuple).

Source code in src/apeGmsh/_kernel/resolvers/_mass_resolver.py
def resolve_point_lumped(
    self,
    defn: PointMassDef,
    node_set: set[int],
) -> list[MassRecord]:
    """Apply the same point mass to every node in *node_set*.

    Honors ``defn.dofs`` (translational mask) and
    ``defn.rotational`` (rotational inertia tuple).
    """
    vec = _build_vec6(float(defn.mass), defn.dofs, defn.rotational)
    accum: dict[int, ndarray] = {}
    for nid in node_set:
        _accumulate(accum, nid, vec)
    return _accum_to_records(accum, name=defn.name)

resolve_line_lumped

resolve_line_lumped(defn: LineMassDef, edges: list[tuple[int, int]]) -> list[MassRecord]

Distribute line mass by tributary length.

Each node receives ρₗ × Σ(adjacent_edge_length / 2) of translational mass on the DOFs listed in defn.dofs (default: 1, 2, 3). defn.rotational (if set) attaches a fixed rotational inertia to every receiving node.

Source code in src/apeGmsh/_kernel/resolvers/_mass_resolver.py
def resolve_line_lumped(
    self,
    defn: LineMassDef,
    edges: list[tuple[int, int]],
) -> list[MassRecord]:
    """Distribute line mass by tributary length.

    Each node receives ``ρₗ × Σ(adjacent_edge_length / 2)`` of
    translational mass on the DOFs listed in ``defn.dofs``
    (default: 1, 2, 3).  ``defn.rotational`` (if set) attaches
    a fixed rotational inertia to every receiving node.
    """
    rho_l = float(defn.linear_density)
    dofs = defn.dofs
    rot = defn.rotational
    accum: dict[int, ndarray] = {}
    for n1, n2 in edges:
        half_L = 0.5 * self.edge_length(n1, n2)
        m_share = rho_l * half_L
        vec = _build_vec6(m_share, dofs, rot)
        _accumulate(accum, n1, vec)
        _accumulate(accum, n2, vec)
    return _accum_to_records(accum, name=defn.name)

resolve_surface_lumped

resolve_surface_lumped(defn: SurfaceMassDef, faces: list[list[int]]) -> list[MassRecord]

Distribute surface mass by tributary area.

Honors defn.dofs (translational mask) and defn.rotational (rotational inertia tuple).

Source code in src/apeGmsh/_kernel/resolvers/_mass_resolver.py
def resolve_surface_lumped(
    self,
    defn: SurfaceMassDef,
    faces: list[list[int]],
) -> list[MassRecord]:
    """Distribute surface mass by tributary area.

    Honors ``defn.dofs`` (translational mask) and
    ``defn.rotational`` (rotational inertia tuple).
    """
    rho_a = float(defn.areal_density)
    dofs = defn.dofs
    rot = defn.rotational
    accum: dict[int, ndarray] = {}
    for face in faces:
        A = self.face_area(face)
        if A <= 0:
            continue
        m_share = rho_a * A / len(face)
        vec = _build_vec6(m_share, dofs, rot)
        for nid in face:
            _accumulate(accum, int(nid), vec)
    return _accum_to_records(accum, name=defn.name)

resolve_volume_lumped

resolve_volume_lumped(defn: VolumeMassDef, elements: list[ndarray]) -> list[MassRecord]

Distribute volume mass equally to element nodes (lumped).

Honors defn.dofs (translational mask) and defn.rotational (rotational inertia tuple).

Source code in src/apeGmsh/_kernel/resolvers/_mass_resolver.py
def resolve_volume_lumped(
    self,
    defn: VolumeMassDef,
    elements: list[ndarray],
) -> list[MassRecord]:
    """Distribute volume mass equally to element nodes (lumped).

    Honors ``defn.dofs`` (translational mask) and
    ``defn.rotational`` (rotational inertia tuple).
    """
    rho = float(defn.density)
    dofs = defn.dofs
    rot = defn.rotational
    accum: dict[int, ndarray] = {}
    for conn_row in elements:
        V = self.element_volume(conn_row)
        if V <= 0:
            continue
        m_share = rho * V / len(conn_row)
        vec = _build_vec6(m_share, dofs, rot)
        for nid in conn_row:
            _accumulate(accum, int(nid), vec)
    return _accum_to_records(accum, name=defn.name)

resolve_point_consistent

resolve_point_consistent(defn: PointMassDef, node_set: set[int]) -> list[MassRecord]

Point mass is unambiguous — same as lumped.

Source code in src/apeGmsh/_kernel/resolvers/_mass_resolver.py
def resolve_point_consistent(
    self,
    defn: PointMassDef,
    node_set: set[int],
) -> list[MassRecord]:
    """Point mass is unambiguous — same as lumped."""
    return self.resolve_point_lumped(defn, node_set)

resolve_line_consistent

resolve_line_consistent(defn: LineMassDef, edges: list[tuple[int, int]]) -> list[MassRecord]

Consistent line mass for 2-node line element.

For a line element with constant linear density ρₗ and length L::

M_consistent = ρₗL/6 · [[2, 1], [1, 2]]

Each node receives the diagonal entry (2ρₗL/6) plus a coupling term to the other node (ρₗL/6). Lumped is M = ρₗL/2 · I, so the diagonal sum is the same — but consistent has off-diagonal coupling.

Per-node accumulation matches lumped at the diagonal level (2ρₗL/6 + ρₗL/6 = ρₗL/2). This implementation emits the diagonal-equivalent (which matches lumped numerically for 2-node lines) — true off-diagonal coupling requires the solver to support consistent element mass via -cMass or equivalent.

Source code in src/apeGmsh/_kernel/resolvers/_mass_resolver.py
def resolve_line_consistent(
    self,
    defn: LineMassDef,
    edges: list[tuple[int, int]],
) -> list[MassRecord]:
    """Consistent line mass for 2-node line element.

    For a line element with constant linear density ρₗ and length L::

        M_consistent = ρₗL/6 · [[2, 1], [1, 2]]

    Each node receives the diagonal entry (2ρₗL/6) plus a coupling
    term to the other node (ρₗL/6).  Lumped is `M = ρₗL/2 · I`,
    so the diagonal sum is the same — but consistent has
    off-diagonal coupling.

    Per-node accumulation matches lumped at the diagonal level
    (2ρₗL/6 + ρₗL/6 = ρₗL/2).  This implementation emits the
    diagonal-equivalent (which matches lumped numerically for
    2-node lines) — true off-diagonal coupling requires the
    solver to support consistent element mass via ``-cMass``
    or equivalent.
    """
    rho_l = float(defn.linear_density)
    accum: dict[int, ndarray] = {}
    for n1, n2 in edges:
        L = self.edge_length(n1, n2)
        # Consistent diagonal contribution: (2/6 + 1/6)·ρₗL = ρₗL/2
        m_share = 0.5 * rho_l * L
        vec = np.array([m_share, m_share, m_share, 0.0, 0.0, 0.0])
        _accumulate(accum, n1, vec)
        _accumulate(accum, n2, vec)
    return _accum_to_records(accum, name=defn.name)

resolve_surface_consistent

resolve_surface_consistent(defn: SurfaceMassDef, faces: list[list[int]]) -> list[MassRecord]

HRZ-lumped surface mass.

Distributes ρ_a · A per face to its nodes using HRZ weights keyed by node count (tri3/quad4 → equal-share by construction; tri6/quad8/quad9 → mid-side nodes correctly weighted). Honors defn.dofs and defn.rotational.

Source code in src/apeGmsh/_kernel/resolvers/_mass_resolver.py
def resolve_surface_consistent(
    self,
    defn: SurfaceMassDef,
    faces: list[list[int]],
) -> list[MassRecord]:
    """HRZ-lumped surface mass.

    Distributes ``ρ_a · A`` per face to its nodes using HRZ
    weights keyed by node count (tri3/quad4 → equal-share by
    construction; tri6/quad8/quad9 → mid-side nodes correctly
    weighted).  Honors ``defn.dofs`` and ``defn.rotational``.
    """
    rho_a = float(defn.areal_density)
    dofs = defn.dofs
    rot = defn.rotational
    derive = getattr(defn, "derive_rotational", False)
    accum: dict[int, ndarray] = {}
    for face in faces:
        A = self.face_area(face)
        if A <= 0:
            continue
        code = surface_code(len(face))
        if derive:
            self._hrz_distribute_derived(
                accum, face, rho_a, rho_a * A, code, dofs,
            )
        else:
            self._hrz_distribute(
                accum, face, rho_a * A, code, dofs, rot,
            )
    return _accum_to_records(accum, name=defn.name)

resolve_volume_consistent

resolve_volume_consistent(defn: VolumeMassDef, elements: list[ndarray]) -> list[MassRecord]

HRZ-lumped volume mass.

Distributes ρ · V per element to its nodes using HRZ weights keyed by node count (tet4/hex8/wedge6 → equal-share by construction; tet10/hex20/hex27 → higher-order nodes correctly weighted). Honors defn.dofs and defn.rotational.

The element volume comes from :meth:element_volume, now isoparametric (∫|J|dξ) for every catalog element type — exact for tet4 / hex8 / affine higher-order, high-accuracy for curved higher-order. The HRZ distribution is correct regardless.

Source code in src/apeGmsh/_kernel/resolvers/_mass_resolver.py
def resolve_volume_consistent(
    self,
    defn: VolumeMassDef,
    elements: list[ndarray],
) -> list[MassRecord]:
    """HRZ-lumped volume mass.

    Distributes ``ρ · V`` per element to its nodes using HRZ
    weights keyed by node count (tet4/hex8/wedge6 → equal-share by
    construction; tet10/hex20/hex27 → higher-order nodes correctly
    weighted).  Honors ``defn.dofs`` and ``defn.rotational``.

    The element volume comes from :meth:`element_volume`, now
    isoparametric (``∫|J|dξ``) for every catalog element type —
    exact for tet4 / hex8 / affine higher-order, high-accuracy for
    curved higher-order.  The HRZ *distribution* is correct
    regardless.
    """
    rho = float(defn.density)
    dofs = defn.dofs
    rot = defn.rotational
    derive = getattr(defn, "derive_rotational", False)
    accum: dict[int, ndarray] = {}
    for conn_row in elements:
        V = self.element_volume(conn_row)
        if V <= 0:
            continue
        code = volume_code(len(conn_row))
        if derive:
            self._hrz_distribute_derived(
                accum, conn_row, rho, rho * V, code, dofs,
            )
        else:
            self._hrz_distribute(
                accum, conn_row, rho * V, code, dofs, rot,
            )
    return _accum_to_records(accum, name=defn.name)