delaunay 0.7.6

D-dimensional Delaunay triangulations and convex hulls in Rust, with exact predicates, multi-level validation, and bistellar flips
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
# API Design: Two-Track Approach

This document explains the dual API design for working with Delaunay triangulations:
the **Builder API** for constructing and maintaining Delaunay triangulations, and
the **Edit API** for explicit topology editing via bistellar flips.

## Overview

The library provides two distinct APIs for different use cases:

1. **Builder API** (`DelaunayTriangulation::insert` / `::remove_vertex`)
   - High-level construction and maintenance of Delaunay triangulations
   - Automatically maintains the Delaunay property (empty circumsphere)
   - Designed for building triangulations from point sets
   - Uses cavity-based insertion and fan retriangulation

2. **Edit API** (`prelude::triangulation::flips::BistellarFlips` trait)
   - Low-level topology editing via bistellar (Pachner) flips
   - Explicit control over individual topology operations
   - Does **not** automatically restore the Delaunay property
   - Designed for topology manipulation, research, and custom algorithms

## When to Use Each API

### Use the Builder API when

- Building a Delaunay triangulation from a set of points
- Adding/removing vertices while maintaining the Delaunay property
- You need automatic geometric property preservation
- Working with standard computational geometry workflows

**Example use cases:**

- Computing convex hulls
- Nearest-neighbor queries
- Voronoi diagram construction
- Mesh generation
- Scientific simulations requiring Delaunay meshes

### Use the Edit API when

- Implementing custom topological algorithms
- Researching bistellar flip sequences
- Building non-Delaunay triangulations with specific properties
- Experimenting with topology transformations
- You need explicit control over topology changes

**Example use cases:**

- Implementing custom Delaunay repair strategies
- Topology optimization algorithms
- Research on triangulation properties
- Building triangulations with non-Delaunay constraints
- Educational demonstrations of bistellar flip theory

## Builder API Reference

### Simple Construction: `DelaunayTriangulation::new()`

For most use cases, the simple constructor is sufficient:

```rust
use delaunay::prelude::triangulation::*;

// Simple construction from vertices (Euclidean space, default options)
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<_, (), (), 3> = 
    DelaunayTriangulation::new(&vertices).unwrap();

// Incremental insertion (maintains Delaunay property)
let new_vertex = vertex!([0.5, 0.5, 0.5]);
dt.insert(new_vertex).unwrap();

// Vertex removal (maintains Delaunay property via fan retriangulation)
let vertex_key = dt.vertices().next().unwrap().0;
dt.remove_vertex(vertex_key).unwrap();
```

### Advanced Construction: `DelaunayTriangulationBuilder`

For advanced configuration (toroidal topology, custom validation policies, etc.),
use `DelaunayTriangulationBuilder`:

```rust
use delaunay::prelude::triangulation::*;
use delaunay::core::triangulation::{TopologyGuarantee, ValidationPolicy};

// Toroidal (periodic) triangulation in 2D
let vertices = vec![
    vertex!([0.1, 0.1]),
    vertex!([0.9, 0.9]),
    vertex!([0.5, 0.5]),
];

let mut dt = DelaunayTriangulationBuilder::new(&vertices)
    .toroidal([1.0, 1.0]) // Phase 1: canonicalized toroidal construction
    .topology_guarantee(TopologyGuarantee::PLManifoldStrict)
    .build::<()>()
    .unwrap();

dt.set_validation_policy(ValidationPolicy::Always);

// Works like any other DelaunayTriangulation
dt.insert(vertex!([0.25, 0.75])).unwrap();
```

**When to use the Builder:**

- **Toroidal/periodic triangulations**: Use `.toroidal()` (Phase 1 canonicalized) or
  `.toroidal_periodic()` (Phase 2 periodic image-point) with explicit domain periods
- **Custom topology guarantees**: Set stricter or more relaxed manifold checks
- **Custom validation policies**: Configure `ValidationPolicy` via
  `dt.set_validation_policy(...)` before or after build; the active policy controls
  automatic topology validation for subsequent insert/remove and other modification
  operations
- **Custom repair policies**: Configure Delaunay repair behavior

See `docs/topology.md` for more on toroidal triangulations and `docs/validation.md`
for topology guarantee and validation policy details.

### Key Characteristics

- **Automatic property preservation**: The Delaunay empty-circumsphere property is maintained automatically
- **Cavity-based insertion**: New vertices are inserted by identifying conflicting cells, removing them, and filling the cavity
- **Fan retriangulation**: Vertex removal uses fan-based retriangulation of the vertex star
- **Auxiliary data**: Vertices and cells carry optional user data (`U` / `V`). Read via `vertex.data()` /
  `cell.data()`, write via `dt.set_vertex_data(key, data)` / `dt.set_cell_data(key, data)` (O(1),
  invariant-preserving). See [`workflows.md`]workflows.md for examples.
