delaunay 0.7.7

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
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
//! Triangulation editing operations (bistellar flips).
//!
//! This module exposes **high-level** flip methods for explicit triangulation editing.
//! These operations do **not** automatically restore the Delaunay property.
//! For Delaunay construction/removal, use
//! [`crate::triangulation::delaunay::DelaunayTriangulation::insert`] and
//! [`crate::triangulation::delaunay::DelaunayTriangulation::remove_vertex`].

#![forbid(unsafe_code)]

pub use crate::core::algorithms::flips::{
    BistellarFlipKind, FlipContextError, FlipDirection, FlipEdgeAdjacencyError, FlipError,
    FlipInfo, FlipMutationError, FlipNeighborWiringError, FlipPredicateError,
    FlipPredicateOperation, FlipTriangleAdjacencyError, FlipVertexAdjacencyError, RidgeHandle,
    TriangleHandle,
};
pub use crate::core::edge::EdgeKey;
pub use crate::core::facet::FacetHandle;
pub use crate::core::tds::{CellKey, VertexKey};

use crate::core::algorithms::flips::{
    apply_bistellar_flip_dynamic, apply_bistellar_flip_k1, apply_bistellar_flip_k1_inverse,
    apply_bistellar_flip_k2, apply_bistellar_flip_k3, build_k2_flip_context,
    build_k2_flip_context_from_edge, build_k3_flip_context, build_k3_flip_context_from_triangle,
};
use crate::core::traits::data_type::DataType;
use crate::core::triangulation::Triangulation;
use crate::core::vertex::Vertex;
use crate::geometry::kernel::Kernel;
use crate::triangulation::delaunay::DelaunayTriangulation;
/// High-level triangulation editing operations via bistellar flips.
///
/// # Example
///
/// ```rust
/// use delaunay::prelude::triangulation::construction::TopologyGuarantee;
/// use delaunay::prelude::triangulation::flips::*;
///
/// 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_with_topology_guarantee(
///         &vertices,
///         TopologyGuarantee::PLManifold,
///     )
///     .unwrap();
/// let cell_key = dt.cells().next().unwrap().0;
///
/// // Split a cell by inserting a vertex (k=1 move).
/// let _info = dt
///     .flip_k1_insert(cell_key, vertex!([0.1, 0.1, 0.1]))
///     .unwrap();
/// ```
pub trait BistellarFlips<K, U, const D: usize>
where
    K: Kernel<D>,
    U: DataType,
{
    /// Apply a forward k=1 move (cell split) by inserting a vertex into a cell.
    ///
    /// # Errors
    ///
    /// Returns [`FlipError`] if the cell is missing, the vertex cannot be inserted,
    /// or the flip would create invalid topology.
    ///
    /// # Example
    ///
    /// ```rust
    /// use delaunay::prelude::triangulation::construction::TopologyGuarantee;
    /// use delaunay::prelude::triangulation::flips::*;
    ///
    /// 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_with_topology_guarantee(
    ///         &vertices,
    ///         TopologyGuarantee::PLManifold,
    ///     )
    ///     .unwrap();
    /// let cell_key = dt.cells().next().unwrap().0;
    ///
    /// // Insert a vertex into the cell
    /// let info = dt.flip_k1_insert(cell_key, vertex!([0.25, 0.25, 0.25])).unwrap();
    /// assert!(!info.new_cells.is_empty());
    /// ```
    fn flip_k1_insert(
        &mut self,
        cell_key: CellKey,
        vertex: Vertex<K::Scalar, U, D>,
    ) -> Result<FlipInfo<D>, FlipError>;

    /// Apply an inverse k=1 move (vertex collapse).
    ///
    /// # Errors
    ///
    /// Returns [`FlipError`] if the vertex star is not collapsible or the flip
    /// would create invalid topology.
    ///
    /// # Example
    ///
    /// ```rust
    /// use delaunay::prelude::triangulation::construction::TopologyGuarantee;
    /// use delaunay::prelude::triangulation::flips::*;
    ///
    /// 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_with_topology_guarantee(
    ///         &vertices,
    ///         TopologyGuarantee::PLManifold,
    ///     )
    ///     .unwrap();
    /// let cell_key = dt.cells().next().unwrap().0;
    /// let inserted = dt.flip_k1_insert(cell_key, vertex!([0.25, 0.25, 0.25])).unwrap();
    /// let inserted_vertex = inserted.inserted_face_vertices[0];
    ///
    /// // Remove the inserted vertex
    /// let info = dt.flip_k1_remove(inserted_vertex).unwrap();
    /// assert!(!info.removed_cells.is_empty());
    /// ```
    fn flip_k1_remove(&mut self, vertex_key: VertexKey) -> Result<FlipInfo<D>, FlipError>;

    /// Apply a k=2 facet flip (forward).
    ///
    /// # Errors
    ///
    /// Returns [`FlipError`] if the facet is invalid, the flip would be degenerate,
    /// or the resulting topology would be non-manifold.
    ///
    /// # Example
    ///
    /// ```rust
    /// use delaunay::prelude::triangulation::flips::*;
    ///
    /// 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]),
    ///     vertex!([0.5, 0.5, 0.3]),
    /// ];
    /// let mut dt: DelaunayTriangulation<_, (), (), 3> =
    ///     DelaunayTriangulation::new(&vertices).unwrap();
    ///
    /// // Find an interior facet and attempt a k=2 flip
    /// // Note: k=2 flips require specific geometric conditions
    /// let cell_key = dt.cells().next().map(|(k, _)| k);
    /// if let Some(key) = cell_key {
    ///     let has_neighbor = dt.tds().cell(key)
    ///         .and_then(|cell| cell.neighbors())
    ///         .map(|neighbors| neighbors.iter().any(|n| n.is_some()))
    ///         .unwrap_or(false);
    ///     if has_neighbor {
    ///         let facet = FacetHandle::new(key, 0);
    ///         let _ = dt.flip_k2(facet);  // May succeed or fail depending on configuration
    ///     }
    /// }
    /// ```
    fn flip_k2(&mut self, facet: FacetHandle) -> Result<FlipInfo<D>, FlipError>;

    /// Apply a k=3 ridge flip (forward).
    ///
    /// # Errors
    ///
    /// Returns [`FlipError`] if the ridge is invalid, the flip would be degenerate,
    /// or the resulting topology would be non-manifold.
    ///
    /// # Example
    ///
    /// ```rust
    /// use delaunay::prelude::triangulation::flips::*;
    ///
    /// 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]),
    ///     vertex!([1.0, 1.0, 1.0]),
    /// ];
    /// let mut dt: DelaunayTriangulation<_, (), (), 3> =
    ///     DelaunayTriangulation::new(&vertices).unwrap();
    ///
    /// // k=3 flips require specific ridge configurations in 3D and above
    /// // This is an illustrative example; actual ridge selection depends on topology
    /// let _ = dt;  // Use dt to prevent unused variable warning
    /// ```
    fn flip_k3(&mut self, ridge: RidgeHandle) -> Result<FlipInfo<D>, FlipError>;

    /// Apply an inverse k=2 flip from an edge star (D >= 3).
    ///
    /// # Errors
    ///
    /// Returns [`FlipError`] if the edge star is invalid or the inverse flip
    /// would create invalid topology.
    fn flip_k2_inverse_from_edge(&mut self, edge: EdgeKey) -> Result<FlipInfo<D>, FlipError>;

    /// Apply an inverse k=3 flip from a triangle star (D >= 4).
    ///
    /// If `D < 4`, this returns [`FlipError::UnsupportedDimension`].
    ///
    /// # Errors
    ///
    /// Returns [`FlipError`] if the triangle star is invalid or the inverse flip
    /// would create invalid topology.
    fn flip_k3_inverse_from_triangle(
        &mut self,
        triangle: TriangleHandle,
    ) -> Result<FlipInfo<D>, FlipError>;
}

