Crate delaunay

Crate delaunay 

Source
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 implementing CoordinateScalar)
  • Copy-able data types associated with vertices and cells (see DataType for constraints)
  • Serialization/Deserialization with serde

§Basic Usage

This library handles arbitrary dimensions (subject to numerical issues). Here’s a 4D triangulation example:

use delaunay::prelude::*;

// Create a 4D Delaunay triangulation (4-dimensional space!)
let vertices = 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]),  // 5 vertices (D+1) creates first 4-simplex
    vertex!([0.2, 0.2, 0.2, 0.2]),  // Additional vertex uses incremental insertion
];

let dt: DelaunayTriangulation<_, (), (), 4> =
    DelaunayTriangulation::new(&vertices).unwrap();

assert_eq!(dt.number_of_vertices(), 6);
assert_eq!(dt.dim(), 4);                    // Full 4D triangulation
assert!(dt.number_of_cells() > 1);          // Cavity-based insertion creates additional 4-simplices

Key insight: The triangulation uses efficient incremental cavity-based insertion after building an initial simplex from the first D+1 vertices (5 vertices for 4D).

§Convex Hull Extraction

Extract d-dimensional convex hulls from Delaunay triangulations:

use delaunay::prelude::*;

// Create two tetrahedrons sharing a triangular facet (double tetrahedron)
let vertices: Vec<_> = 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 dt: DelaunayTriangulation<_, (), (), 3> =
    DelaunayTriangulation::new(&vertices).unwrap();

// Extract the convex hull (boundary facets of the triangulation)
let hull = ConvexHull::from_triangulation(dt.triangulation()).unwrap();

println!("Convex hull has {} facets in {}D", hull.number_of_facets(), 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, dt.triangulation()).unwrap());  // Inside the hull
assert!(hull.is_point_outside(&outside_point, dt.triangulation()).unwrap());   // Outside the hull

// Find visible facets from an external point (useful for incremental construction)
let visible_facets = hull.find_visible_facets(&outside_point, dt.triangulation()).unwrap();
println!("Point sees {} out of {} facets", visible_facets.len(), hull.number_of_facets());

// Works in any dimension!
let vertices_4d: Vec<_> = 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 dt_4d: DelaunayTriangulation<_, (), (), 4> =
    DelaunayTriangulation::new(&vertices_4d).unwrap();
let hull_4d = ConvexHull::from_triangulation(dt_4d.triangulation()).unwrap();

assert_eq!(hull_4d.number_of_facets(), 5);  // 4-simplex has 5 boundary facets
assert_eq!(hull_4d.dimension(), 4);     // 4D convex hull

§Triangulation invariants and validation hierarchy