- **Error handling**: Operations fail gracefully if they would violate invariants (see
  [`invariants.md`]invariants.md).
- **Validation**: The active `ValidationPolicy` (set with
  `dt.set_validation_policy(...)`) governs automatic topology validation for
  subsequent construction/modification operations

## Edit API Reference

The Edit API is exposed through the `BistellarFlips` trait in `prelude::triangulation::flips`:

```rust
use delaunay::prelude::triangulation::*;
use delaunay::prelude::triangulation::flips::*;

// Start with a valid 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<_, (), (), 3> = 
    DelaunayTriangulation::new(&vertices).unwrap();

// k=1 move: Insert a vertex into a cell (splits cell into D+1 cells)
let cell_key = dt.cells().next().unwrap().0;
let info = dt.flip_k1_insert(cell_key, vertex!([0.25, 0.25, 0.25])).unwrap();

// k=1 inverse: Remove a vertex (collapses its star)
let vertex_key = info.inserted_face_vertices[0];
dt.flip_k1_remove(vertex_key).unwrap();

// k=2 move: Flip a facet (2 cells ↔ D cells)
let facet = /* FacetHandle */;
let info = dt.flip_k2(facet).unwrap();

// k=2 inverse: Flip from an edge star (D cells ↔ 2 cells)
let edge = EdgeKey::new(info.inserted_face_vertices[0], info.inserted_face_vertices[1]);
dt.flip_k2_inverse_from_edge(edge).unwrap();

// k=3 move: Flip a ridge (3 cells ↔ D-1 cells, requires D ≥ 3)
let ridge = /* RidgeHandle */;
let info = dt.flip_k3(ridge).unwrap();

// k=3 inverse: Flip from a triangle star (D-1 cells ↔ 3 cells)
let triangle = TriangleHandle::new(
    info.inserted_face_vertices[0],
    info.inserted_face_vertices[1],
    info.inserted_face_vertices[2],
);
dt.flip_k3_inverse_from_triangle(triangle).unwrap();
```

### Available Flip Operations

#### k=1 Moves (Cell Split/Merge)

- **Forward (`flip_k1_insert`)**: Insert a vertex into a cell, splitting it into D+1 cells
  - Valid for D ≥ 1
  - Replaces 1 cell with D+1 cells
  - Removed face: the entire cell (D-simplex)
  - Inserted face: the new vertex (0-simplex)

- **Inverse (`flip_k1_remove`)**: Remove a vertex, collapsing its star
  - Requires the vertex star to be collapsible (star of D+1 cells forming a ball)
  - Replaces D+1 cells with 1 cell

#### k=2 Moves (Facet Flip)

- **Forward (`flip_k2`)**: Flip a facet shared by 2 cells
  - Valid for D ≥ 2
  - Replaces 2 cells with D cells
  - Removed face: the shared facet ((D-1)-simplex)
  - Inserted face: an edge (1-simplex)

- **Inverse (`flip_k2_inverse_from_edge`)**: Flip from an edge star
  - Requires an edge with star of D cells
  - Replaces D cells with 2 cells

#### k=3 Moves (Ridge Flip)

- **Forward (`flip_k3`)**: Flip a ridge
  - Valid for D ≥ 3
  - Replaces 3 cells with D-1 cells
  - Removed face: a ridge ((D-2)-simplex)
  - Inserted face: a triangle (2-simplex)

- **Inverse (`flip_k3_inverse_from_triangle`)**: Flip from a triangle star
  - Requires a triangle with star of D-1 cells
  - Replaces D-1 cells with 3 cells

### Key Characteristics

- **Explicit control**: You specify exactly which flip to perform
- **No automatic property preservation**: The Delaunay property is **not** maintained automatically
- **Reversible**: Each forward move has a corresponding inverse
- **Geometric validation**: Flips check for degeneracy and manifold preservation
- **Flexible**: Can be used to build custom repair or optimization algorithms

### Important Caveats

⚠️ **The Edit API does not preserve the Delaunay property automatically.**

After applying flips, you should:

1. Manually verify the Delaunay property if needed:

   ```rust
   dt.is_valid().unwrap();  // Check Level 4 (Delaunay property)
   ```

2. Consider running a repair pass if you need the Delaunay property again (requires `K: ExactPredicates`):
   - `dt.repair_delaunay_with_flips()` (flip-based repair)
   - `dt.repair_delaunay_with_flips_advanced(DelaunayRepairHeuristicConfig::default())` (includes a heuristic rebuild fallback)

## Combining Both APIs

You can mix both APIs in the same workflow:

