Skip to main content

Crate delaunay

Crate delaunay 

Source
Expand description

Β§delaunay

DOI Crates.io Downloads License Docs.rs CI rust-clippy analyze codecov Audit dependencies Codacy Badge

D-dimensional Delaunay triangulations and convex hulls in Rust, with exact predicates, multi-level validation, and bistellar flips inspired by CGAL.

Β§πŸ“ Introduction

This library implements d-dimensional Delaunay triangulations in Rust for computational geometry and scientific workflows. It is inspired by CGAL, which is a C++ library for computational geometry, and Spade, a Rust library that implements 2D Delaunay triangulations, Constrained Delaunay triangulations, and Voronoi diagrams. The goal is to provide a lightweight but rigorous alternative to CGAL in the Rust ecosystem, with explicit topology settings, validation levels, and repair behavior.

Β§πŸ§ͺ Scientific Basis

This crate models triangulations of finite point sets in R^d as oriented simplicial complexes with explicit combinatorial and geometric checks.

  • Core operational invariant for editing/repair: coherent orientation + PL-manifold validity
  • Local move system: bistellar flips (k = 1, 2, 3 and inverses), providing the supported Pachner moves set in dimensions ≀ 5
  • Geometric convergence basis in finite-point workflows: Herbert Edelsbrunner and Nina R. Shah, Incremental Topological Flipping Works for Regular Triangulations, Discrete & Computational Geometry (1996), https://doi.org/10.1007/BF01975867
  • Scope of claims: guarantees apply to supported library workflows (construction, flip-based repair, and validation APIs), not arbitrary external mutation of internal structures

§✨ Features

  • Copyable data types associated with vertices and cells (integers, floats, chars, custom enums)
  • d-dimensional Delaunay triangulations
  • d-dimensional Convex hulls
  • Toroidal (periodic) triangulations via DelaunayTriangulationBuilder with Phase 1 (canonicalization) and Phase 2 (image-point method) support
  • Geometry quality metrics for simplices: radius ratio and normalized volume (dimension-agnostic)
  • Serialization/deserialization of all data structures to/from JSON
  • Tested for 2-, 3-, 4-, and 5-dimensional triangulations
  • Configurable predicate kernels: AdaptiveKernel (default; exact arithmetic + SoS), RobustKernel (exact, preserves degeneracy signals), FastKernel (raw f64, 2D only)
  • Bulk insertion ordering (InsertionOrderStrategy): Hilbert curve (default) or input order
  • Batch construction options (ConstructionOptions): optional deduplication and deterministic retries
  • Incremental construction APIs: insertion plus vertex removal (remove_vertex)
  • 4-level validation hierarchy (element validity β†’ TDS structural validity β†’ manifold topology β†’ Delaunay property), including full diagnostics via validation_report
  • Local topology validation (PL-manifold default, Pseudomanifold opt-out)
  • Coherent combinatorial orientation validation/normalization for cells, maintaining oriented simplicial complexes
  • The complete set of Pachner moves up to 5D implemented as bistellar k-flips for k = 1, 2, 3 plus inverse moves
  • Delaunay repair using bistellar flips for k=2/k=3 with inverse edge/triangle queues in 4D/5D
  • Safe Rust: #![forbid(unsafe_code)]

See CHANGELOG.md for details. Older releases are archived under docs/archive/changelog/.

§🟒 Minimal Construction Example

The construction API has two entry points:

  • DelaunayTriangulation::new(&vertices) - simple constructor for the common case
  • DelaunayTriangulationBuilder - Advanced configuration (custom options, toroidal topology)

Add the library to your crate:

cargo add delaunay
use delaunay::prelude::triangulation::*;

// Create a 4D Delaunay triangulation from a set of vertices (uses AdaptiveKernel by default).
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]), // Adding this vertex creates new simplices via flips
];

let dt = DelaunayTriangulation::new(&vertices).unwrap();

assert_eq!(dt.dim(), 4);
assert_eq!(dt.number_of_vertices(), 6);
assert_eq!(dt.number_of_cells(), 5); // 1 simplex from first 5 vertices + 4 new simplices from last vertex

// Optional verification:
// - `dt.is_valid()` checks Level 4 only (Delaunay property).
// - `dt.validate()` checks Levels 1–4 (elements + structure + topology + Delaunay).
assert!(dt.is_valid().is_ok());

Β§Toroidal (Periodic) Triangulations

For periodic boundary conditions, use DelaunayTriangulationBuilder:

use delaunay::prelude::triangulation::*;

// Phase 1: Canonicalization (wraps coordinates into [0, 1)Β²)
let vertices = vec![
    vertex!([0.1, 0.2]),
    vertex!([0.8, 0.3]),
    vertex!([0.5, 0.7]),
    vertex!([1.2, 0.4]),  // Wraps to [0.2, 0.4]
];

