Skip to content

compas_cgal.booleans ¤

Functions¤

boolean_chain ¤

boolean_chain(
    meshes: Iterable[VerticesFaces],
    operations: Iterable[Literal["union", "difference", "intersection", "xor"]],
    hybrid: bool = False,
) -> VerticesFacesNumpy

Run a chain of boolean operations entirely in C++ without round-tripping intermediates.

Computes result = meshes[0]; result = result OP_i meshes[i+1] for each operation in order. The full collection of input meshes is sent to C++ in a single call; intermediate meshes never leave C++; only the final (V, F) is returned to Python.

Parameters:

Name Type Description Default
meshes iterable of (V, F)

Triangle meshes. Length must be len(operations) + 1.

required
operations iterable of {"union", "difference", "intersection", "xor"}

Per-step operation. "difference" is result - meshes[i+1]. "xor" is the symmetric difference.

required
hybrid bool

If False, runs the chain in CGAL's Exact_predicates_exact_constructions_kernel (EPECK). Constructions are lazy-exact, so the chain handles geometrically degenerate input — e.g., three cylinders meeting at the origin — without any shifts.

If True, runs on an EPICK Surface_mesh augmented with an EPECK vertex property map passed to corefinement via parameters::vertex_point_map() — the pattern from CGAL's "consecutive boolean operations with exact point maps" example. Storage stays at native double precision while every intersection vertex is constructed exactly: same robustness as full EPECK without lazy-exact storage cost.

False

Returns:

Type Description
(V, F) : VerticesFacesNumpy

The final mesh.

boolean_chain_with_face_source ¤

boolean_chain_with_face_source(
    meshes: Iterable[VerticesFaces],
    operations: Iterable[Literal["union", "difference", "intersection"]],
    hybrid: bool = False,
) -> VerticesFacesSourceNumpy

Boolean chain that also tracks, for every output face, which input mesh and which input face produced it. Returns (V, F, S) where S[i] = [mesh_id, face_id]. mesh_id is the position of the source mesh in meshes and face_id is its row in that mesh's input face array.

Tracking is done via a CGAL corefinement visitor that propagates per-face tags through subface creations and face copies — the technique used in the Cockroach project for CGAL face-color tracking through booleans.

See :func:boolean_chain for the meaning of the hybrid flag.

"xor" is not supported here.

boolean_difference_mesh_mesh ¤

boolean_difference_mesh_mesh(A: VerticesFaces, B: VerticesFaces) -> VerticesFacesNumpy

Boolean difference of two meshes.

Parameters:

Name Type Description Default
A VerticesFaces

Mesh A.

required
B VerticesFaces

Mesh B.

required

Returns:

Type Description
VerticesFacesNumpy

Examples:

>>> from compas.geometry import Box, Sphere, Polyhedron
>>> from compas_cgal.booleans import boolean_difference_mesh_mesh
>>> box = Box(1)
>>> sphere = Sphere(0.5, point=[1, 1, 1])
>>> A = box.to_vertices_and_faces(triangulated=True)
>>> B = sphere.to_vertices_and_faces(u=32, v=32, triangulated=True)
>>> C = boolean_difference_mesh_mesh(A, B)
>>> shape = Polyhedron(*C)

boolean_difference_mesh_mesh_with_edges ¤

boolean_difference_mesh_mesh_with_edges(
    A: VerticesFaces, B: VerticesFaces
) -> VerticesFacesEdgesNumpy

Boolean difference returning (V, F, E). See boolean_union_mesh_mesh_with_edges.

boolean_difference_mesh_mesh_with_face_source ¤

boolean_difference_mesh_mesh_with_face_source(
    A: VerticesFaces, B: VerticesFaces
) -> VerticesFacesSourceNumpy

Boolean difference returning (V, F, S). See boolean_union_mesh_mesh_with_face_source.

boolean_difference_mesh_meshes ¤

boolean_difference_mesh_meshes(
    A: VerticesFaces, Bs: Iterable[VerticesFaces]
) -> VerticesFacesNumpy

Subtract many meshes from A in a single corefinement.

Sequential A = A - B_i chains accumulate subdivision and round-off, which can crash CGAL's corefinement on long sequences (see CGAL issue #9282). Merging the cutters into one disjoint operand and doing a single difference avoids that.

The cutters in Bs should not overlap each other; if they do, merge them first via boolean_union_mesh_mesh and pass the union as a single operand.

boolean_intersection_mesh_mesh ¤

boolean_intersection_mesh_mesh(
    A: VerticesFaces, B: VerticesFaces
) -> VerticesFacesNumpy

Boolean intersection of two meshes.

Parameters:

Name Type Description Default
A VerticesFaces

Mesh A.

required
B VerticesFaces

