Quick Start Guide
This guide will walk you through the essential features of MeshMash with real, working examples.
Creating a Mesh
MeshMash works with meshes represented as tuples of vertices and faces:
import numpy as np
# Define vertices (N x 3 array)
vertices = np.array([
[ 0.0 , 0.0 , 0.0 ],
[ 1.0 , 0.0 , 0.0 ],
[ 0.0 , 1.0 , 0.0 ],
[ 0.0 , 0.0 , 1.0 ]
], dtype = np.float64)
# Define faces (M x 3 array of vertex indices)
faces = np.array([
[ 0 , 1 , 2 ],
[ 0 , 1 , 3 ],
[ 0 , 2 , 3 ],
[ 1 , 2 , 3 ]
], dtype = np.int32)
# Create mesh tuple
mesh = (vertices, faces)
MeshMash also accepts objects with vertices and faces attributes, such as PyVista PolyData.
Basic Operations
Computing the Laplacian
The cotangent Laplacian is fundamental to spectral mesh processing:
from meshmash import cotangent_laplacian
# Compute cotangent Laplacian and area matrix
L, M = cotangent_laplacian(mesh)
print ( f "Laplacian shape: { L.shape } " ) # (n_vertices, n_vertices)
print ( f "Area matrix shape: { M.shape } " ) # (n_vertices, n_vertices)
For non-manifold meshes, use the robust formulation: L, M = cotangent_laplacian(mesh, robust = True , mollify_factor = 1e-5 )
Spectral Decomposition
Decompose the mesh into its spectral components:
from meshmash import decompose_mesh
# Compute first 50 eigenpairs
eigenvalues, eigenvectors = decompose_mesh(
mesh,
n_components = 50 ,
robust = True
)
print ( f "Eigenvalues shape: { eigenvalues.shape } " ) # (50,)
print ( f "Eigenvectors shape: { eigenvectors.shape } " ) # (n_vertices, 50)
Advanced Features
Heat Kernel Signatures (HKS)
Compute multi-scale geometric features using HKS:
Import the function
from meshmash import compute_hks
Set parameters
max_eigenvalue = 1e-8 # Frequency cutoff
n_components = 32 # Number of scales
Compute HKS features
hks_features = compute_hks(
mesh,
max_eigenvalue = max_eigenvalue,
n_components = n_components,
robust = True ,
verbose = True
)
print ( f "HKS features shape: { hks_features.shape } " ) # (n_vertices, 32)
HKS features are invariant to isometric deformations, making them ideal for shape matching and analysis.
Mesh Utilities
Extract Largest Component
Subset by Indices
Convert to Graph
Combine Meshes
from meshmash import largest_mesh_component
# Get only the largest connected component
largest = largest_mesh_component(mesh)
vertices, faces = largest
Spectral Filtering
Custom Spectral Filters
Apply custom filters in the spectral domain:
from meshmash import spectral_geometry_filter, construct_bspline_filter
# Create a B-spline filter
filter_func = construct_bspline_filter(
e_min = 0.0 ,
e_max = 1e-8 ,
n_components = 32
)
# Apply filter to mesh geometry
filtered_features = spectral_geometry_filter(
mesh,
filter = filter_func,
max_eigenvalue = 1e-8 ,
band_size = 50 ,
robust = True ,
verbose = True
)
Geometry Vectors
Compute geometric descriptors using spectral methods:
from meshmash import compute_geometry_vectors
geometry_features = compute_geometry_vectors(
mesh,
max_eigenvalue = 1e-8 ,
n_components = 32 ,
robust = True
)
print ( f "Shape: { geometry_features.shape } " ) # (n_vertices, 32)
Working with PyVista
MeshMash integrates seamlessly with PyVista:
import pyvista as pv
from meshmash import mesh_to_poly, poly_to_mesh, compute_hks
# Load mesh from file
poly = pv.read( 'model.ply' )
# Convert to MeshMash format
mesh = poly_to_mesh(poly)
# Process with MeshMash
hks = compute_hks(mesh, max_eigenvalue = 1e-8 , n_components = 16 )
# Add features to PyVista mesh for visualization
poly[ 'hks_0' ] = hks[:, 0 ]
poly.plot( scalars = 'hks_0' , cmap = 'viridis' )
Hierarchical Mesh Splitting
For large meshes, use spectral splitting:
from meshmash import fit_mesh_split, apply_mesh_split
# Split mesh hierarchically
split_info = fit_mesh_split(
mesh,
max_vertex_threshold = 20000 , # Maximum vertices per submesh
min_vertex_threshold = 100 , # Minimum vertices per submesh
verbose = True
)
# Access submeshes
for submesh_idx, submesh in enumerate (split_info[ 'submeshes' ]):
vertices, faces = submesh
print ( f "Submesh { submesh_idx } : { len (vertices) } vertices" )
Mesh splitting is stochastic. Results may vary slightly between runs due to the spectral decomposition.
Connected Components
Analyze mesh topology:
from meshmash import mesh_connected_components, threshold_mesh_by_component_size
# Iterate over components larger than threshold
for component in mesh_connected_components(mesh, size_threshold = 100 ):
vertices, faces = component
print ( f "Component with { len (vertices) } vertices" )
# Filter mesh by component size
filtered_mesh, kept_indices = threshold_mesh_by_component_size(
mesh,
size_threshold = 1000
)
Label Propagation
Propagate labels across mesh surfaces:
from meshmash import label_propagation
import numpy as np
# Create initial labels (most unlabeled = -1)
labels = np.full( len (vertices), - 1 , dtype = np.int32)
labels[ 0 ] = 0 # Label first vertex as class 0
labels[ 10 ] = 1 # Label 11th vertex as class 1
# Propagate labels
propagated_labels = label_propagation(mesh, labels)
Next Steps
Mesh Splitting Learn how to split large meshes for parallel processing
Spectral Filtering Apply custom spectral filters to meshes
Heat Kernel Signature Compute HKS features for shape analysis
GitHub View source code and contribute