let dt = DelaunayTriangulationBuilder::new(&vertices)
    .toroidal([1.0, 1.0])  // Fundamental domain periods
    .build::<()>()
    .unwrap();

assert_eq!(dt.topology_kind(), TopologyKind::Toroidal);

For the full periodic image-point method (Phase 2), see the DelaunayTriangulationBuilder documentation.

Β§Need more control?

Β§βœ… Validation and Guarantees

LevelWhat is validatedPrimary API
1Element validity (vertex/cell primitives)dt.validate() / dt.validation_report()
2TDS structural validity (keys, incidences, neighbors)dt.tds().is_valid()
3Manifold topology (link checks, Euler/topological consistency)dt.as_triangulation().is_valid()
4Delaunay property (empty-circumsphere via local predicates)dt.is_valid()
1-4Cumulative checks with diagnosticsdt.validate() or dt.validation_report()

TopologyGuarantee controls which Level 3 manifold constraints are enforced, and ValidationPolicy controls when Level 3 checks run automatically during incremental insertion.

Β§πŸ”¬ Reproducibility

The construction pipeline exposes deterministic controls for experiments and regression testing:

  • Deterministic insertion ordering via InsertionOrderStrategy: Hilbert (default) or Input (use Input to preserve caller-provided order exactly)
  • Deterministic preprocessing via DedupPolicy
  • Deterministic retry behavior via RetryPolicy (including seeded shuffled retries) or RetryPolicy::Disabled
  • Explicit topology/validation configuration via TopologyGuarantee and ValidationPolicy
use delaunay::core::delaunay_triangulation::{
    ConstructionOptions, DedupPolicy, InsertionOrderStrategy, RetryPolicy,
};
use delaunay::core::triangulation::{TopologyGuarantee, ValidationPolicy};
use delaunay::prelude::triangulation::*;

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

let options = ConstructionOptions::default()
    .with_insertion_order(InsertionOrderStrategy::Input)
    .with_dedup_policy(DedupPolicy::Exact)
    .with_retry_policy(RetryPolicy::Disabled);

let mut dt = DelaunayTriangulationBuilder::new(&vertices)
    .topology_guarantee(TopologyGuarantee::PLManifold)
    .construction_options(options)
    .build::<()>()
    .unwrap();

dt.set_validation_policy(ValidationPolicy::Always);
assert!(dt.validate().is_ok());

For reproducible checks in CI/local runs, use just check, just test, just doc-check, or just ci.

§⚠️ Limitations

  • Dimension coverage: CI and property-test coverage targets 2D–5D.
  • Exact predicates (D ≀ 5): orientation and insphere predicates use provably correct exact arithmetic (det_sign_exact) for D ≀ 5. For D β‰₯ 6, insphere falls back to symbolic perturbation (the (D+2)Γ—(D+2) matrix exceeds the stack-allocation limit).
  • Flip convergence: at larger scales (~130+ points in 3D, ~40+ in 4D), flip-based Delaunay repair can encounter cycles. Since v0.7.3, Simulation of Simplicity (SoS) via AdaptiveKernel breaks predicate-degeneracy ties deterministically; remaining convergence issues at scale are typically cavity/topology interactions rather than predicate ambiguity.
  • 4D+ bulk construction: see Known Issues for details.
  • Validation/repair guarantees assume the library-managed construction/editing pipeline.

§🚧 Project History

This crate was originally maintained at https://github.com/oovm/shape-rs through version 0.1.0. The original implementation provided basic Delaunay triangulation functionality.

Starting with version 0.3.4, maintenance transferred to this repository, which hosts a completely rewritten d-dimensional implementation focused on computational geometry research applications.

§🀝 How to Contribute

We welcome contributions! Here’s a quickstart:

# Clone and setup
git clone https://github.com/acgetchell/delaunay.git
cd delaunay

# Setup development environment (installs tools, builds project)
cargo install just
just setup            # Installs all development tools and dependencies

# Development workflow
just fix              # Apply formatters/auto-fixes (mutating)
just check            # Lint/validators (non-mutating)
just ci               # Full CI run (checks + all tests + examples + bench compile)
just ci-slow          # CI + slow tests (100+ vertices)
just --list           # See all available commands
just help-workflows   # Show common workflow patterns

Try the examples:

just examples         # Run all examples
# Or run specific examples:
cargo run --release --example triangulation_3d_100_points
cargo run --release --example convex_hull_3d_100_points

Β§πŸ“‹ Examples

