Expand description
§delaunay
This is a library for computing the Delaunay triangulation of a set of n-dimensional points in a simplicial complex inspired by CGAL.
§Features
- d-dimensional Delaunay triangulations
- d-dimensional convex hulls
- Generic floating-point coordinate types (supports
f32,f64, and other types implementingCoordinateScalar) - Copy-able data types associated with vertices and cells (see
DataTypefor constraints) - Serialization/Deserialization with serde
§Basic Usage
This library handles arbitrary dimensions (subject to numerical issues). Here’s a 4D triangulation example:
use delaunay::core::triangulation_data_structure::Tds;
use delaunay::vertex;
// Create a 4D triangulation (4-dimensional space!)
let mut tds: Tds<f64, Option<()>, Option<()>, 4> = Tds::empty();
// Add vertices incrementally - triangulation evolves automatically
tds.add(vertex!([0.0, 0.0, 0.0, 0.0])).unwrap(); // 1 vertex, 0 cells
tds.add(vertex!([1.0, 0.0, 0.0, 0.0])).unwrap(); // 2 vertices, 0 cells
tds.add(vertex!([0.0, 1.0, 0.0, 0.0])).unwrap(); // 3 vertices, 0 cells
tds.add(vertex!([0.0, 0.0, 1.0, 0.0])).unwrap(); // 4 vertices, 0 cells
assert_eq!(tds.number_of_cells(), 0);
tds.add(vertex!([0.0, 0.0, 0.0, 1.0])).unwrap(); // 5 vertices, 1 cell (first 4-simplex!)
tds.add(vertex!([0.2, 0.2, 0.2, 0.2])).unwrap(); // 6 vertices, multiple cells
assert_eq!(tds.number_of_vertices(), 6);
assert_eq!(tds.dim(), 4); // Full 4D triangulation
assert!(tds.number_of_cells() > 1); // Bowyer-Watson creates additional 4-simplices
assert!(tds.is_valid().is_ok()); // Maintains Delaunay property in 4DKey insight: The transition happens at D+1 vertices (5 vertices for 4D), where the first 4-simplex (5-vertex cell) is created. Additional vertices trigger the Bowyer-Watson algorithm to maintain the 4D Delaunay triangulation.
§Convex Hull Extraction
Extract d-dimensional convex hulls from Delaunay triangulations:
use delaunay::core::triangulation_data_structure::Tds;
use delaunay::geometry::algorithms::convex_hull::ConvexHull;
use delaunay::geometry::point::Point;
use delaunay::geometry::traits::coordinate::Coordinate;
use delaunay::vertex;
// Create two tetrahedrons sharing a triangular facet (double tetrahedron)
let vertices = vec![
// Shared triangular facet vertices (forms base of both tetrahedrons)
vertex!([0.0, 0.0, 0.0]), // Shared vertex A
vertex!([2.0, 0.0, 0.0]), // Shared vertex B
vertex!([1.0, 2.0, 0.0]), // Shared vertex C
// Apex of first tetrahedron (above the shared facet)
vertex!([1.0, 0.7, 1.5]), // First tet apex
// Apex of second tetrahedron (below the shared facet)
vertex!([1.0, 0.7, -1.5]), // Second tet apex
];
let tds: Tds<f64, Option<()>, Option<()>, 3> = Tds::new(&vertices).unwrap();
// Extract the convex hull (boundary facets of the triangulation)
let hull: ConvexHull<f64, Option<()>, Option<()>, 3> =
ConvexHull::from_triangulation(&tds).unwrap();
println!("Convex hull has {} facets in {}D", hull.facet_count(), hull.dimension());
// Test point containment
let inside_point = Point::new([1.0, 0.5, 0.5]);
let outside_point = Point::new([3.0, 3.0, 3.0]);
assert!(!hull.is_point_outside(&inside_point, &tds).unwrap()); // Inside the hull
assert!(hull.is_point_outside(&outside_point, &tds).unwrap()); // Outside the hull
// Find visible facets from an external point (useful for incremental construction)
let visible_facets = hull.find_visible_facets(&outside_point, &tds).unwrap();
println!("Point sees {} out of {} facets", visible_facets.len(), hull.facet_count());
// Works in any dimension!
let vertices_4d = vec![
vertex!([0.0, 0.0, 0.0, 0.0]),
vertex!([1.0, 0.0, 0.0, 0.0]),
vertex!([0.0, 1.0, 0.0, 0.0]),
vertex!([0.0, 0.0, 1.0, 0.0]),
vertex!([0.0, 0.0, 0.0, 1.0]),
];
let tds_4d: Tds<f64, Option<()>, Option<()>, 4> = Tds::new(&vertices_4d).unwrap();
let hull_4d: ConvexHull<f64, Option<()>, Option<()>, 4> =
ConvexHull::from_triangulation(&tds_4d).unwrap();
assert_eq!(hull_4d.facet_count(), 5); // 4-simplex has 5 boundary facets
assert_eq!(hull_4d.dimension(), 4); // 4D convex hull§Triangulation Invariants
The library maintains several critical triangulation invariants that ensure geometric and topological correctness:
§Invariant Enforcement
| Invariant Type | Enforcement Location | Method | Timing |
|---|---|---|---|
| Delaunay Property | bowyer_watson::find_bad_cells() | Empty circumsphere test using insphere() | Proactive (during construction) |
| Facet Sharing | validate_facet_sharing() | Each facet shared by ≤ 2 cells | Reactive (via validation) |
| No Duplicate Cells | validate_no_duplicate_cells() | No cells with identical vertex sets | Reactive (via validation) |
| Neighbor Consistency | validate_neighbors_internal() | Mutual neighbor relationships | Reactive (via validation) |
| Cell Validity | CellBuilder::validate() (vertex count) + cell.is_valid() (comprehensive) | Construction + runtime validation | Both (construction + validation) |
| Vertex Validity | Point::from() (coordinates) + UUID auto-gen + vertex.is_valid() | Construction + runtime validation | Both (construction + validation) |
The Delaunay property (empty circumsphere) is enforced proactively during construction by removing violating cells, while structural invariants are enforced reactively through validation methods.
For detailed information, see:
core::algorithms::bowyer_watson- Primary invariant enforcement during triangulation constructioncore::triangulation_data_structure::Tds::is_valid- Comprehensive validation of all invariants
§References
The algorithms and geometric predicates in this library are based on established computational geometry literature. For a comprehensive list of academic references and bibliographic citations, see the REFERENCES.md file in the repository.
§Project History
Versions ≤ 0.1.0 were maintained at old repo. Versions ≥ 0.3.4 are maintained here.
See https://docs.rs/delaunay/0.1.0 for historical documentation. See https://docs.rs/delaunay for the latest documentation.
Modules§
- core
- The
coremodule contains the primary data structures and algorithms for building and manipulating Delaunay triangulations. - geometry
- Contains geometric types including the
Pointstruct and geometry predicates. - prelude
- A prelude module that re-exports commonly used types and macros. This makes it easier to import the most commonly used items from the crate.
Macros§
- cell
Deprecated - Convenience macro for creating cells with less boilerplate.
- vertex
- Convenience macro for creating vertices with less boilerplate.
Functions§
- is_
normal - The function
is_normalchecks that structs implementautotraits. Traits are checked at compile time, so this function is only used for testing.