Skip to content

Commit 35b7383

Browse files
perf(multires meshes): vastly faster mesh chunking (#214)
* wip: try to use new zmesh.chunk_mesh * wip: meshes appearing, but of low quality * fix: use a more reasonable radius for merging close vertices * fix: grid_origin was incorrect for these meshes, use minpt instead * fix: use grid offsets * fix: remove unused parameter * install: set minimum zmesh version as 1.9.0
1 parent 82a7a80 commit 35b7383

File tree

2 files changed

+27
-62
lines changed

2 files changed

+27
-62
lines changed

igneous/tasks/mesh/multires.py

Lines changed: 26 additions & 61 deletions
Original file line numberDiff line numberDiff line change
@@ -109,7 +109,11 @@ def process_mesh(
109109
return (None, None)
110110

111111
lods = [
112-
create_octree_level_from_mesh(lods[lod], chunk_shape, lod, len(lods), grid_origin, mesh_shape)
112+
create_octree_level_from_mesh(
113+
lods[lod], chunk_shape,
114+
lod, len(lods),
115+
grid_origin, mesh_shape
116+
)
113117
for lod in range(len(lods))
114118
]
115119
fragment_positions = [ nodes for submeshes, nodes in lods ]
@@ -132,6 +136,9 @@ def process_mesh(
132136
mesh_binaries = []
133137
for lod, submeshes in enumerate(lods):
134138
for frag_no, submesh in enumerate(submeshes):
139+
if submesh.empty():
140+
continue
141+
135142
submesh.vertices = to_stored_model_space(
136143
submesh.vertices, manifest,
137144
lod=lod,
@@ -222,7 +229,8 @@ def MultiResShardedMeshMergeTask(
222229
# important to iterate this way to avoid
223230
# creating a copy of meshes vs. { ... for in }
224231
for label in labels:
225-
meshes[label] = Mesh.concatenate(*meshes[label])
232+
mesh = Mesh.concatenate(*meshes[label], segid=label)
233+
meshes[label] = mesh.consolidate()
226234
del labels
227235

228236
fname, shard = create_mesh_shard(
@@ -324,7 +332,7 @@ def generate_lods(
324332
)
325333

326334
lods.append(
327-
Mesh(*simplifier.getMesh())
335+
Mesh(*simplifier.getMesh(), segid=label)
328336
)
329337

330338
return lods
@@ -493,7 +501,7 @@ def less_msb(x: int, y: int) -> bool:
493501
return lhs[msd] - rhs[msd]
494502

495503
def determine_mesh_shape_from_lods(
496-
lods: list[trimesh.Trimesh],
504+
lods: list[Mesh],
497505
) -> Tuple[npt.NDArray[np.float32], npt.NDArray[np.int_]]:
498506
mesh_starts = [np.min(lod.vertices, axis=0) for lod in lods]
499507
mesh_ends = [np.max(lod.vertices, axis=0) for lod in lods]
@@ -503,54 +511,17 @@ def determine_mesh_shape_from_lods(
503511

504512
return grid_origin, mesh_shape
505513

506-
def generate_gridded_submeshes(
507-
mesh: trimesh.Trimesh,
508-
offset: np.ndarray,
509-
grid_size: Vec,
510-
scale: Vec,
511-
) -> Iterator[Tuple[trimesh.Trimesh, Tuple[int, int, int]]]:
512-
nx, ny, nz = np.eye(3)
513-
ox, oy, oz = offset * np.eye(3)
514-
515-
for x in range(0, grid_size.x):
516-
# list(...) required b/c it doesn't like Vec classes
517-
mesh_x = trimesh.intersections.slice_mesh_plane(mesh, plane_normal=nx, plane_origin=list(nx*x*scale.x+ox))
518-
mesh_x = trimesh.intersections.slice_mesh_plane(mesh_x, plane_normal=-nx, plane_origin=list(nx*(x+1)*scale.x+ox))
519-
for y in range(0, grid_size.y):
520-
mesh_y = trimesh.intersections.slice_mesh_plane(mesh_x, plane_normal=ny, plane_origin=list(ny*y*scale.y+oy))
521-
mesh_y = trimesh.intersections.slice_mesh_plane(mesh_y, plane_normal=-ny, plane_origin=list(ny*(y+1)*scale.y+oy))
522-
for z in range(0, grid_size.z):
523-
mesh_z = trimesh.intersections.slice_mesh_plane(mesh_y, plane_normal=nz, plane_origin=list(nz*z*scale.z+oz))
524-
mesh_z = trimesh.intersections.slice_mesh_plane(mesh_z, plane_normal=-nz, plane_origin=list(nz*(z+1)*scale.z+oz))
525-
526-
if len(mesh_z.vertices) == 0:
527-
continue
528-
529-
# test for totally degenerate meshes by checking if
530-
# all of two axes match, meaning the mesh must be a
531-
# point or a line.
532-
if np.sum([ np.all(mesh_z.vertices[:,i] == mesh_z.vertices[0,i]) for i in range(3) ]) >= 2:
533-
continue
534-
535-
yield mesh_z, (x, y, z)
536-
537514
def retriangulate_mesh(
538-
mesh: trimesh.Trimesh,
515+
mesh: Mesh,
539516
offset: np.ndarray,
540-
grid_size: Vec,
541517
scale: Vec,
542-
) -> trimesh.Trimesh:
518+
) -> Mesh:
543519
"""
544520
Retriangulate the input mesh to avoid any cases where the boundaries of a triangle are split across the boundaries of the submeshes
545521
"""
546-
new_mesh = trimesh.Trimesh()
547-
548-
for submesh, _ in generate_gridded_submeshes(
549-
mesh, offset, grid_size, scale
550-
):
551-
new_mesh = trimesh.util.concatenate(new_mesh, submesh)
552-
553-
return new_mesh
522+
chunks = zmesh.chunk_mesh(mesh, scale, offset)
523+
new_mesh = zmesh.Mesh.concatenate(*chunks.values(), id=mesh.segid)
524+
return new_mesh.merge_close_vertices(radius=1e-5)
554525

555526
def create_octree_level_from_mesh(mesh, chunk_shape, lod, num_lods, offset, grid_length):
556527
"""
@@ -559,7 +530,6 @@ def create_octree_level_from_mesh(mesh, chunk_shape, lod, num_lods, offset, grid
559530
560531
This creates (2^lod)^3 submeshes.
561532
"""
562-
mesh = trimesh.Trimesh(vertices=mesh.vertices, faces=mesh.faces)
563533
scale = Vec(*(np.array(chunk_shape) * (2**lod)))
564534
grid_size = Vec(*(np.ceil(grid_length / scale)), dtype=int)
565535

@@ -568,26 +538,21 @@ def create_octree_level_from_mesh(mesh, chunk_shape, lod, num_lods, offset, grid
568538
# at the higher level of the octree
569539
if lod > 0:
570540
upper_grid_scale = Vec(*(np.array(chunk_shape) * (2 ** (lod - 1))))
571-
upper_grid_shape = Vec(*np.ceil(grid_length / upper_grid_scale), dtype=int)
572-
mesh = retriangulate_mesh(mesh, offset, upper_grid_shape, upper_grid_scale)
541+
mesh = retriangulate_mesh(mesh, offset, upper_grid_scale)
573542

574543
if lod == num_lods - 1:
575544
return ([Mesh(mesh.vertices, mesh.faces)], ((0, 0, 0),))
576545

577-
submeshes = []
578-
nodes = []
579-
for submesh, node in generate_gridded_submeshes(
580-
mesh, offset, grid_size, scale
581-
):
582-
submeshes.append(submesh)
583-
nodes.append(node)
546+
grid = zmesh.chunk_mesh(mesh, scale, offset)
584547

585548
# Sort in Z-curve order
586-
submeshes, nodes = zip(
587-
*sorted(zip(submeshes, nodes),
588-
key=functools.cmp_to_key(lambda x, y: cmp_zorder(x[1], y[1])))
549+
nodes, submeshes = zip(
550+
*sorted(
551+
grid.items(),
552+
key=functools.cmp_to_key(lambda x, y: cmp_zorder(x[0], y[0]))
553+
)
589554
)
590-
# convert back from trimesh to CV Mesh class
591-
submeshes = [ Mesh(m.vertices, m.faces) for m in submeshes ]
555+
# convert back from zmesh.Mesh to CV Mesh class
556+
submeshes = [ Mesh(m.vertices, m.faces, segid=mesh.segid) for m in submeshes ]
592557

593558
return (submeshes, nodes)

requirements.txt

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -25,4 +25,4 @@ task-queue>=2.4.0
2525
tqdm
2626
trimesh[easy]
2727
xs3d>=1.11.0
28-
zmesh>=1.4,<2.0
28+
zmesh>=1.9.0,<2.0

0 commit comments

Comments
 (0)