The examples/ directory contains several demonstrations:

  • convex_hull_3d_100_points: 3D convex hull extraction and analysis on the same 100-point configuration
  • into_from_conversions: Demonstrates Into/From trait conversions and utilities
  • memory_analysis: Memory usage analysis for triangulations across dimensions with allocation tracking
  • pachner_roundtrip_4d: 4D Pachner move (k=1,2,3) roundtrip checks (flip + inverse preserves the triangulation)
  • point_comparison_and_hashing: Demonstrates point comparison and hashing behavior
  • topology_editing_2d_3d: Builder API vs Edit API in 2D/3D (bistellar flips and Delaunay preservation)
  • triangulation_3d_100_points: 3D Delaunay triangulation with a stable 100-point random configuration
  • zero_allocation_iterator_demo: Performance comparison between allocation and zero-allocation iterators

For detailed documentation, sample output, and usage instructions for each example, see examples/README.md.

For comprehensive guidelines on development environment setup, testing, benchmarking, performance analysis, and development workflow, please see CONTRIBUTING.md.

This includes information about:

  • Building and testing the library
  • Running benchmarks and performance analysis
  • Code style and standards
  • Submitting changes and pull requests
  • Project structure and development tools

Β§πŸ“– Documentation

  • API Design - Builder vs Edit API design (explicit bistellar flips)
  • Code Organization - Project structure and module patterns
  • Invariants - Theoretical background and rationale for the topological and geometric invariants
  • Numerical Robustness Guide - Robustness strategies, kernels, and retry/repair behavior
  • Property Testing Summary - Property-based testing with proptest (where tests live, how to run)
  • Known Issues - Dimensional and large-scale limitations (exact predicate bounds, flip convergence, 4D+)
  • Releasing - Release workflow (changelog + benchmarks + publish)
  • Topology - Level 3 topology validation (manifoldness + Euler characteristic) and module overview
  • Validation Guide - Comprehensive 4-level validation hierarchy guide (element β†’ structural β†’ manifold β†’ Delaunay)
  • Workflows - Happy-path construction plus practical Builder/Edit recipes (stats, repairs, and minimal flips)

Β§πŸ“Ž How to Cite

If you use this software in academic work, cite the Zenodo DOI and include the software metadata from CITATION.cff.