impl<K, U, V, const D: usize> BistellarFlips<K, U, D> for Triangulation<K, U, V, D>
where
    K: Kernel<D>,
    U: DataType,
    V: DataType,
{
    fn flip_k1_insert(
        &mut self,
        cell_key: CellKey,
        vertex: Vertex<K::Scalar, U, D>,
    ) -> Result<FlipInfo<D>, FlipError> {
        apply_bistellar_flip_k1(&mut self.tds, cell_key, vertex)
    }

    fn flip_k1_remove(&mut self, vertex_key: VertexKey) -> Result<FlipInfo<D>, FlipError> {
        apply_bistellar_flip_k1_inverse(&mut self.tds, vertex_key)
    }

    fn flip_k2(&mut self, facet: FacetHandle) -> Result<FlipInfo<D>, FlipError> {
        let context = build_k2_flip_context(&self.tds, facet)?;
        apply_bistellar_flip_k2(&mut self.tds, &context)
    }

    fn flip_k3(&mut self, ridge: RidgeHandle) -> Result<FlipInfo<D>, FlipError> {
        let context = build_k3_flip_context(&self.tds, ridge)?;
        apply_bistellar_flip_k3(&mut self.tds, &context)
    }

    fn flip_k2_inverse_from_edge(&mut self, edge: EdgeKey) -> Result<FlipInfo<D>, FlipError> {
        let context = build_k2_flip_context_from_edge(&self.tds, edge)?;
        apply_bistellar_flip_dynamic(&mut self.tds, D, &context)
    }

    fn flip_k3_inverse_from_triangle(
        &mut self,
        triangle: TriangleHandle,
    ) -> Result<FlipInfo<D>, FlipError> {
        if D < 4 {
            return Err(FlipError::UnsupportedDimension { dimension: D });
        }

        let context = build_k3_flip_context_from_triangle(&self.tds, triangle)?;

        // Avoid const-eval underflow for invalid instantiations (e.g. D=0), even though
        // the public contract for this method requires D>=4.
        let k_move = D
            .checked_sub(1)
            .ok_or(FlipError::UnsupportedDimension { dimension: D })?;

        apply_bistellar_flip_dynamic(&mut self.tds, k_move, &context)
    }
}