Mesh B.

required

Returns:

Type Description
VerticesFacesNumpy

Examples:

>>> from compas.geometry import Box, Sphere, Polyhedron
>>> from compas_cgal.booleans import boolean_intersection_mesh_mesh
>>> box = Box(1)
>>> sphere = Sphere(0.5, point=[1, 1, 1])
>>> A = box.to_vertices_and_faces(triangulated=True)
>>> B = sphere.to_vertices_and_faces(u=32, v=32, triangulated=True)
>>> C = boolean_intersection_mesh_mesh(A, B)
>>> shape = Polyhedron(*C)

boolean_intersection_mesh_mesh_with_edges ¤

boolean_intersection_mesh_mesh_with_edges(
    A: VerticesFaces, B: VerticesFaces
) -> VerticesFacesEdgesNumpy

Boolean intersection returning (V, F, E). See boolean_union_mesh_mesh_with_edges.

boolean_intersection_mesh_mesh_with_face_source ¤

boolean_intersection_mesh_mesh_with_face_source(
    A: VerticesFaces, B: VerticesFaces
) -> VerticesFacesSourceNumpy

Boolean intersection returning (V, F, S). See boolean_union_mesh_mesh_with_face_source.

boolean_union_mesh_mesh ¤

boolean_union_mesh_mesh(A: VerticesFaces, B: VerticesFaces) -> VerticesFacesNumpy

Boolean union of two meshes.

Parameters:

Name Type Description Default
A VerticesFaces

Mesh A.

required
B VerticesFaces

Mesh B.

required

Returns:

Type Description
VerticesFacesNumpy

Examples:

>>> from compas.geometry import Box, Sphere, Polyhedron
>>> from compas_cgal.booleans import boolean_union_mesh_mesh
>>> box = Box(1)
>>> sphere = Sphere(0.5, point=[1, 1, 1])
>>> A = box.to_vertices_and_faces(triangulated=True)
>>> B = sphere.to_vertices_and_faces(u=32, v=32, triangulated=True)
>>> C = boolean_union_mesh_mesh(A, B)
>>> shape = Polyhedron(*C)

boolean_union_mesh_mesh_with_edges ¤

boolean_union_mesh_mesh_with_edges(
    A: VerticesFaces, B: VerticesFaces
) -> VerticesFacesEdgesNumpy

Boolean union returning (V, F, E) where E lists vertex-index pairs of intersection-curve edges.

The edges are the corefinement intersection curve in the output mesh and can be used to recover clean polygonal outlines / face groupings without remeshing. See CGAL discussion #8711.

boolean_union_mesh_mesh_with_face_source ¤

boolean_union_mesh_mesh_with_face_source(
    A: VerticesFaces, B: VerticesFaces
) -> VerticesFacesSourceNumpy

Boolean union returning (V, F, S) where S[i] = [mesh_id, face_id] of the input face that produced output face i.

mesh_id is 0 for A and 1 for B; face_id is the original face index into FA / FB. Tracking is done through a CGAL corefinement visitor so face splits during the boolean are correctly propagated to all sub-faces.

split_by_source ¤

split_by_source(
    V: ndarray, F: ndarray, S: ndarray
) -> dict[int, tuple[ndarray, ndarray]]

Split a face-source-tagged boolean result into one mesh per source.

Given (V, F, S) from :func:boolean_chain_with_face_source, returns a dict mapping each mesh_id present in S[:, 0] to its own (V_sub, F_sub) pair. V_sub contains only the vertices referenced by that submesh's faces and F_sub is reindexed accordingly.

Note: vertices on the boundary between two source regions are duplicated across the resulting submeshes (each submesh gets its own copy), so the submeshes are no longer connected at cut boundaries. The original (V, F, S) remains the canonical output — use this helper only when a per-source mesh layout is convenient (e.g. assigning a single colour or material to a viewer scene object).

split_mesh_mesh ¤

split_mesh_mesh(A: VerticesFaces, B: VerticesFaces) -> VerticesFacesNumpy

Split one mesh with another.

Parameters:

Name Type Description Default
A VerticesFaces

Mesh A.

required
B VerticesFaces

Mesh B.

required

Returns:

Type Description
VerticesFacesNumpy

Examples:

>>> from compas.geometry import Box, Sphere, Polyhedron
>>> from compas_cgal.booleans import split_mesh_mesh
>>> box = Box(1)
>>> sphere = Sphere(0.5, point=[1, 1, 1])
>>> A = box.to_vertices_and_faces(triangulated=True)
>>> B = sphere.to_vertices_and_faces(u=32, v=32, triangulated=True)
>>> V, F = split_mesh_mesh(A, B)
>>> shape = Polyhedron(V.tolist(), F.tolist())