Mesh#
Meshes#
- fluxfem.mesh.BaseMesh#
alias of
BaseMeshClosure
- class fluxfem.mesh.BaseMeshPytree(coords: jax.numpy.ndarray, conn: jax.numpy.ndarray, cell_tags: jax.numpy.ndarray | None = None, node_tags: jax.numpy.ndarray | None = None)#
BaseMesh variant that registers as a JAX pytree.
- class fluxfem.mesh.HexMesh(coords: jax.numpy.ndarray, conn: jax.numpy.ndarray, cell_tags: jax.numpy.ndarray | None = None, node_tags: jax.numpy.ndarray | None = None)#
Structured / unstructured hex mesh (8-node linear hex elements).
coords: (n_nodes, 3) float32 conn: (n_elems, 8) int64 # node indices of each element
- class fluxfem.mesh.HexMeshPytree(*args: Any, **kwargs: Any)#
Hex mesh registered as a JAX pytree.
- class fluxfem.mesh.StructuredHexBox(nx: int, ny: int, nz: int, lx: float = 1.0, ly: float = 1.0, lz: float = 1.0, origin: tuple[float, float, float] = (0.0, 0.0, 0.0), order: int = 1)#
Uniform hex mesh generator on a box, returned as an unstructured HexMesh.
- fluxfem.mesh.tag_axis_minmax_facets(mesh: HexMesh, axis: int = 0, dirichlet_tag: int = 1, neumann_tag: int = 2, tol: float = 1e-08) Tuple[jax.numpy.ndarray, jax.numpy.ndarray]#
Tag boundary facets on min/max of the given axis.
- Returns:
(n_facets, 4) int64, quad node ids facet_tags: (n_facets,) int64, dirichlet_tag on min side, neumann_tag on max side
- Return type:
facets
- class fluxfem.mesh.TetMesh(coords: jax.numpy.ndarray, conn: jax.numpy.ndarray, cell_tags: jax.numpy.ndarray | None = None, node_tags: jax.numpy.ndarray | None = None)#
Unstructured tetra mesh.
- class fluxfem.mesh.TetMeshPytree(*args: Any, **kwargs: Any)#
Unstructured tetra mesh (pytree).
- class fluxfem.mesh.StructuredTetBox(nx: int, ny: int, nz: int, lx: float = 1.0, ly: float = 1.0, lz: float = 1.0, origin: tuple[float, float, float] = (0.0, 0.0, 0.0), order: int = 1)#
Regular grid subdivided into 5 tets per cube (simplex) in a structured layout.
- class fluxfem.mesh.StructuredTetTensorBox(nx: int, ny: int, nz: int, lx: float = 1.0, ly: float = 1.0, lz: float = 1.0, origin: tuple[float, float, float] = (0.0, 0.0, 0.0), order: int = 1)#
Regular grid subdivided into 6 tets per cube (matches skfem MeshTet.init_tensor).
Surface#
- class fluxfem.mesh.SurfaceMesh(coords: jax.numpy.ndarray, conn: jax.numpy.ndarray, cell_tags: jax.numpy.ndarray | None = None, node_tags: jax.numpy.ndarray | None = None, facet_tags: jax.numpy.ndarray | None = None)#
Simple boundary mesh made of facets (tri/quad) that live in the volume mesh nodes. Uses BaseMesh.conn to store facets.
- class fluxfem.mesh.SurfaceMeshPytree(*args: Any, **kwargs: Any)#
Simple boundary mesh made of facets (tri/quad) that live in the volume mesh nodes. Uses BaseMesh.conn to store facets.
IO#
- fluxfem.mesh.load_gmsh_mesh(path: str)#
Load a Gmsh .msh (v2/v4) file containing hex or tet elements (and optional boundary facets).
- Returns:
HexMesh or TetMesh facets: (n_facets, 3 or 4) or None facet_tags: (n_facets,) or None (gmsh physical tags if present)
- Return type:
mesh
- fluxfem.mesh.load_gmsh_hex_mesh(path: str)#
Load a Gmsh mesh and return a HexMesh with optional facets/tags.
- fluxfem.mesh.load_gmsh_tet_mesh(path: str)#
Load a Gmsh mesh and return a TetMesh with optional facets/tags.
- fluxfem.mesh.make_surface_from_facets(coords: ndarray, facets: ndarray, tags=None) SurfaceMesh#
Helper to build SurfaceMesh from raw coords/facets (optional tags).
Predicates#
- fluxfem.mesh.bbox_predicate(mins: ndarray, maxs: ndarray, tol: float = 1e-08)#
Return predicate True for nodes on the axis-aligned bounding box.
- fluxfem.mesh.plane_predicate(axis: int, value: float, tol: float = 1e-08)#
Return predicate True when all nodes lie on plane x[axis]=value (within tol).
- fluxfem.mesh.axis_plane_predicate(axis: int, value: float, tol: float = 1e-08)#
Alias of plane_predicate for readability.
- fluxfem.mesh.slab_predicate(axis: int, min_val: float, max_val: float, tol: float = 1e-08)#
Return predicate True when nodes lie in a slab min<=x[axis]<=max (within tol).
Contact#
For new public-facing contact setup, prefer ContactSpaces,
ContactGroupSpaces, and OneSidedContactSpaces. The concrete contact
space classes remain available as lower-level constructors.
- class fluxfem.mesh.ContactSide(surface: 'SurfaceMesh', elem_conn: 'np.ndarray | None', value_dim: 'int', space: 'object | None' = None)#
- class fluxfem.mesh.ContactSpaces(master: ContactSide, slave: ContactSide, field_master: str = 'a', field_slave: str = 'b')
Public spec that binds contact roles to contact sides.
- class fluxfem.mesh.ContactGroupSpaces(master: ContactSide, slaves: Sequence[ContactSide], field_master: str = 'master', field_slave: str = 'slave')
Public spec that binds one-master/many-slave contact roles.
- class fluxfem.mesh.OneSidedContactSpaces(side: ContactSide, surface_master: SurfaceMesh | None = None, elem_conn_master: ndarray | None = None, facet_to_elem_master: ndarray | None = None)
Public spec that binds one-sided contact roles to a contact side.
- class fluxfem.mesh.ContactSurfaceSpace(surface_master: SurfaceMesh, surface_slave: SurfaceMesh, supermesh_coords: ndarray, supermesh_conn: ndarray, source_facets_master: ndarray, source_facets_slave: ndarray, elem_conn_master: ndarray | None, elem_conn_slave: ndarray | None, facet_to_elem_master: ndarray | None, facet_to_elem_slave: ndarray | None, field_master: str = 'a', field_slave: str = 'b', value_dim_master: int = 1, value_dim_slave: int = 1, space_mode_master: str = 'nodal', space_mode_slave: str = 'nodal', facet_dofs_master: ndarray | None = None, facet_dofs_slave: ndarray | None = None, trial_value_dim_master: int | None = None, trial_value_dim_slave: int | None = None, trial_space_mode_master: str | None = None, trial_space_mode_slave: str | None = None, trial_facet_dofs_master: ndarray | None = None, trial_facet_dofs_slave: ndarray | None = None, quad_order: int = 0, normal_sign: float | None = None, tol: float = 1e-08, backend: str = 'jax', batch_jac: bool | None = None, supermesh_quad_cache: Any | None = None)#
Surface interface wrapper for contact assembly on a supermesh.
- class fluxfem.mesh.OneToManyContactSurfaceSpace(contacts: tuple[ContactSurfaceSpace, ...], field_master: str = 'master', field_slave: str = 'slave')#
One-master/multi-slave wrapper built from pairwise ContactSurfaceSpace objects.
- class fluxfem.mesh.OneSidedContactSurfaceSpace(surface_slave: SurfaceMesh, elem_conn_slave: ndarray, facet_to_elem_slave: ndarray, value_dim: int = 1, quad_order: int = 2, normal_sign: float = 1.0, tol: float = 1e-08, surface_master: SurfaceMesh | None = None, elem_conn_master: ndarray | None = None, facet_to_elem_master: ndarray | None = None)#
Surface wrapper for one-sided (Dirichlet) contact assembly.
- class fluxfem.mesh.ContactOperators(enforcement: str, law: str | None = None, formulation: str | None = None, coupling_aa: Any | None = None, coupling_ab: Any | None = None, B_a: Any | None = None, B_b: Any | None = None, B: Any | None = None, Kuu: Any | None = None, residual: Any | None = None, jacobian: Any | None = None, facet_conn_master: ndarray | None = None, rho: float | None = None, multiplier: Any | None = None)#
Container for assembled contact operators.
- fluxfem.mesh.assemble_contact_constraint_operators(contact, *, law: str | None = None, formulation: str | None = None, rho: float = 0.0, multiplier: ContactMultiplierSpace, backend: str | None = None, weak_form: MixedSurfaceResidualForm | None = None, state: Mapping[str, npt.ArrayLike] | Sequence[Any] | None = None, res_form: MixedSurfaceResidualForm | None = None, u: Mapping[str, npt.ArrayLike] | Sequence[Any] | None = None, params: 'WeakParams' | None = None, normal_source: str = 'master', sparse: bool = False, batch_jac: bool | None = None) ContactOperators#
- fluxfem.mesh.assemble_contact_penalty_operators(contact, *, law: str | None = None, formulation: str | None = None, backend: str | None = None, weak_form: MixedSurfaceResidualForm | None = None, state: Mapping[str, npt.ArrayLike] | Sequence[Any] | None = None, res_form: MixedSurfaceResidualForm | None = None, u: Mapping[str, npt.ArrayLike] | Sequence[Any] | None = None, params: 'WeakParams' | None = None, normal_source: str = 'master', sparse: bool = False, batch_jac: bool | None = None) ContactOperators#
- fluxfem.mesh.assemble_contact_interface_residual(*args, **kwargs)#
Assemble residual on a contact interface supermesh.
- fluxfem.mesh.assemble_contact_interface_jacobian(*args, **kwargs)#
Assemble Jacobian on a contact interface supermesh.
- fluxfem.mesh.assemble_contact_coupling_matrices(*args, **kwargs)#
Assemble coupling matrices for contact interface constraints.
- fluxfem.mesh.assemble_contact_kkt(coupling_aa, coupling_ab, *, rho: float = 0.0, multiplier: ContactMultiplierSpace, facet_conn_master: ndarray | None = None, backend: str | None = None, format: str = 'fluxsparse', return_blocks: bool = False)#
- fluxfem.mesh.solve_contact_kkt(kkt_matrix, rhs, *, backend: str | None = None, diagonal_shift: float = 0.0, config: ContactKKTSolveConfig | None = None)#
Solve KKT linear system
KKT * x = rhs.config is the preferred control surface. backend=None auto-selects from
kkt_matrix/rhswhen no explicit config is provided.
- fluxfem.mesh.facet_gap_values(coords: ndarray, facets: ndarray, u: ndarray, n: ndarray, c: float, *, value_dim: int | None = None, reduce: str = 'min') tuple[ndarray, float]#
Compute per-facet gap values for a one-sided contact plane.
Returns (g_f, min_g_all) where g_f is reduced per facet and min_g_all is the global minimum node gap.
- fluxfem.mesh.active_contact_facets(coords: ndarray, facets: ndarray, u: ndarray, n: ndarray, c: float, *, value_dim: int | None = None, reduce: str = 'min', threshold: float = 0.0) tuple[ndarray, float]#
Return active facet indices and global minimum gap for one-sided contact.
Remote Constraints#
- fluxfem.mesh.assemble_rbe2_constraint_matrix(ref_point: ndarray, slave_coords: ndarray, *, backend: str | None = None)#
Assemble 3D RBE2-style rigid kinematic constraints.
- Unknown ordering:
q = [u_ref(3), omega_ref(3), u_slave_0(3), …, u_slave_{n-1}(3)]
- Constraint for each slave node i:
u_slave_i - u_ref - (omega_ref x (x_i - x_ref)) = 0
- fluxfem.mesh.assemble_rbe3_constraint_matrix(ref_point: ndarray, slave_coords: ndarray, *, weights: ndarray | None = None, normalize_weights: bool = True, backend: str | None = None)#
Assemble a weighted 3D RBE3-style interpolation constraint.
- Unknown ordering:
q = [u_ref(3), omega_ref(3), u_slave_0(3), …, u_slave_{n-1}(3)]
The constraints are formed from weighted rigid-body reconstruction in normal- equation form:
(sum_i w_i B_i^T B_i) q_ref - sum_i w_i B_i^T u_i = 0
where
B_i = [I, -[r_i]_x]andr_i = x_i - x_ref.This yields a 6 x (6 + 3*n_slave) matrix. Repeated use of this helper allows multiple user-defined RBE3 couplings to be added to one system.
- fluxfem.mesh.build_rbe3_weights(ref_point: ndarray, slave_coords: ndarray, *, method: str = 'equal', surface: SurfaceMesh | None = None, power: float = 2.0, eps: float = 1e-12, normalize: bool = True) ndarray#
Build convenience weights for RBE3-style interpolation.
Supported methods#
equal:Uniform node weights.
distance:Inverse-distance^power weights from the remote point.
facet_area:Lump each facet area equally to its nodes, then normalize per node. Requires
surfacewhose node numbering matchesslave_coordsorder.
- fluxfem.mesh.build_rbe3_remote_resultant(ref_point: ndarray, slave_coords: ndarray, *, surface: SurfaceMesh, load: ArrayLike | None = None, pressure: float | ArrayLike | None = None, outward_from: ArrayLike | None = None) ndarray#
Build the equivalent remote-point resultant for an RBE3-supported surface.
The returned 6-vector is ordered as
[force(3), moment(3)]and is compatible with a 6-DOF remote-point field ordered as[u_ref(3), omega_ref(3)].Exactly one of
loadorpressuremust be provided:load: constant vector load per unit area with shape(3,)or(n_facets, 3)pressure: scalar normal traction with shape()or(n_facets,)