GPU-Accelerated Mesh Processing for Physical Simulation and Scientific Visualization in Any Dimension
"It's not just a bag of triangles -- it's a fast bag of triangles!"
The word "mesh" means different things to different communities:
- CFD/FEM engineers think "volume mesh" (3D tetrahedra filling a 3D domain)
- Graphics programmers think "surface mesh" (2D triangles in 3D space)
- Computer vision researchers think "point cloud" (0D vertices in 3D space)
- Robotics engineers think "curves" (1D edges in 2D or 3D space)
TorchMesh handles all of these in a unified, dimensionally-generic framework. More precisely, TorchMesh operates on arbitrary-dimensional pure simplicial complexes embedded in arbitrary-dimensional Euclidean spaces.
This means you can work with:
- 2D triangles in 2D space (planar meshes for 2D simulations)
- 2D triangles in 3D space (surface meshes for graphics/CFD)
- 3D tetrahedra in 3D space (volume meshes for FEM/CFD)
- 1D edges in 3D space (curve meshes for path planning)
- Any other n-dimensional manifold in m-dimensional space (where n ≤ m)
all with the same API. TorchMesh's API design takes heavy inspiration from PyVista, but it is designed to be a) end-to-end GPU-accelerated, b) dimensionally generic, c) autograd-differentiable where possible, d) allow arbitrary-rank tensor fields as data, and e) support nested field data. TorchMesh is extremely fast and lightweight, depending only on PyTorch and TensorDict (an official PyTorch data structure).
The only restriction: meshes must be simplicial (composed of points, line segments, triangles, tetrahedra, and higher-dimensional n-simplices). This enables a plethora of rigorous differential geometry + discrete calculus computations, as well as significant performance benefits.
Core Capabilities:
- GPU-Accelerated: All operations vectorized with PyTorch, run natively on CUDA
- Dimensionally Generic: Works with n-D manifolds embedded in m-D spaces
- TensorDict Integration: Structured data management with TensorDict and automatic device handling
- Differentiable: Most features offer seamless integration with PyTorch autograd
Mathematical Operations:
- Discrete Calculus: Gradient, divergence, curl, Laplace-Beltrami operator (Note: these are all the core ingredients required for a high-performance manifold PDE solver for many PDEs of industrial interest.)
- Both DEC (Discrete Exterior Calculus) and LSQ (Least-Squares) methods
- Intrinsic (tangent space) and extrinsic (ambient space) derivatives
- Differential Geometry: Gaussian curvature, mean curvature, normals, tangent spaces
- Curvature Analysis: Angle defect (intrinsic) and cotangent Laplacian (extrinsic) methods
Mesh Operations:
- Subdivision: Linear, Loop (C²), and Butterfly (interpolating) schemes
- Smoothing: Laplacian smoothing with feature preservation
- Remeshing: Uniform remeshing via clustering (dimension-agnostic)
- Repair: Remove duplicates, fix orientation, fill holes, clean topology
Analysis Tools:
- Topology: Boundary detection, watertight/manifold checking
- Neighbors: Point-to-point, point-to-cell, cell-to-cell, cell-to-point adjacency, computed and stored efficiently
- Quality Metrics: Aspect ratio, edge lengths, angles, quality scores
- Spatial Queries: BVH-accelerated point containment and nearest-cell search
pip install torchmeshThis installs TorchMesh with its core dependencies: PyTorch and TensorDict.
pip install torchmesh[io]This adds PyVista, matplotlib, and numpy, which are soft dependencies used only for:
- Loading/saving mesh files (STL, OBJ, VTK/VTU/VTP/VTM, PLY, etc.)
- Interactive 3D visualization
- Uniform remeshing (via pyacvd)
All other TorchMesh operations (calculus, curvature, subdivision, etc.) work without these dependencies and can stay entirely on GPU.
git clone https://github.com/peterdsharpe/torchmesh.git
cd torchmesh
pip install -e ".[dev]"import torch
from torchmesh import Mesh
# Create a triangle mesh in 2D
points = torch.tensor([
[0.0, 0.0],
[1.0, 0.0],
[0.5, 1.0]
])
cells = torch.tensor([[0, 1, 2]]) # Indicates that points 0, 1, and 2 form a cell (in 2D, a triangle)
mesh = Mesh(points=points, cells=cells)
print(mesh)Output: (dimensionality is inferred from points and cells shapes)
Mesh(manifold_dim=2, spatial_dim=2, n_points=3, n_cells=1)
point_data : {}
cell_data : {}
global_data: {}
# Scalar data (shape: n_points or n_cells)
mesh.point_data["temperature"] = torch.tensor([300.0, 350.0, 325.0])
mesh.cell_data["pressure"] = torch.tensor([101.3])
# Vector data (shape: n_points × n_spatial_dims or n_cells × n_spatial_dims)
mesh.point_data["velocity"] = torch.tensor([[1.0, 0.5], [0.8, 1.2], [0.0, 0.9]])
# Tensor data (shape: n_cells × n_spatial_dims × n_spatial_dims)
# Reynolds stress tensor: symmetric 2×2 tensor for each cell
mesh.cell_data["reynolds_stress"] = torch.tensor([[[2.1, 0.3], [0.3, 1.8]]])
print(mesh)Output: (data fields show the trailing dimensions)
Mesh(manifold_dim=2, spatial_dim=2, n_points=3, n_cells=1)
point_data : {temperature: (), velocity: (2,)}
cell_data : {pressure: (), reynolds_stress: (2, 2)}
global_data: {}
from torchmesh.io import from_pyvista
import pyvista as pv
# Load any mesh format PyVista supports
pv_mesh = pv.examples.load_airplane() # See pyvista.org for more datasets
mesh = from_pyvista(pv_mesh)
# Or, equivalently:
from torchmesh.examples.pyvista_datasets.airplane import load
mesh = load()
print(mesh)Output:
Mesh(manifold_dim=2, spatial_dim=3, n_points=1335, n_cells=2452)
point_data : {}
cell_data : {}
global_data: {}
This is a 2D surface mesh (triangles) embedded in 3D space - a typical graphics/CAD mesh.
Then, with mesh.draw(), you can visualize the mesh:
Starting with the airplane mesh, we can compute its surface Gaussian curvature:
mesh = mesh.subdivide(levels=2, filter="loop")
mesh.point_data["gaussian_curvature"] = mesh.gaussian_curvature_vertices
mesh.draw(
point_scalars="gaussian_curvature",
show_edges=False,
...
)Warmer colors indicate positive Gaussian curvature (convex regions), cooler colors indicate negative Gaussian curvature (concave regions).
Or, compute the mean curvature:
mesh.point_data["mean_curvature"] = mesh.mean_curvature_vertices
mesh.draw(
point_scalars="mean_curvature",
show_edges=False,
...
)Warmer colors indicate positive mean curvature (convex regions), cooler colors indicate negative mean curvature (concave regions).
# Create scalar field: T = x + 2y
mesh.point_data["temperature"] = mesh.points[:, 0] + 2 * mesh.points[:, 1]
# Compute gradient using least-squares reconstruction
mesh_with_grad = mesh.compute_point_derivatives(keys="temperature", method="lsq")
grad_T = mesh_with_grad.point_data["temperature_gradient"]
print(f"Gradient shape: {grad_T.shape}") # (n_points, n_spatial_dims)
print(f"∇T = {grad_T[0]}") # tensor([1.0000, 2.0000])# Move entire mesh and all data to GPU
mesh_gpu = mesh.to("cuda")
# Compute on GPU
K_gpu = mesh_gpu.gaussian_curvature_vertices
# Move back to CPU
mesh_cpu = mesh_gpu.to("cpu")Comprehensive overview of TorchMesh capabilities:
| Feature | Status | Notes |
|---|---|---|
| Core Operations | ||
| Mesh creation & manipulation | ✅ | n-dimensional simplicial meshes |
| Point/cell/global data | ✅ | TensorDict-based (including nested data) |
| GPU acceleration | ✅ | Full CUDA support |
| Merge multiple meshes | ✅ | |
| Device management (CPU/GPU) | ✅ | |
| Calculus | ||
| Gradient/Jacobian (LSQ) | ✅ | Weighted least-squares reconstruction |
| Gradient/Jacobian (DEC) | ✅ | Via sharp operator |
| Divergence (LSQ) | ✅ | Component-wise gradients |
| Divergence (DEC) | ✅ | Explicit dual volume formula |
| Curl (LSQ, 3D only) | ✅ | Antisymmetric Jacobian |
| Laplace-Beltrami (DEC) | ❌ | Work in progress |
| Intrinsic derivatives | ✅ | Tangent space projection |
| Extrinsic derivatives | ✅ | Ambient space |
| Geometry | ||
| Cell centroids | ✅ | Arithmetic mean of vertices |
| Cell areas/volumes | ✅ | Gram determinant method |
| Cell normals | ✅ | Generalized cross product |
| Point normals | ✅ | Area-weighted from adjacent cells |
| Facet extraction | ✅ | Extract all (n-1)-dimensional simplices |
| Boundary detection and extraction | ✅ | Extract only the boundary (n-1)-dimensional simplices |
| Curvature | ||
| Gaussian curvature (vertices) | ✅ | Angle defect method |
| Gaussian curvature (cells) | ✅ | |
| Mean curvature | ✅ | Cotangent Laplacian |
| Subdivision | ||
| Linear | ✅ | Midpoint subdivision |
| Loop | ✅ | C² smooth, approximating |
| Butterfly | ✅ | Interpolating |
| Smoothing | ||
| Laplacian smoothing | ✅ | |
| Remeshing | ||
| Uniform remeshing | ✅ | Clustering-based |
| Spatial Queries | ||
| BVH construction | ✅ | |
| Point containment | ✅ | |
| Nearest cell search | ✅ | |
| Data interpolation | ✅ | Barycentric coordinates |
| Sampling | ||
| Random points on cells | ✅ | Dirichlet distribution |
| Data sampling at points | ✅ | |
| Transformations | ||
| Translation | ✅ | |
| Rotation | ✅ | In 2D or 3D (angle-axis); for higher dimensions rotation is ill-defined, use transform() instead |
| Scaling | ✅ | Uniform or anisotropic |
| Arbitrary matrix transform | ✅ | |
| Extrusion | ✅ | Manifold → higher dimension |
| Projection / Intersection | ❌ | Manifold → lower dimension; work in progress |
| Neighbors & Adjacency | ||
| Point-to-points | ✅ | Graph edges |
| Point-to-cells | ✅ | Vertex star |
| Cell-to-cells | ✅ | Shared facets |
| Cells-to-points | ✅ | Cell vertices |
| Ragged array format | ✅ | Efficient sparse encoding via offsets + indices |
| Geometry | ||
| Delaunay triangulation from points | ❌ | Work in progress |
| Voronoi areas | ✅ | |
| Convex hulls | ❌ | Work in progress |
| Topology & Repair | ||
| Watertight detection | ✅ | |
| Manifold detection | ✅ | |
| Remove duplicate vertices | ✅ | |
| Remove duplicate cells | ✅ | |
| Remove degenerate cells | ✅ | |
| Remove isolated vertices | ✅ | |
| Fix orientation | ✅ | |
| Fill holes | ✅ | |
| Clean mesh (all-in-one) | ✅ | |
| Validation & Analysis | ||
| Quality metrics | ✅ | Aspect ratio, angles, edge ratios |
| Mesh statistics | ✅ | |
| I/O | ||
| PyVista integration | ✅ | All PyVista-supported formats |
| Visualization | ||
| Matplotlib backend | ✅ | 2D/3D plotting |
| PyVista backend | ✅ | Interactive 3D |
| Scalar colormapping on points or cells | ✅ | Auto L2-norm for vectors |
Legend: ✅ Complete | ❌ Not Implemented, but planned
import numpy as np
mesh_translated = mesh.translate([1.0, 0.0, 0.0])
mesh_rotated = mesh.rotate(axis=[0, 0, 1], angle=np.pi/4)
mesh_scaled = mesh.scale(2.0) # Or [2.0, 1.0, 0.5] for anisotropicrefined = mesh.subdivide(levels=2, filter="linear") # Topology only
smooth = mesh.subdivide(levels=2, filter="loop") # C² continuous
interp = mesh.subdivide(levels=2, filter="butterfly") # Interpolatingfrom torchmesh.calculus import compute_divergence_points_lsq, compute_curl_points_lsq
# Gradient
mesh.point_data["pressure"] = (mesh.points ** 2).sum(dim=-1)
mesh_grad = mesh.compute_point_derivatives(keys="pressure", method="lsq")
grad_p = mesh_grad.point_data["pressure_gradient"] # (n_points, n_spatial_dims)
# Divergence and curl
mesh.point_data["velocity"] = mesh.points.clone()
div_v = compute_divergence_points_lsq(mesh, mesh.point_data["velocity"])
curl_v = compute_curl_points_lsq(mesh, mesh.point_data["velocity"]) # 3D only
# For surfaces: intrinsic (tangent space) vs extrinsic (ambient space)
grad_intrinsic = mesh.compute_point_derivatives(keys="T", gradient_type="intrinsic")
grad_extrinsic = mesh.compute_point_derivatives(keys="T", gradient_type="extrinsic")# Boundary detection
boundary = mesh.get_boundary_mesh()
facets = mesh.get_facet_mesh()
is_watertight = mesh.is_watertight()
# Repair
clean = mesh.clean() # All-in-onefrom torchmesh.spatial import BVH
bvh = BVH.from_mesh(mesh)
cell_candidates = bvh.find_candidate_cells(torch.rand(1000, 3))
sampled = mesh.sample_data_at_points(query_points, data_source="points")Note that these use an efficient sparse (indices, offsets) encoding of the adjacency relationships, which is used internally for all computations. (See the dedicated torchmesh.neighbors._adjacency.py module.) You can convert these to a typical ragged list-of-lists representation with .to_list(), which is useful for debugging or interoperability, at the cost of performance:
point_neighbors = mesh.get_point_to_points_adjacency().to_list()
cell_neighbors = mesh.get_cell_to_cells_adjacency().to_list()The Mesh class is a tensorclass with five components:
Mesh(
points: torch.Tensor, # (n_points, n_spatial_dims)
cells: torch.Tensor, # (n_cells, n_manifold_dims + 1), integer dtype
point_data: TensorDict, # Per-vertex data
cell_data: TensorDict, # Per-cell data
global_data: TensorDict, # Mesh-level data
)All data moves together with .to("cuda") or .to("cpu"). Expensive computations (centroids, normals, curvature) are automatically cached in the data dictionaries with keys starting with _.
TorchMesh is built on three principles:
- Correctness First: Rigorous mathematical foundations, extensive validation
- Performance Second: Fully vectorized GPU operations, no Python loops over mesh elements
- Usability Third: Clean APIs that don't sacrifice power for simplicity
Key design decisions enable these principles:
- Simplicial meshes only (enables rigorous discrete exterior calculus)
- Explicit dimensionality (
n_spatial_dims,n_manifold_dimsas first-class concepts) - Fail loudly with helpful error messages (no silent failures)
- Examples: See
examples/directory for runnable demonstrations - Tests: See
test/directory for comprehensive test suite showing usage patterns - Source: Explore
src/torchmesh/for implementation details
Module Organization:
torchmesh.calculus- Discrete differential operatorstorchmesh.curvature- Gaussian and mean curvaturetorchmesh.subdivision- Mesh refinement schemestorchmesh.boundaries- Boundary detection and facet extractiontorchmesh.neighbors- Adjacency computationstorchmesh.spatial- BVH and spatial queriestorchmesh.sampling- Point sampling and interpolationtorchmesh.transformations- Geometric operationstorchmesh.repair- Mesh cleaning and topology repairtorchmesh.validation- Quality metrics and statisticstorchmesh.visualization- Matplotlib and PyVista backendstorchmesh.io- PyVista import/exporttorchmesh.examples- Example mesh generators
If you use TorchMesh in your research, please cite:
@software{torchmesh2025,
author = {Sharpe, Peter},
title = {TorchMesh: GPU-Accelerated Mesh Processing for Scientific Computing},
year = {2025},
url = {https://github.com/peterdsharpe/torchmesh},
version = {0.1.0}
}TorchMesh is licensed under the Apache License 2.0. See LICENSE for details.
TorchMesh draws inspiration for its API design and mathematical foundation from:
- PyTorch team for the foundational deep learning framework
- PyVista team for the excellent 3D visualization and I/O library
- Discrete Exterior Calculus: Desbrun, Hirani, Leok, Marsden (2005) - arXiv:math/0508341, and the eponymous dissertation on Discrete Exterior Calculus by Hirani (2003)
- Discrete Differential Operators: Meyer, Desbrun, Schröder, Barr (2003) - Discrete Differential-Geometry Operators for Triangulated 2-Manifolds
- Loop Subdivision: Loop (1987) - Smooth Subdivision Surfaces Based on Triangles
- Butterfly Subdivision: Zorin, Schröder, Sweldens (1996) - Interpolating Subdivision for Meshes with Arbitrary Topology
Questions? Issues? Feature requests?
Open an issue on GitHub or start a discussion!