The crate is organized as a small validation stack, where each layer adds additional invariants on top of the preceding one:

  • Tds (Triangulation Data Structure) stores the combinatorial / structural representation. Level 2 (structural) validation checks invariants such as:

    • Vertex mappings – every vertex UUID has a corresponding key and vice versa.
    • Cell mappings – every cell UUID has a corresponding key and vice versa.
    • No duplicate cells – no two maximal cells share the same vertex set.
    • Facet sharing – each facet is shared by at most 2 cells (1 on the boundary, 2 in the interior).
    • Neighbor consistency – neighbor relationships are mutual and reference a shared facet.

    These checks are surfaced via Tds::is_valid (structural only) and Tds::validate (Levels 1–2, elements + structural). For cumulative diagnostics across the full stack, use DelaunayTriangulation::validation_report.

  • Triangulation builds on the TDS and validates manifold topology. Level 3 (topology) validation is performed by Triangulation::is_valid (Level 3 only) and Triangulation::validate (Levels 1–3), which:

    • Strengthens facet sharing to the manifold facet property: each facet belongs to exactly 1 cell (boundary) or exactly 2 cells (interior).
    • Checks the Euler characteristic of the triangulation (using the topology module).
  • DelaunayTriangulation builds on Triangulation and validates the geometric Delaunay condition. Level 4 (Delaunay property) validation is performed by DelaunayTriangulation::is_valid (Level 4 only) and DelaunayTriangulation::validate (Levels 1–4). Construction is designed to satisfy the Delaunay property, but in rare cases it may be violated for near-degenerate inputs (see Issue #120).

§Validation

The crate exposes four validation levels (element → structural → manifold → Delaunay). The canonical guide (when to use each level, complexity, examples, troubleshooting) lives in docs/validation.md: https://github.com/acgetchell/delaunay/blob/main/docs/validation.md

In brief:

  • Level 2 (structural / Tds): dt.tds().is_valid() for a quick check, or dt.tds().validate() for Levels 1–2.
  • Level 3 (topology / Triangulation): dt.triangulation().is_valid() for topology-only checks, or dt.triangulation().validate() for Levels 1–3.
  • Level 4 (Delaunay / DelaunayTriangulation): dt.is_valid() for the empty-circumsphere property, or dt.validate() for Levels 1–4.
  • Full diagnostics: dt.validation_report() returns all violated invariants across Levels 1–4.

§Automatic topology validation during insertion (ValidationPolicy)

In addition to explicit validation calls, incremental construction (new() / insert*()) can run an automatic Level 3 topology validation pass after insertion, controlled by ValidationPolicy.

The default is ValidationPolicy::OnSuspicion: Level 3 validation runs only when insertion takes a suspicious path (e.g. perturbation retries, repair loops, or neighbor-pointer repairs that actually changed pointers).

This automatic pass only runs Level 3 (Triangulation::is_valid()). It does not run Level 4.

use delaunay::prelude::*;
let mut dt: DelaunayTriangulation<_, (), (), 3> = DelaunayTriangulation::new(&vertices).unwrap();

// Default:
assert_eq!(dt.validation_policy(), ValidationPolicy::OnSuspicion);

// Tests/debugging:
dt.set_validation_policy(ValidationPolicy::Always);

// Max performance (you can still validate explicitly when desired):
dt.set_validation_policy(ValidationPolicy::Never);
use delaunay::prelude::*;

let vertices = vec![
    vertex!([0.0, 0.0, 0.0]),
    vertex!([1.0, 0.0, 0.0]),
    vertex!([0.0, 1.0, 0.0]),
    vertex!([0.0, 0.0, 1.0]),
];
let dt = DelaunayTriangulation::new(&vertices).unwrap();
assert!(dt.tds().is_valid().is_ok());
assert!(dt.triangulation().is_valid().is_ok());
assert!(dt.is_valid().is_ok());
assert!(dt.validate().is_ok());

For implementation details on invariant enforcement, see core::algorithms::incremental_insertion.

§Correctness Guarantees and Limitations

The library provides strong correctness guarantees for vertex insertion operations while being transparent about edge cases and limitations.

§Guarantees

When using DelaunayTriangulation::insert() or the underlying insertion algorithms:

  1. Successful insertions are designed to maintain all invariants - If insertion succeeds (Ok(_)), the triangulation is expected to satisfy all structural and topological invariants. The incremental cavity-based insertion algorithm is designed to maintain these invariants. For applications requiring strict guarantees, use DelaunayTriangulation::validate() (Levels 1–4) or DelaunayTriangulation::is_valid() (Level 4 only) to verify the Delaunay property.

  2. Failed insertions leave triangulation in valid state - If insertion fails (Err(_)), the triangulation remains in a valid state with all invariants maintained. No partial or corrupted state is possible.

  3. Clear error messages - Insertion failures include detailed error messages specifying which constraint or invariant was violated, along with context about what went wrong.

  4. No silent failures - The library never silently produces incorrect triangulations. Operations either succeed with guarantees or fail with explicit errors.

  5. Duplicate vertex detection - Duplicate and near-duplicate vertices (within 1e-10 epsilon) are automatically detected and rejected with InsertionError::DuplicateCoordinates or InsertionError::DuplicateUuid, preventing numerical instabilities.

When constructing a triangulation from a batch of vertices using DelaunayTriangulation::new:

  • Successful construction yields a triangulation that is designed to satisfy the Delaunay property. Use dt.validate() (Levels 1–4) for cumulative verification.
  • Duplicate coordinates are automatically detected and rejected.

Incremental construction via DelaunayTriangulation::insert follows the same invariant rules on each insertion: on success the triangulation remains structurally valid; on failure the data structure is rolled back to its previous state. Use DelaunayTriangulation::is_valid() (Level 4) if you need explicit verification of the Delaunay property.

§Incremental insertion algorithm

Triangulations are built using an efficient incremental cavity-based insertion algorithm:

  • Initial simplex construction - The first D+1 affinely independent vertices are used to create an initial valid simplex using robust orientation predicates. If no non-degenerate simplex can be formed, construction fails with TriangulationConstructionError::GeometricDegeneracy.

  • Incremental insertion - Each subsequent vertex is inserted using a cavity-based algorithm that:

    1. Locates the vertex using efficient point location
    2. Identifies conflicting cells (those whose circumsphere contains the new vertex)
    3. Removes conflicting cells to create a cavity
    4. Fills the cavity with new cells connecting the cavity boundary to the new vertex
    5. Wires neighbor relationships locally without global recomputation

The incremental insertion algorithm maintains all structural invariants throughout construction, and is designed to satisfy the Delaunay (empty-circumsphere) property. Vertices are only rejected if they would violate fundamental geometric constraints (duplicates, near-duplicates, or degenerate configurations).

§Delaunay validation

The incremental insertion algorithm is designed to maintain the Delaunay property, aiming to ensure that the empty circumsphere property holds after each insertion. Global Delaunay validation can be performed explicitly using DelaunayTriangulation::is_valid when verification is needed (see Issue #120 for rare edge cases where validation may be necessary).

For construction from a batch of vertices using DelaunayTriangulation::new, the resulting triangulation is constructed to satisfy the Delaunay property. Call DelaunayTriangulation::is_valid() if you need explicit verification.

§Error handling

The incremental insertion algorithm provides clear error reporting for vertices that cannot be inserted:

  • Duplicate detection - Exact and near-duplicate vertices are detected and rejected with InsertionError::DuplicateCoordinates or InsertionError::DuplicateUuid
  • Geometric failures - Degenerate configurations that would violate the Delaunay property are rejected with appropriate error messages
  • Validation failures - If insertion would break structural invariants, the operation fails and the triangulation is left in its previous valid state
use delaunay::prelude::*;

let vertices = vec![
    vertex!([0.0, 0.0, 0.0]),
    vertex!([1.0, 0.0, 0.0]),
    vertex!([0.0, 1.0, 0.0]),
    vertex!([0.0, 0.0, 1.0]),
];

let dt: DelaunayTriangulation<_, (), (), 3> =
    DelaunayTriangulation::new(&vertices).unwrap();

assert_eq!(dt.number_of_vertices(), 4);
assert!(dt.validate().is_ok());

§Degenerate input handling

When the input vertices cannot form a non-degenerate simplex (for example, when all points are collinear in 2D), construction fails during initial simplex construction with TriangulationConstructionError::GeometricDegeneracy. This occurs because degenerate simplices (collinear in 2D, coplanar in 3D, etc.) are detected early using robust orientation predicates before any topology is built.

use delaunay::prelude::*;

// All points lie on a line in 2D: no non-degenerate simplex exists.
let degenerate = vec![
    vertex!([0.0, 0.0]),
    vertex!([1.0, 0.0]),
    vertex!([2.0, 0.0]),
    vertex!([3.0, 0.0]),
];

let result: Result<DelaunayTriangulation<_, (), (), 2>, _> =
    DelaunayTriangulation::new(&degenerate);

// Collinear points fail during initial simplex construction due to degeneracy
assert!(matches!(
    result,
    Err(DelaunayTriangulationConstructionError::Triangulation(
        TriangulationConstructionError::GeometricDegeneracy { .. },
    ))
));

§Limitations

  1. Degenerate geometry in higher dimensions - Highly degenerate point configurations (e.g., many nearly collinear or coplanar points) in 4D and 5D may cause insertion to fail gracefully with InsertionError. This is a known limitation of incremental algorithms in high-dimensional spaces with degenerate inputs.

  2. Iterative refinement constraints - The cavity-based insertion algorithm uses iterative refinement to maintain the Delaunay property. In rare cases with complex geometries, refinement may hit topological constraints and fail gracefully rather than producing an invalid triangulation.

  3. Numerical precision - Like all computational geometry libraries, numerical precision can affect results near floating-point boundaries. The library uses robust predicates to minimize these issues, but extreme coordinate values or ill-conditioned point sets may still cause problems.

§Simple API Usage

use delaunay::prelude::*;

// Create 4D triangulation - uses fast predicates by default (f64)
let vertices = 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 dt: DelaunayTriangulation<_, (), (), 4> =
    DelaunayTriangulation::new(&vertices).unwrap();

assert_eq!(dt.number_of_vertices(), 5);
assert_eq!(dt.dim(), 4);
assert_eq!(dt.number_of_cells(), 1);  // Single 4-simplex

// Also works in 2D
let vertices_2d = vec![
    vertex!([0.0, 0.0]),
    vertex!([1.0, 0.0]),
    vertex!([0.5, 1.0]),
];

let dt_2d: DelaunayTriangulation<_, (), (), 2> =
    DelaunayTriangulation::new(&vertices_2d).unwrap();

assert_eq!(dt_2d.number_of_vertices(), 3);
assert_eq!(dt_2d.dim(), 2);

For implementation details on invariant validation and error handling, see core::algorithms::incremental_insertion.

§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 core module contains the primary data structures and algorithms for building and manipulating Delaunay triangulations.
geometry
Contains geometric types including the Point struct 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.
topology
Topology analysis and validation for triangulated spaces.

Macros§

assert_jaccard_gte
Assert that the Jaccard index between two sets meets or exceeds a threshold.
vertex
Convenience macro for creating vertices with less boilerplate.

Functions§

is_normal
The function is_normal checks that structs implement auto traits. Traits are checked at compile time, so this function is only used for testing.