```rust
use delaunay::prelude::triangulation::*;
use delaunay::prelude::triangulation::flips::*;

// 1. Build initial triangulation (Builder API)
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<_, (), (), 3> = 
    DelaunayTriangulation::new(&vertices).unwrap();

// 2. Add vertices using Builder API (maintains Delaunay)
dt.insert(vertex!([0.5, 0.5, 0.5])).unwrap();

// 3. Make custom topology edits (Edit API)
let facet = /* ... */;
dt.flip_k2(facet).unwrap();

// 4. Verify Delaunay property if needed
if let Err(e) = dt.is_valid() {
    eprintln!("Warning: Delaunay property violated after manual edit: {}", e);
    // Optionally restore using Builder API or custom repair
}
```

## Validation and Guarantees

Both APIs work with the same validation framework but have different guarantees:

### Builder API Guarantees

- ✅ Maintains **structural invariants** (Level 1-2)
- ✅ Maintains **manifold topology** (Level 3, controlled by `TopologyGuarantee`)
- ✅ Designed to maintain **Delaunay property** (Level 4)
- ✅ Fails gracefully if invariants cannot be maintained

### Edit API Guarantees

- ✅ Maintains **structural invariants** (Level 1-2)
- ✅ Checks **geometric degeneracy** (prevents degenerate flips)
- ⚠️ Does **not** automatically maintain Delaunay property
- ⚠️ User is responsible for property preservation

### Validation Levels

Use the appropriate validation level for your needs:

```rust
// Level 2: Structural only (fast)
dt.tds().is_valid().unwrap();

// Level 3: + Manifold topology
dt.as_triangulation().is_valid().unwrap();

// Level 4: + Delaunay property (most comprehensive)
dt.is_valid().unwrap();

// Full diagnostic report
let report = dt.validation_report();
```

## Implementation Details

### Internal Organization

- **Builder API**: Implemented in `triangulation::delaunay`, `triangulation::builder`,
  and `core::algorithms::incremental_insertion`
- **Edit API**: Implemented in `triangulation::flips` (public trait) and `core::algorithms::flips` (internal implementation)
- **Low-level primitives**: Context builders and flip application functions are `pub(crate)` in `core::algorithms::flips`

### Design Rationale

The separation serves several purposes:

1. **Clear contracts**: Builder API guarantees Delaunay property; Edit API does not
2. **Safety**: Low-level flip primitives are not exposed to prevent accidental misuse
3. **Flexibility**: Edit API enables research and custom algorithms without restricting the design
4. **Documentation**: Clear distinction between "construction" and "manipulation" workflows

## Examples

See the following examples for practical demonstrations:

- `examples/topology_editing_2d_3d.rs` - 2D+3D example showing both APIs
- `examples/pachner_roundtrip_4d.rs` - Advanced 4D example with all flip types
- `examples/triangulation_3d_100_points.rs` - Builder API usage for construction
- `examples/delaunayize_repair.rs` - Delaunayize workflow (2D/3D/4D, flip-then-repair, custom config)

## Delaunayize Workflow

The `triangulation::delaunayize` module provides a single entrypoint for the
common "repair topology then restore Delaunay" workflow:

```rust
use delaunay::prelude::triangulation::delaunayize::*;

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<_, (), (), 3> =
    DelaunayTriangulation::new(&vertices).unwrap();

let outcome = delaunayize_by_flips(&mut dt, DelaunayizeConfig::default()).unwrap();
assert!(outcome.topology_repair.succeeded);
```

### Steps

1. **PL-manifold topology repair** — bounded deterministic removal of cells
   that cause facet over-sharing (codimension-1 facet degree > 2).
2. **Delaunay flip repair** — k=2/k=3 bistellar flips to restore the
   empty-circumsphere property.
3. **Optional fallback rebuild** — rebuilds from the vertex set when both
   repair passes fail (`DelaunayizeConfig { fallback_rebuild: true, .. }`).

### Configuration

`DelaunayizeConfig` controls:

- `topology_max_iterations` (default 64): max repair iterations.
- `topology_max_cells_removed` (default 10,000): max cells removed.
- `fallback_rebuild` (default false): rebuild from vertices on failure,
  restoring cell data for rebuilt cells whose sorted vertex UUID set still
  matches exactly one original cell.
- `delaunay_max_flips` (default `None`): optional per-attempt flip budget.

### Data Preservation

`PlManifoldRepairStats` carries `removed_cells` and `removed_vertices`
(identified by UUID) so callers can recover user data from entities removed
during topology repair. The fallback rebuild path also preserves cell payloads
when a rebuilt cell has the same vertex UUID set as exactly one original cell;
changed or ambiguous cells receive no payload.

### Explicitly Deferred

- Dedicated targeted repair stages for boundary-ridge multiplicity,
  ridge-link manifoldness, and vertex-link manifoldness (#304).

## Further Reading

- **Bistellar flip theory**: See `REFERENCES.md` for academic citations
- **Validation framework**: See `docs/validation.md` for detailed validation guide
- **Invariant rationale**: See [`invariants.md`]invariants.md for theory and implementation pointers
- **Topology analysis**: See `docs/topology.md` for topological concepts
- **API implementation**: See `triangulation::flips` module documentation