Loads — g.loads¶
Solver-agnostic load definitions, records, and resolver. Loads are
declared on geometry (with optional pattern grouping) and
resolved on the mesh by g.mesh.queries.get_fem_data.
Two-stage pipeline¶
Stage 1 — declare before meshing. The factory methods on
g.loads (point, point_closest, line, surface, gravity,
body, face_load, face_sp) store
LoadDef dataclasses describing
intent at the geometry level. The active
pattern
context tags every def created inside it.
Stage 2 — resolve after meshing.
LoadResolver converts each
def to a list of resolved records. Records land on the FEM broker
according to type:
| Record family | Lives on | Emitted by |
|---|---|---|
NodalLoadRecord |
fem.nodes.loads |
tributary / consistent reductions |
ElementLoadRecord |
fem.elements.loads |
target_form="element" (eleLoad style) |
SPRecord |
fem.nodes.sp |
face_sp only |
Patterns¶
Loads (and only loads — not constraints, not masses) are grouped
under named patterns via the
pattern
context manager:
with g.loads.pattern("Dead"):
g.loads.gravity("Slab", density=2400)
g.loads.line("BeamEdge", magnitude=-15e3)
with g.loads.pattern("Live"):
g.loads.surface("Slab", magnitude=-2.5e3, normal=False)
Defs declared outside any pattern block belong to the implicit
"default" pattern. Downstream solvers emit one
timeSeries/pattern block per group.
Reduction & emission form¶
Three of the distributed-load factories (line, surface,
gravity, body) accept two orthogonal flags that change how the
load is converted to records:
reduction |
target_form |
Effect |
|---|---|---|
"tributary" |
"nodal" |
Default. Length/area/volume-weighted nodal lumping — one NodalLoadRecord per node |
"consistent" |
"nodal" |
Shape-function (Gauss-quadrature) integration — required for higher-order elements |
"tributary" |
"element" |
Skip nodal lumping entirely; emit one ElementLoadRecord per element |
"consistent" |
"element" |
Same as above; the solver's element handles the integration |
Use element form for beam-element line loads
(eleLoad -beamUniform), shell pressures handled inside the
element, or any solver-side load that you don't want decomposed at
the apeGmsh layer.
Target identification¶
All factory methods accept a flexible positional target argument
plus three explicit keyword overrides (pg=, label=, tag=) that
pin the lookup source. The auto path tries, in order: raw
(dim, tag) list → mesh selection → label → physical group → part
label. The first match wins. See the
LoadsComposite
class docstring for the full disambiguation rules.
Worked example¶
from apeGmsh import apeGmsh
with apeGmsh(model_name="frame") as g:
# ... geometry + Parts already imported ...
with g.loads.pattern("Dead"):
g.loads.gravity("Slab", density=2400) # body load
g.loads.line("BeamEdge", magnitude=-15e3, # distributed
direction=(0, 0, -1), # line load
reduction="tributary")
with g.loads.pattern("Push"):
g.loads.point_closest( # snaps to
xyz=(5.0, 2.5, 3.0), within="Slab", # nearest
force_xyz=(120e3, 0.0, 0.0), # mesh node
)
g.mesh.generation.generate(dim=3)
fem = g.mesh.queries.get_fem_data(dim=3)
# Pattern-by-pattern emission into OpenSees
for pat in g.loads.patterns():
ops.timeSeries("Linear", pat_tag(pat))
ops.pattern("Plain", pat_tag(pat), pat_tag(pat))
for r in fem.nodes.loads.by_pattern(pat):
ops.load(r.node_id, *(r.force_xyz or (0,)*3))
Composite¶
apeGmsh.core.LoadsComposite.LoadsComposite ¶
Loads composite — define + resolve loads.
Target resolution¶
All factory methods (point, line, surface, gravity,
body, face_load, face_sp) accept a flexible positional
target argument plus three explicit keyword overrides::
g.loads.point("my_pt", force_xyz=(0, 0, -1)) # auto
g.loads.point(pg="my_pg", force_xyz=(0, 0, -1)) # force PG
g.loads.point(label="top", force_xyz=(0, 0, -1)) # force label
g.loads.point(tag=[(0, 7)], force_xyz=(0, 0, -1)) # raw DimTag
When the caller passes target=... (the auto path),
:meth:_resolve_target tries each of these in order until one
matches:
=== ======================== =============================
Source Provided by¶
=== ======================== =============================
1 raw list[(dim, tag)] the caller
2 mesh selection name g.mesh_selection
3 label (Tier 1, prefixed) _label: physical groups
4 physical group (Tier 2) user-authored PGs
5 part label g.parts._instances
=== ======================== =============================
The first match wins. If two namespaces share a name (e.g. a label
and a PG both called "top"), label wins because it is checked
first. To bypass auto resolution and pin a specific source use the
keyword form: pg= skips straight to step 4, label= to step
3, tag= to step 1.
A KeyError is raised if auto resolution exhausts all five
sources without finding the name.
Patterns¶
All load definitions inherit the pattern of the active
:meth:pattern context (default "default"). Group loads into
named patterns so downstream solvers can emit one timeSeries /
pattern block per group.
Source code in src/apeGmsh/core/LoadsComposite.py
pattern ¶
Group subsequent load definitions under a named pattern.
Example¶
::
with g.loads.pattern("dead"):
g.loads.gravity("concrete", g=(0, 0, -9.81), density=2400)
g.loads.line("beams", magnitude=-2e3, direction="z")
with g.loads.pattern("live"):
g.loads.surface("slabs", magnitude=-3e3)
Source code in src/apeGmsh/core/LoadsComposite.py
point ¶
point(target=None, *, pg=None, label=None, tag=None, force_xyz=None, moment_xyz=None, name=None) -> PointLoadDef
Concentrated force and/or moment applied at every node of target.
Each node in the resolved target receives the same force
and moment vectors. Use this when the load point lives on a
named entity (a physical group, label, part, mesh selection,
or raw (dim, tag) list); use :meth:point_closest instead
when you only have world coordinates.
Resolution emits one
:class:~apeGmsh.solvers.Loads.NodalLoadRecord per
targeted node onto fem.nodes.loads. Both force_xyz
and moment_xyz may be supplied (or either alone), but
at least one of the two must be non-None for the load
to do anything useful.
Parameters¶
target : str or list of (dim, tag), optional
Auto-resolved positional target — see the
:class:LoadsComposite docstring for the lookup order.
Pass pg=, label=, or tag= to bypass
auto-resolution.
pg, label, tag :
Explicit-source overrides. See the class docstring.
force_xyz : (Fx, Fy, Fz), optional
Concentrated force vector applied at each targeted
node, in model force units.
moment_xyz : (Mx, My, Mz), optional
Concentrated moment vector. For 2-D models pass a
length-1 tuple (Mz,) — the resolver will accept it.
name : str, optional
Friendly name for :meth:summary and the viewer.
Returns¶
PointLoadDef
The stored definition (also appended to
self.load_defs).
Raises¶
KeyError
If target is a string that doesn't resolve to any
of label, physical group, part, or mesh selection.
ValueError
If neither target nor an explicit-source kwarg is
given.
See Also¶
point_closest : Coordinate-driven variant — snap to the nearest mesh node. face_load : Apply a centroidal force/moment to a whole face without rigidising it.
Examples¶
with g.loads.pattern("Lateral"): ... g.loads.point( ... "ColTop", ... force_xyz=(120e3, 0.0, 0.0), ... )
Source code in src/apeGmsh/core/LoadsComposite.py
point_closest ¶
point_closest(xyz, *, within=None, pg=None, label=None, tag=None, force_xyz=None, moment_xyz=None, tol=None, name=None) -> PointClosestLoadDef
Concentrated load at the mesh node closest to xyz.
Coordinate-driven targeting — useful when the load point doesn't
live on a named PG/label. The snap happens at :meth:resolve,
and the snap distance is recorded back on the def.
Parameters¶
xyz : (x, y, z)
World-coordinate target.
within : str | list, optional
Restrict the snap pool to nodes inside this PG/label/part/
DimTag list. pg=/label=/tag= force the source.
Default = global (search every mesh node).
tol : float, optional
If given, every node within tol of xyz receives the
load. Default None = single nearest node.
Source code in src/apeGmsh/core/LoadsComposite.py
line ¶
line(target=None, *, pg=None, label=None, tag=None, magnitude=None, direction=(0.0, 0.0, -1.0), q_xyz=None, normal=False, away_from=None, reduction='tributary', target_form='nodal', name=None) -> LineLoadDef
Distributed load (force per unit length) along the curve(s) of target.
Three ways to specify the load vector:
magnitude+direction: scalar magnitude along a fixed unit vector (or axis name"x"/"y"/"z").q_xyz: explicit(qx, qy, qz)force-per-length vector.normal=True+away_from: edge-by-edge in-plane pressure. The pressure direction is the edge normal lying in the plane of the loaded curves (any plane — XY, XZ, YZ, or arbitrary; the plane is fitted from the curve geometry), sign-flipped per edge so it points away fromaway_from(a reference point representing the source of the load — e.g. the centre of an arched cavity loaded by internal pressure). Positivemagnitudethen pushes into the structure.
For normal=True without away_from, apeGmsh consults
the parent surface's Gmsh-oriented normal + boundary loop to
decide which side is "into the structure" (also plane-
general). If the curve has no adjacent surface, or bounds
more than one, the resolver raises ValueError —
disambiguate by passing away_from, or fall back to
direction/q_xyz. If the loaded curves are collinear
or non-planar the in-plane normal is undefined and the
resolver raises ValueError.
Reduction and emission form¶
reduction="tributary"(default): split each edge's length-weighted load equally between its two end nodes. Emits :class:NodalLoadRecordonfem.nodes.loads.reduction="consistent": shape-function integration (line2 / line3) — equivalent to the FEM consistent load vector. Required for higher-order elements where simple tributary lumping is wrong.target_form="element": skip nodal lumping entirely and emit oneElementLoadRecordper beam element withload_type="beamUniform"— the solver's element formulation handles the integration.
Parameters¶
target : str or list of (dim, tag), optional
Curve(s) to load.
pg, label, tag :
Explicit-source overrides. See class docstring.
magnitude : float or callable, optional
Scalar force per unit length. Required if q_xyz is
None. Required when normal=True.
May also be a **callable** ``q(xyz) -> float`` that
receives the ``(x, y, z)`` coordinate as a length-3
array and returns the local force-per-length — for a
spatially varying load such as a depth-dependent ground
/ convergence pressure, e.g.
``magnitude=lambda p: gamma * (z_top - p[2])``. Works
with ``normal=True`` and with ``direction=``, in every
``target_form``; mutually exclusive with ``q_xyz``.
Accuracy of the varying field depends on ``reduction``:
* ``reduction="consistent"`` — the field is integrated
against the shape functions at the element **Gauss
points** (:func:`integrate_edge_scaled`), i.e. the
exact consistent load vector to quadrature order. No
over/undershoot of the resultant; mesh-converged even
on a coarse mesh.
* ``reduction="tributary"`` (default) — sampled once at
each edge **midpoint** and lumped (the tributary model
is itself a lumping approximation). This is the
midpoint rule: ``O(h^2)`` and exact for a linear field
on straight edges, but it can over/undershoot the true
∫q on curved edges or steep gradients — pass
``reduction="consistent"`` or refine the mesh.
direction : tuple or {"x", "y", "z"}, default (0, 0, -1)
Unit direction for magnitude. Ignored when
q_xyz or normal=True is given.
q_xyz : (qx, qy, qz), optional
Explicit force-per-length vector — overrides
magnitude × direction.
normal : bool, default False
If True, treat the load as a pressure normal to each
edge, acting in the plane of the loaded curves (fitted
from the curve geometry — any plane, not just XY).
away_from : (x, y, z), optional
Reference point for normal=True direction
disambiguation.
reduction : "tributary" or "consistent", default
"tributary"
How distributed loads are reduced to nodal records.
target_form : "nodal" or "element", default
"nodal"
Output record type. "element" skips nodal lumping
and emits eleLoad-style records.
name : str, optional
Friendly name.
Returns¶
LineLoadDef
Raises¶
ValueError
If neither magnitude nor q_xyz is supplied, or
normal=True is set without magnitude.
KeyError
If target doesn't resolve.
Examples¶
Uniform vertical line load on a beam edge::
g.loads.line(
"BeamEdge",
magnitude=-15e3,
direction=(0, 0, -1),
)
Internal pressure on a curved 2-D arch::
g.loads.line(
"InnerArc",
magnitude=p_int,
normal=True,
away_from=(0.0, 0.0, 0.0),
)
Element-form output for a beam carrying its own eleLoad
per element::
g.loads.line(
"Girder",
magnitude=-25e3,
direction=(0, 0, -1),
target_form="element",
)
Source code in src/apeGmsh/core/LoadsComposite.py
280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 | |
surface ¶
surface(target=None, *, pg=None, label=None, tag=None, magnitude=0.0, normal=True, direction=(0.0, 0.0, -1.0), reduction='tributary', target_form='nodal', name=None) -> SurfaceLoadDef
Pressure or traction on the surface(s) of target.
Two regimes selected by normal:
normal=True(default): scalar pressure normal to each face. The face normal is computed at resolution time from the mesh; positivemagnitudepushes into the face (i.e. acts opposite to the outward normal).normal=False: vector traction alongdirection, independent of face orientation.
Reduction and emission form¶
reduction="tributary"(default): split each face's area-weighted load equally among its corner nodes (tri3 / quad4 corner mass).reduction="consistent": shape-function integration via Gauss quadrature on the curved face — required for tri6, quad8, quad9. Fornormal=True, the curved normal at each Gauss point is used.target_form="element": emit oneElementLoadRecordper face withload_type="surfacePressure"and let the solver's element handle integration.
Parameters¶
target : str or list of (dim, tag)
Surface(s) to load.
pg, label, tag :
Explicit-source overrides.
magnitude : float, default 0.0
Pressure (force per unit area) when normal=True,
or traction magnitude along direction otherwise.
normal : bool, default True
True → normal pressure; False → vector
traction.
direction : (dx, dy, dz), default (0, 0, -1)
Unit traction direction. Ignored when normal=True.
reduction : "tributary" or "consistent", default
"tributary"
Lumping scheme.
target_form : "nodal" or "element", default
"nodal"
Output record type.
name : str, optional
Friendly name.
Returns¶
SurfaceLoadDef
Raises¶
KeyError
If target doesn't resolve.
Examples¶
Wind pressure on a vertical façade (positive into the face)::
g.loads.surface(
"Facade",
magnitude=1.2e3,
normal=True,
)
Vertical live load on a slab (vector traction, not pressure)::
g.loads.surface(
"Slab",
magnitude=2.5e3,
normal=False,
direction=(0, 0, -1),
)
Higher-order pressure with consistent reduction on a quad8 mesh::
g.loads.surface(
"CurvedShell",
magnitude=p,
normal=True,
reduction="consistent",
)
Source code in src/apeGmsh/core/LoadsComposite.py
445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 | |
gravity ¶
gravity(target=None, *, pg=None, label=None, tag=None, g=(0.0, 0.0, -9.81), density=None, reduction='tributary', target_form='nodal', name=None) -> GravityLoadDef
Body weight (ρ · g) over the volume(s) of target.
Convenience wrapper over :meth:body for the common case of
gravity loading. The total per-element load is
density × element_volume × g_vec, distributed to the
element's nodes.
Reduction and emission form¶
reduction="tributary"(default): split each element's weight equally among its corner nodes. Requiresdensity.reduction="consistent": for tet4 / hex8 with constant density, reduces to the same per-node share as tributary (so behaviourally equivalent today, but the path is kept separate for higher-order extensions).target_form="element": emit oneElementLoadRecordper volume element withload_type="bodyForce"carryingganddensity; the solver's element formulation handles integration.density=Noneis allowed in this form — the solver reads it from the assigned material.
Parameters¶
target : str or list of (dim, tag)
Volume(s) carrying body weight.
pg, label, tag :
Explicit-source overrides.
g : (gx, gy, gz), default (0, 0, -9.81)
Gravitational acceleration vector. Unit-sensitive —
use (0, 0, -9810) for mm models with kg-mm-s units,
etc.
density : float, optional
Material density (mass per unit volume). Required when
target_form="nodal"; optional in element form.
reduction : "tributary" or "consistent", default
"tributary"
Lumping scheme.
target_form : "nodal" or "element", default
"nodal"
Output record type.
name : str, optional
Friendly name.
Returns¶
GravityLoadDef
Raises¶
ValueError
If density is missing for target_form="nodal".
KeyError
If target doesn't resolve.
See Also¶
body : Generic per-volume body force vector.
masses.volume : Add the same density as nodal mass for
inertial response (don't double-count if the OpenSees
material already carries rho).
Examples¶
Self-weight of a concrete slab (kg-m-s, ρ = 2400 kg/m³)::
with g.loads.pattern("Dead"):
g.loads.gravity("Slab", density=2400)
Element-form gravity reading density from the material::
g.loads.gravity(
"ConcreteBlock",
target_form="element",
)
Source code in src/apeGmsh/core/LoadsComposite.py
545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 | |
body ¶
body(target=None, *, pg=None, label=None, tag=None, force_per_volume=(0.0, 0.0, 0.0), reduction='tributary', target_form='nodal', name=None) -> BodyLoadDef
Generic per-volume body force on the volume(s) of target.
General sibling of :meth:gravity — accepts an arbitrary
force-per-volume vector. The total per-element load is
force_per_volume × element_volume, distributed to the
element's nodes.
Use cases beyond gravity:
- Centrifugal / rotational body force.
- Magnetic body force in coupled-physics models.
- Thermal expansion modelled as an equivalent body force.
- Any prescribed loading proportional to volume.
Parameters¶
target : str or list of (dim, tag)
Volume(s) to load.
pg, label, tag :
Explicit-source overrides.
force_per_volume : (bx, by, bz), default (0, 0, 0)
Body force vector in force per unit volume.
reduction : "tributary" or "consistent", default
"tributary"
Lumping scheme. "consistent" falls back to
tributary for tet4/hex8 (same per-node share for
constant body force).
target_form : "nodal" or "element", default
"nodal"
Output record type. "element" emits one
ElementLoadRecord per volume element with
load_type="bodyForce" and params={"bf": ...}.
name : str, optional
Friendly name.
Returns¶
BodyLoadDef
See Also¶
gravity : Convenience wrapper for ρ · g body force.
Examples¶
Centrifugal body force ρ · ω² · r evaluated as a
constant approximation::
g.loads.body(
"Rotor",
force_per_volume=(omega**2 * rho * r_cg, 0, 0),
)
Source code in src/apeGmsh/core/LoadsComposite.py
face_load ¶
face_load(target=None, *, pg=None, label=None, tag=None, force_xyz=None, moment_xyz=None, magnitude=0.0, normal=False, direction=None, name=None) -> FaceLoadDef
Concentrated force/moment at face centroid, distributed to nodes.
force_xyz is split equally among all face nodes.
moment_xyz is converted to statically equivalent nodal
forces via least-norm distribution.
magnitude (a scalar in Newtons, not pressure) combined
with normal=True or an explicit direction produces the
equivalent total force without manually computing the face
normal. Sign convention: total = magnitude * unit_direction,
where unit_direction is +n_avg for normal=True or
the normalised direction vector otherwise. Pass
magnitude=-F for an "into-face" (pressure-like) load.
Combining magnitude with force_xyz is an error;
combining with moment_xyz is fine.
Use this instead of a reference node when you only need to apply a load to a face without structural coupling.
Parameters¶
target : str or list[(dim, tag)]
Surface to load (PG name, label, part, or raw DimTag list).
force_xyz : tuple, optional
Concentrated force (Fx, Fy, Fz) at the face centroid.
moment_xyz : tuple, optional
Concentrated moment (Mx, My, Mz) about the face centroid.
magnitude : float, default 0.0
Total scalar force in Newtons. Routed by normal/
direction to produce the equivalent force_xyz.
normal : bool, default False
When True, the area-weighted average face normal
supplies the direction; positive magnitude pushes into
the face.
direction : (dx, dy, dz), optional
Explicit unit-direction override (auto-normalised) for the
magnitude path; mutually exclusive with normal=True.
Examples¶
Symmetric pull on the two faces of an embedded crack —
normal=True resolves a per-face physical outward via
adjacent-tet centroids, so the same negative magnitude on
both coincident entities pulls each face away from its own
bonded body (opening the crack)::
with m.loads.pattern("Open"):
m.loads.face_load("Crack_normal", magnitude=-1e3, normal=True)
m.loads.face_load("Crack_inverted", magnitude=-1e3, normal=True)
Source code in src/apeGmsh/core/LoadsComposite.py
699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 | |
face_sp ¶
face_sp(target=None, *, pg=None, label=None, tag=None, dofs=None, disp_xyz=None, rot_xyz=None, magnitude=0.0, normal=False, direction=None, name=None) -> FaceSPDef
Prescribed displacement/rotation at face centroid, mapped to nodes.
Each face node receives u_i = disp_xyz + rot_xyz x r_i
where r_i is the arm from the face centroid to node i.
When disp_xyz, rot_xyz, and magnitude are all
zero / None, the result is a homogeneous fix (equivalent
to fix()).
magnitude (a scalar centroid translation, in mesh length
units) combined with normal=True or an explicit
direction produces the equivalent disp_xyz without
manually computing the face normal. Sign convention matches
:meth:face_load: total = magnitude * unit_direction,
along +n_avg for normal=True. Combining magnitude
with disp_xyz is an error; combining with rot_xyz is
fine.
Parameters¶
target : str or list[(dim, tag)]
Surface to constrain.
dofs : list[int], optional
Restraint mask (1 = constrained, 0 = free).
Defaults to [1, 1, 1].
disp_xyz : tuple, optional
Prescribed translation (ux, uy, uz) at centroid.
rot_xyz : tuple, optional
Prescribed rotation (θx, θy, θz) about centroid.
magnitude : float, default 0.0
Scalar centroid translation routed by normal/direction.
normal : bool, default False
When True, area-weighted face normal supplies the direction.
direction : (dx, dy, dz), optional
Explicit unit-direction override for the magnitude path.
Source code in src/apeGmsh/core/LoadsComposite.py
validate_pre_mesh ¶
Validate every registered load's target can be resolved.
Called by :meth:Mesh.generate before meshing so typos fail
fast instead of after minutes of meshing. Raw (dim, tag)
lists are skipped — only string targets are looked up.
Source code in src/apeGmsh/core/LoadsComposite.py
resolve ¶
resolve(node_tags, node_coords, elem_tags=None, connectivity=None, *, node_map=None, face_map=None) -> LoadSet
Resolve all stored LoadDefs into a :class:LoadSet.
Source code in src/apeGmsh/core/LoadsComposite.py
summary ¶
DataFrame of the declared load intent — one row per def.
Columns: kind, name, pattern, target, source, reduction,
target_form, params. params is a short stringified view
of the kind-specific fields (force, magnitude, direction, ...).
Source code in src/apeGmsh/core/LoadsComposite.py
Base classes¶
apeGmsh._kernel.defs.loads.LoadDef
dataclass
¶
LoadDef(kind: str, target: object, pattern: str = 'default', name: str | None = None, reduction: str = 'tributary', target_form: str = 'nodal', target_source: str = 'auto')
Base class for all load definitions.
apeGmsh._kernel.records._loads.LoadRecord
dataclass
¶
Base class for all resolved load records.
Concentrated loads¶
Concentrated forces and moments — applied either to nodes that already exist on a named target, or to the mesh node nearest a world coordinate.
apeGmsh._kernel.defs.loads.PointLoadDef
dataclass
¶
PointLoadDef(kind: str, target: object, pattern: str = 'default', name: str | None = None, reduction: str = 'tributary', target_form: str = 'nodal', target_source: str = 'auto', force_xyz: tuple[float, float, float] | None = None, moment_xyz: tuple[float, float, float] | None = None)
Bases: LoadDef
Concentrated force/moment at a node (or set of nodes).
All targeted nodes receive the same force/moment. Use
:meth:PointLoadDef.force_xyz for translational forces and
:meth:PointLoadDef.moment_xyz for moments (3D rotational DOFs).
Either may be None.
apeGmsh._kernel.defs.loads.PointClosestLoadDef
dataclass
¶
PointClosestLoadDef(kind: str, target: object, pattern: str = 'default', name: str | None = None, reduction: str = 'tributary', target_form: str = 'nodal', target_source: str = 'auto', force_xyz: tuple[float, float, float] | None = None, moment_xyz: tuple[float, float, float] | None = None, xyz_request: tuple[float, float, float] = (0.0, 0.0, 0.0), within: object | None = None, within_source: str = 'auto', tol: float | None = None, snap_distance: float | None = None)
Bases: PointLoadDef
Concentrated load at the mesh node(s) closest to a coordinate.
Coordinate-driven targeting (no PG/label required). At resolve time,
the composite snaps xyz_request to the nearest mesh node — or, if
tol is given, to every node within that radius. Pass within
(PG/label/part/DimTag list) to restrict the candidate node pool.
The actual snap distance is written back to snap_distance after
:meth:LoadsComposite.resolve, so it surfaces in summary().
Distributed loads¶
Length-, area-, or volume-distributed loads. All four accept the
reduction × target_form flags described above.
apeGmsh._kernel.defs.loads.LineLoadDef
dataclass
¶
LineLoadDef(kind: str, target: object, pattern: str = 'default', name: str | None = None, reduction: str = 'tributary', target_form: str = 'nodal', target_source: str = 'auto', magnitude: object = 0.0, direction: object = (0.0, 0.0, -1.0), q_xyz: tuple[float, float, float] | None = None, normal: bool = False, away_from: tuple[float, float, float] | None = None)
Bases: LoadDef
Distributed load along a 1-D entity (curve / beam element).
Three ways to specify the load vector:
magnitude+direction— scalar magnitude (force per unit length) and a direction vector or axis name ("x","y","z").q_xyz— explicit(qx, qy, qz)force-per-length vector.normal=True+away_from=(x0, y0, z0)— pressure perpendicular to each edge, in the plane of the loaded curves (any plane; fitted from the geometry). The in-plane normal is sign-flipped per edge so it points away fromaway_from;magnitudeis then force per unit length along that normal.
magnitude may be a constant float or a callable
q(xyz) -> float evaluated per edge midpoint (spatially varying
line load); the resolver in :mod:apeGmsh.mesh._load_resolver
handles both.
apeGmsh._kernel.defs.loads.SurfaceLoadDef
dataclass
¶
SurfaceLoadDef(kind: str, target: object, pattern: str = 'default', name: str | None = None, reduction: str = 'tributary', target_form: str = 'nodal', target_source: str = 'auto', magnitude: float = 0.0, normal: bool = True, direction: tuple[float, float, float] = (0.0, 0.0, -1.0))
Bases: LoadDef
Pressure or traction on a 2-D entity.
normal=True: scalar pressure perpendicular to the face (positive into the face).normal=False: vector traction in the given direction.
apeGmsh._kernel.defs.loads.GravityLoadDef
dataclass
¶
GravityLoadDef(kind: str, target: object, pattern: str = 'default', name: str | None = None, reduction: str = 'tributary', target_form: str = 'nodal', target_source: str = 'auto', g: tuple[float, float, float] = (0.0, 0.0, -9.81), density: float | None = None)
Bases: LoadDef
Body load from gravity = ρ·g over a volume.
If density is None, the solver bridge is expected to read
density from the assigned material/section.
apeGmsh._kernel.defs.loads.BodyLoadDef
dataclass
¶
Face load and face SP¶
Face-centroid versions used when you want to apply a centroidal force/moment or prescribed motion to a whole face without introducing a reference node and a coupling constraint.
apeGmsh._kernel.defs.loads.FaceLoadDef
dataclass
¶
FaceLoadDef(kind: str, target: object, pattern: str = 'default', name: str | None = None, reduction: str = 'tributary', target_form: str = 'nodal', target_source: str = 'auto', force_xyz: tuple[float, float, float] | None = None, moment_xyz: tuple[float, float, float] | None = None, magnitude: float = 0.0, normal: bool = False, direction: tuple[float, float, float] | None = None)
Bases: LoadDef
Concentrated force/moment at face centroid, distributed to face nodes.
force_xyz is split equally among all face nodes (F / N).
moment_xyz is converted to statically equivalent nodal forces
via a least-norm distribution such that Sum(r_i x f_i) = M and
Sum(f_i) = 0.
A scalar magnitude (total Newtons, NOT pressure) can be combined
with either normal=True or an explicit direction to produce
the equivalent force_xyz without manually computing the face
normal. Sign convention: magnitude * direction_unit always —
i.e. +magnitude with normal=True acts along the
area-weighted average outward normal +n_avg (and
-magnitude flips it, matching :class:SurfaceLoadDef's
"into-face" pressure when desired). Composes with moment_xyz;
combining with force_xyz is an error.
Use this instead of a reference node + coupling when you only need to apply a load to a face without structural coupling to another element.
apeGmsh._kernel.defs.loads.FaceSPDef
dataclass
¶
FaceSPDef(kind: str, target: object, pattern: str = 'default', name: str | None = None, reduction: str = 'tributary', target_form: str = 'nodal', target_source: str = 'auto', dofs: list[int] = (lambda: [1, 1, 1])(), disp_xyz: tuple[float, float, float] | None = None, rot_xyz: tuple[float, float, float] | None = None, magnitude: float = 0.0, normal: bool = False, direction: tuple[float, float, float] | None = None)
Bases: LoadDef
Prescribed displacement/rotation at face centroid, mapped to face nodes.
Maps a rigid-body motion at the face centroid to per-node
displacements using u_i = disp_xyz + rot_xyz x r_i.
When disp_xyz, rot_xyz, and magnitude are all None /
zero, the result is a homogeneous fix.
A scalar magnitude (displacement, in mesh length units) can be
combined with normal=True or an explicit direction to
derive the centroid translation without computing the face normal
by hand. Sign convention matches :class:FaceLoadDef: total =
magnitude * unit_direction along +n_avg (or the normalised
direction). Composes with rot_xyz; combining with
disp_xyz is an error.
Parameters¶
dofs : list[int]
Restraint mask — 1 for constrained DOFs, 0 for free.
disp_xyz : tuple or None
Prescribed translation at the face centroid.
rot_xyz : tuple or None
Prescribed rotation about the face centroid.
magnitude : float
Scalar centroid translation, routed via normal/direction.
normal : bool
When True, use the area-weighted face normal as the direction.
direction : tuple or None
Explicit unit direction (auto-normalised); mutually exclusive
with normal=True.
Resolved records¶
What ends up on the FEM broker after meshing.
apeGmsh._kernel.records._loads.NodalLoadRecord
dataclass
¶
NodalLoadRecord(kind: str, pattern: str = 'default', name: str | None = None, node_id: int = 0, force_xyz: tuple[float, float, float] | None = None, moment_xyz: tuple[float, float, float] | None = None)
Bases: LoadRecord
Force and/or moment at a single node.
force_xyz and moment_xyz are pure 3D spatial vectors (or
None when absent). The record is DOF-agnostic — mapping onto
a solver's DOF space is the caller's responsibility.
apeGmsh._kernel.records._loads.ElementLoadRecord
dataclass
¶
apeGmsh._kernel.records._loads.SPRecord
dataclass
¶
SPRecord(kind: str, pattern: str = 'default', name: str | None = None, node_id: int = 0, dof: int = 1, value: float = 0.0, is_homogeneous: bool = True)
Bases: LoadRecord
Single-point constraint: prescribed displacement or homogeneous fix.
One record per DOF per node. When is_homogeneous is True
the downstream emitter can use ops.fix(); otherwise it must use
ops.sp(node, dof, value).
Resolver¶
apeGmsh._kernel.resolvers._load_resolver.LoadResolver ¶
LoadResolver(node_tags: ndarray, node_coords: ndarray, elem_tags: ndarray | None = None, connectivity: ndarray | None = None)
Convert :class:LoadDef instances to :class:LoadRecord lists.
Pure mesh math — receives raw arrays and a DOF context, returns record lists. No Gmsh queries (the composite handles target resolution before calling here).
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
face_area ¶
Polygonal face area via fan triangulation from node[0].
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
face_normal ¶
Outward normal estimate from the first three nodes.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
element_volume ¶
Approximate element volume from its node coordinates.
Tet4: actual volume. Hex8: approximate via 6 tets. For other element types, returns the convex hull bbox volume as a coarse fallback.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_point ¶
Apply the same force/moment to every node in node_set.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_line_tributary ¶
Distribute a line load by length-weighted nodal share.
edges is a list of (node_a, node_b) pairs covering the loaded
curve. Each node receives magnitude * Σ(adjacent_edge_len/2)
in the direction vector.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_line_per_edge_tributary ¶
resolve_line_per_edge_tributary(defn: LineLoadDef, items: list[tuple[int, int, ndarray]]) -> list[NodalLoadRecord]
Tributary line-load reduction with a per-edge force-per-length.
items is a list of (n1, n2, q_xyz) triples; each q_xyz
is the force-per-length vector applied to that single edge.
Used by the composite for normal=True loads where the
direction varies edge-by-edge.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_surface_tributary ¶
resolve_surface_tributary(defn: SurfaceLoadDef, faces: list[list[int]], outwards: list[ndarray] | None = None) -> list[NodalLoadRecord]
Distribute a surface load by tributary area.
faces is a list of node-id lists (one per face element).
For each face, total = magnitude * area and is split
equally among the face's nodes. normal=True projects
along the face normal; otherwise the explicit direction vector.
outwards, when given, supplies a per-face physical
outward unit normal that overrides the connectivity-derived
:meth:face_normal. This is needed for embedded crack
faces (whose connectivity normal can disagree with physical
outward) and for tilted faces with unpredictable connectivity
orientation; the composite layer fills it in via
:meth:LoadsComposite._face_outward_normals. When None,
the connectivity normal is used (preserving backward compat
for direct callers that don't go through the composite).
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_gravity_tributary ¶
Distribute body weight equally to element nodes.
elements is a list of connectivity rows (each is an array
of node IDs). Each element contributes ρ·V·g total,
split equally among its nodes.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_body_tributary ¶
Distribute a per-volume body force equally to element nodes.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_line_consistent ¶
Consistent line-load reduction via shape-function integration.
edges is a list of node-id sequences; each sequence's length determines element order:
2 nodes -> line2 (linear), 2-pt Gauss
3 nodes -> line3 (quadratic), 3-pt Gauss
Any other node count raises :class:NotImplementedError rather
than silently producing wrong numbers.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_line_per_edge_consistent ¶
resolve_line_per_edge_consistent(defn: LineLoadDef, items: list[tuple[list[int], ndarray]]) -> list[NodalLoadRecord]
Consistent line-load reduction with a per-edge force-per-length.
items is a list of (node_seq, q_xyz) pairs; each
q_xyz is treated as constant along that edge. Shape-
function integration is otherwise identical to
:meth:resolve_line_consistent.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_line_per_edge_consistent_varying ¶
Consistent reduction of a spatially varying line load.
items is a list of (node_seq, dir_vec, scalar_fn) tuples.
For each edge the scalar force-per-length scalar_fn(xyz) is
integrated against the shape functions at the element's Gauss
points (:func:integrate_edge_scaled) and applied along the
per-edge dir_vec (the in-plane normal for normal=True,
or the direction vector otherwise). This is the exact
consistent load vector to quadrature order, so a varying
magnitude does not over/undershoot the way a single midpoint
sample does.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_surface_consistent ¶
Consistent surface load via shape-function integration.
Each face is a node-id sequence whose length determines the face type:
3 -> tri3, 4 -> quad4, 6 -> tri6, 8 -> quad8, 9 -> quad9
For defn.normal=True the pressure follows the curved face
normal evaluated at each Gauss point. Any other node count
raises :class:NotImplementedError.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_gravity_consistent ¶
Consistent gravity reduction.
For tet4 / hex8 with constant density, the consistent vector equals the tributary vector (each node gets V/n × ρ × g).
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_line_element ¶
Emit one ElementLoadRecord per beam element with beamUniform params.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_line_element_varying ¶
resolve_line_element_varying(defn: LineLoadDef, items: list[tuple[int, tuple[float, float, float]]]) -> list[ElementLoadRecord]
Per-element beamUniform for a spatially varying line load.
items is a list of (element_id, (wx, wy, wz)) pairs — the
composite has already sampled the callable magnitude at each
element's midpoint, so each element gets its own constant
beamUniform (the only thing OpenSees eleLoad supports).
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_surface_element ¶
Emit one ElementLoadRecord per face element with surfacePressure.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_gravity_element ¶
Emit one ElementLoadRecord per volume element with bodyForce.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_body_element ¶
Emit one ElementLoadRecord per volume element with bodyForce.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
resolve_face_load ¶
resolve_face_load(defn: FaceLoadDef, face_node_ids: list[int], faces: list[list[int]] | None = None, outwards: list[ndarray] | None = None) -> list[NodalLoadRecord]
Distribute centroidal force/moment to face nodes.
force_xyz: equal share F / N per node.
moment_xyz: least-norm nodal forces satisfying
Sum(f_i) = 0 and Sum(r_i x f_i) = M.
magnitude + normal/direction: equivalent force_xyz
derived from face geometry. Requires faces (per-element
node-id lists) when normal=True so the area-weighted
average normal can be computed.
outwards, when given, supplies a per-face physical outward
unit normal that overrides the connectivity-derived
:meth:face_normal in the area-weighted average. See the
outwards= discussion on
:meth:resolve_surface_tributary.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 | |
resolve_face_sp ¶
resolve_face_sp(defn: FaceSPDef, face_node_ids: list[int], faces: list[list[int]] | None = None, outwards: list[ndarray] | None = None) -> list[SPRecord]
Map centroidal rigid-body motion to per-node SP constraints.
For each constrained DOF d and each node i:
u_i = disp_xyz + rot_xyz x r_i, then emit
SPRecord(node_id=i, dof=d, value=u_i[d-1]).
magnitude + normal/direction derive an additional
translation contribution (along +n_avg for normal=True,
otherwise along the normalised direction). Requires
faces when normal=True.
outwards, when given, supplies a per-face physical outward
unit normal that overrides the connectivity-derived
:meth:face_normal. See
:meth:resolve_surface_tributary.
Source code in src/apeGmsh/_kernel/resolvers/_load_resolver.py
786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 | |