impl<K, U, V, const D: usize> BistellarFlips<K, U, D> for DelaunayTriangulation<K, U, V, D>
where
    K: Kernel<D>,
    U: DataType,
    V: DataType,
{
    fn flip_k1_insert(
        &mut self,
        cell_key: CellKey,
        vertex: Vertex<K::Scalar, U, D>,
    ) -> Result<FlipInfo<D>, FlipError> {
        let result = self.tri.flip_k1_insert(cell_key, vertex);
        if result.is_ok() {
            self.invalidate_repair_caches();
        }
        result
    }

    fn flip_k1_remove(&mut self, vertex_key: VertexKey) -> Result<FlipInfo<D>, FlipError> {
        let result = self.tri.flip_k1_remove(vertex_key);
        if result.is_ok() {
            self.invalidate_repair_caches();
        }
        result
    }

    fn flip_k2(&mut self, facet: FacetHandle) -> Result<FlipInfo<D>, FlipError> {
        let result = self.tri.flip_k2(facet);
        if result.is_ok() {
            self.invalidate_locate_hint_cache();
        }
        result
    }

    fn flip_k3(&mut self, ridge: RidgeHandle) -> Result<FlipInfo<D>, FlipError> {
        let result = self.tri.flip_k3(ridge);
        if result.is_ok() {
            self.invalidate_locate_hint_cache();
        }
        result
    }

    fn flip_k2_inverse_from_edge(&mut self, edge: EdgeKey) -> Result<FlipInfo<D>, FlipError> {
        let result = self.tri.flip_k2_inverse_from_edge(edge);
        if result.is_ok() {
            self.invalidate_locate_hint_cache();
        }
        result
    }

    fn flip_k3_inverse_from_triangle(
        &mut self,
        triangle: TriangleHandle,
    ) -> Result<FlipInfo<D>, FlipError> {
        let result = self.tri.flip_k3_inverse_from_triangle(triangle);
        if result.is_ok() {
            self.invalidate_locate_hint_cache();
        }
        result
    }
}

#[cfg(test)]
mod tests {
    use super::*;

    use crate::geometry::kernel::{AdaptiveKernel, FastKernel};
    use crate::vertex;
    use slotmap::KeyData;

    #[test]
    fn triangulation_flip_k1_insert_and_remove_roundtrip() {
        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_with_topology_guarantee(
                &vertices,
                crate::core::triangulation::TopologyGuarantee::PLManifold,
            )
            .unwrap();
        let mut tri = dt.as_triangulation().clone();
        let cell_key = tri.cells().next().unwrap().0;

        let inserted = tri
            .flip_k1_insert(cell_key, vertex!([0.25, 0.25, 0.25]))
            .unwrap();
        let inserted_vertex = inserted.inserted_face_vertices[0];
        assert!(!inserted.new_cells.is_empty());
        assert!(tri.validate().is_ok());

        let removed = tri.flip_k1_remove(inserted_vertex).unwrap();
        assert!(!removed.removed_cells.is_empty());
        assert!(tri.validate().is_ok());
    }

    #[test]
    fn triangulation_flip_k2_rejects_invalid_facet_index() {
        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_with_topology_guarantee(
                &vertices,
                crate::core::triangulation::TopologyGuarantee::PLManifold,
            )
            .unwrap();
        let mut tri = dt.as_triangulation().clone();
        let cell_key = tri.cells().next().unwrap().0;

        let err = tri
            .flip_k2(FacetHandle::new(cell_key, u8::MAX))
            .unwrap_err();

        assert!(matches!(
            err,
            FlipError::InvalidFacetIndex {
                cell_key: key,
                facet_index: u8::MAX,
                vertex_count: 4,
            } if key == cell_key
        ));
    }

    #[test]
    fn delaunay_flip_k3_inverse_rejects_unsupported_dimension() {
        let mut dt: DelaunayTriangulation<AdaptiveKernel<f64>, (), (), 3> =
            DelaunayTriangulation::with_empty_kernel(AdaptiveKernel::new());
        let a = VertexKey::from(KeyData::from_ffi(1));
        let b = VertexKey::from(KeyData::from_ffi(2));
        let c = VertexKey::from(KeyData::from_ffi(3));

        let err = dt
            .flip_k3_inverse_from_triangle(TriangleHandle::new(a, b, c))
            .unwrap_err();

        assert!(matches!(
            err,
            FlipError::UnsupportedDimension { dimension: 3 }
        ));
    }

    #[test]
    fn triangulation_flip_k3_inverse_rejects_zero_dimension_without_underflow() {
        let mut tri: Triangulation<FastKernel<f64>, (), (), 0> =
            Triangulation::new_empty(FastKernel::new());
        let a = VertexKey::from(KeyData::from_ffi(1));
        let b = VertexKey::from(KeyData::from_ffi(2));
        let c = VertexKey::from(KeyData::from_ffi(3));

        let err = tri
            .flip_k3_inverse_from_triangle(TriangleHandle::new(a, b, c))
            .unwrap_err();

        assert!(matches!(
            err,
            FlipError::UnsupportedDimension { dimension: 0 }
        ));
    }
}