@software{getchell_delaunay,
  author = {Adam Getchell},
  title = {delaunay: A d-dimensional Delaunay triangulation library},
  doi = {10.5281/zenodo.16931097},
  url = {https://github.com/acgetchell/delaunay}
}

For release-specific fields (version, release date, ORCID), prefer CITATION.cff.

Β§πŸ“š References

For a comprehensive list of academic references and bibliographic citations used throughout the library, see REFERENCES.md.

Β§πŸ€– AI Agents

This repository contains an AGENTS.md file, which defines the canonical rules and invariants for all AI coding assistants and autonomous agents working on this codebase.

AI tools (including ChatGPT, Claude, GitHub Copilot, Cursor, Warp, and CI repair agents) are expected to read and follow AGENTS.md when proposing or applying changes.

Portions of this library were developed with the assistance of these AI tools:

All code was written and/or reviewed and validated by the author.


Β§Documentation map

The README above is included verbatim and serves as the user-facing introduction to the crate (overview, features, and quick-start examples).

Everything below this line specifies the semantic and correctness contract of the delaunay crate and is intended for users who need stronger guarantees, deeper understanding of invariants, or who are extending the implementation.

This crate’s documentation is intentionally layered by audience and intent:

  • README.md (included above): User-facing overview, feature list, and quick-start examples.

  • Crate-level documentation (lib.rs) (this document): The programming contract of the library: what invariants are enforced, when validation runs, and what errors mean.

    In particular, this document covers:

    • The validation hierarchy and invariant stack (Levels 1–4)
    • Topological guarantees (TopologyGuarantee) and insertion-time validation policy (ValidationPolicy)
    • High-level error semantics and programming contract (transactional operations, duplicate rejection)
  • docs/workflows.md: Task-oriented, end-to-end usage recipes (Builder API, Edit API, validation, repairs, diagnostics, and statistics).

  • docs/validation.md: Formal definitions of validation Levels 1–4, their costs, and guidance on when each level should be applied.

  • docs/invariants.md: Deeper theoretical discussion of topological and geometric invariants (PL-manifold conditions, ridge/vertex links, ordering heuristics, and convergence assumptions), plus algorithmic background and limitations.

Β§Examples (contract-oriented)

Β§Validation hierarchy (Levels 1–4)

use delaunay::prelude::triangulation::*;

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 = DelaunayTriangulationBuilder::new(&vertices).build::<()>().unwrap();

// Levels 1–2: elements + structural (TDS)
assert!(dt.tds().validate().is_ok());

// Levels 1–3: elements + structural + topology
assert!(dt.as_triangulation().validate().is_ok());

// Level 4 only: Delaunay property (assumes Levels 1–3)
assert!(dt.is_valid().is_ok());

// Levels 1–4: full cumulative validation
assert!(dt.validate().is_ok());

Β§Topology guarantees and insertion-time validation (TopologyGuarantee, ValidationPolicy)

use delaunay::prelude::triangulation::*;

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 mut dt = DelaunayTriangulationBuilder::new(&vertices).build::<()>().unwrap();

assert_eq!(dt.topology_guarantee(), TopologyGuarantee::PLManifold);
assert_eq!(dt.validation_policy(), ValidationPolicy::OnSuspicion);

dt.set_topology_guarantee(TopologyGuarantee::Pseudomanifold);
dt.set_validation_policy(ValidationPolicy::Always);

assert_eq!(dt.topology_guarantee(), TopologyGuarantee::Pseudomanifold);
assert_eq!(dt.validation_policy(), ValidationPolicy::Always);

Β§Transactional operations and duplicate rejection

use delaunay::prelude::triangulation::*;

let vertices = vec![
    vertex!([0.0, 0.0]),
    vertex!([1.0, 0.0]),
    vertex!([0.0, 1.0]),
];
let mut dt = DelaunayTriangulationBuilder::new(&vertices).build::<()>().unwrap();

let before_vertices = dt.number_of_vertices();
let before_cells = dt.number_of_cells();

// Duplicate coordinates are rejected.
let result = dt.insert(vertex!([0.0, 0.0]));
assert!(matches!(result, Err(InsertionError::DuplicateCoordinates { .. })));

// On error, the triangulation is unchanged.
assert_eq!(dt.number_of_vertices(), before_vertices);
assert_eq!(dt.number_of_cells(), before_cells);

Β§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:

  • Vertex and Cell provide element validity checks. Level 1 (elements) validation checks invariants such as:

    • Vertex coordinates – finite (no NaN/∞) and UUID is non-nil.
    • Cell shape – exactly D+1 distinct vertex keys, valid UUID, and neighbor buffer length (if present) is D+1.

    These checks are surfaced via Vertex::is_valid and Cell::is_valid, and are automatically run by Tds::validate (Levels 1–2).

  • 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 1 (elements / Vertex + Cell): Vertex::is_valid() / Cell::is_valid() for element checks, or dt.tds().validate() for Levels 1–2.
  • 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.as_triangulation().is_valid() for topology-only checks, or dt.as_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::triangulation::*;

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 mut dt = DelaunayTriangulation::new(&vertices).unwrap();

// Performance mode: disable insertion-time Level 3 topology validation.
dt.set_validation_policy(ValidationPolicy::Never);

// Do incremental work...
dt.insert(vertex!([0.2, 0.2, 0.2])).unwrap();

// ...then explicitly validate the topology layer when you need a certificate.
assert!(dt.as_triangulation().validate().is_ok());

Β§Choosing Level 3 topology guarantee (TopologyGuarantee)

This section specifies what invariants are enforced. The formal topological definitions and rationale live in docs/invariants.md.

Level 3 topology validation is parameterized by TopologyGuarantee. This is separate from ValidationPolicy: it controls what invariants Level 3 enforces, not when automatic validation runs.

  • TopologyGuarantee::PLManifold (default): enforces manifold facet degree, boundary closure, connectedness, Euler characteristic, and link-based manifold conditions. Ridge-link checks are applied incrementally during insertion, with vertex-link validation performed at construction completion.

    The formal topological definitions, link conditions, and rationale for this validation strategy are documented in docs/invariants.md.

  • TopologyGuarantee::PLManifoldStrict: vertex-link validation after every insertion (slowest, maximum safety).

  • TopologyGuarantee::Pseudomanifold: skips vertex-link validation (may be faster), but bistellar flip convergence is not guaranteed and you may want to validate the Delaunay property explicitly for near-degenerate inputs.

use delaunay::prelude::triangulation::*;

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 = DelaunayTriangulationBuilder::new(&vertices).build::<()>().unwrap();

// For `TopologyGuarantee::PLManifold`, full certification includes a completion-time
// vertex-link validation pass.
assert!(dt.as_triangulation().validate_at_completion().is_ok());
use delaunay::prelude::triangulation::*;

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 = DelaunayTriangulationBuilder::new(&vertices).build::<()>().unwrap();

// `validate()` returns the first violation; `validation_report()` is intended for
// debugging/telemetry where you want the full set of violated invariants.
assert!(dt.validation_report().is_ok());

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

Β§Programming contract (high-level)

  • Transactional mutations: Construction and incremental operations are designed to be all-or-nothing. If an operation returns Err(_), the triangulation is rolled back to its previous state.
  • Duplicate detection: Near-duplicate coordinates are rejected using a small Euclidean tolerance (currently 1e-10 for the default floating-point scalars), returning InsertionError::DuplicateCoordinates. Duplicate UUIDs return InsertionError::DuplicateUuid.
  • Explicit verification: Use dt.validate() for cumulative verification (Levels 1–4), or dt.is_valid() for Level 4 only.

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.
triangulation
Triangulation-facing APIs.

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.