delaunay 0.7.5

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
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
//! End-to-end "repair then delaunayize" workflow.
//!
//! This module provides `delaunayize_by_flips`, a single public entrypoint that
//! takes an existing [`DelaunayTriangulation`], performs bounded deterministic
//! topology repair toward
//! [`TopologyGuarantee::PLManifold`](crate::core::triangulation::TopologyGuarantee::PLManifold),
//! and then applies
//! standard flip-based Delaunay repair.
//!
//! # Workflow
//!
//! 1. **PL-manifold topology repair** — removes cells that cause facet
//!    over-sharing (codimension-1 facet degree > 2) using a bounded,
//!    deterministic pruning algorithm.
//! 2. **Delaunay flip repair** — runs k=2/k=3 bistellar flips to restore the
//!    empty-circumsphere property.
//! 3. **Optional fallback rebuild** — if configured and both repair passes
//!    fail, rebuilds the triangulation from its vertex set.
//!
//! # Example
//!
//! ```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);
//! ```
//!
//! # Explicitly Deferred
//!
//! - Dedicated targeted repair stages for boundary-ridge multiplicity,
//!   ridge-link manifoldness, and vertex-link manifoldness (#304).
//! - Stronger cell-payload preservation in the fallback rebuild path (#305).

#![forbid(unsafe_code)]

// Re-export outcome field types so users can name them without reaching into
// internal modules.
pub use crate::core::algorithms::flips::DelaunayRepairStats;
pub use crate::core::algorithms::pl_manifold_repair::{
    PlManifoldRepairError, PlManifoldRepairStats,
};

use crate::core::algorithms::flips::DelaunayRepairError;
use crate::core::algorithms::pl_manifold_repair::{
    PlManifoldRepairConfig, repair_facet_oversharing,
};
use crate::core::delaunay_triangulation::{DelaunayRepairHeuristicConfig, DelaunayTriangulation};
use crate::core::traits::data_type::DataType;
use crate::core::vertex::Vertex;
use crate::geometry::kernel::{ExactPredicates, Kernel};
use crate::geometry::traits::coordinate::{CoordinateScalar, ScalarAccumulative};
use thiserror::Error;

// =============================================================================
// CONFIGURATION
// =============================================================================

/// Configuration for the [`delaunayize_by_flips`] workflow.
///
/// # Defaults
///
/// - `topology_max_iterations`: 64
/// - `topology_max_cells_removed`: 10,000
/// - `fallback_rebuild`: false
///
/// # Examples
///
/// ```rust
/// use delaunay::triangulation::delaunayize::DelaunayizeConfig;
///
/// let config = DelaunayizeConfig::default();
/// assert_eq!(config.topology_max_iterations, 64);
/// assert_eq!(config.topology_max_cells_removed, 10_000);
/// assert!(!config.fallback_rebuild);
/// assert!(config.delaunay_max_flips.is_none());
/// ```
#[derive(Debug, Clone, Copy, PartialEq, Eq)]
pub struct DelaunayizeConfig {
    /// Maximum number of topology-repair iterations.
    pub topology_max_iterations: usize,
    /// Maximum number of cells that may be removed during topology repair.
    pub topology_max_cells_removed: usize,
    /// If `true`, rebuild the triangulation from the vertex set when both
    /// topology repair and flip-based Delaunay repair fail.
    ///
    /// **Warning:** the fallback rebuild discards cell-level user data (`V`).
    pub fallback_rebuild: bool,
    /// Optional per-attempt flip budget cap for Delaunay repair.
    ///
    /// `None` (default) uses the internal dimension-dependent budget.
    /// Set to `Some(n)` to limit each repair attempt to at most `n` flips.
    pub delaunay_max_flips: Option<usize>,
}

impl Default for DelaunayizeConfig {
    fn default() -> Self {
        Self {
            topology_max_iterations: 64,
            topology_max_cells_removed: 10_000,
            fallback_rebuild: false,
            delaunay_max_flips: None,
        }
    }
}

// =============================================================================
// OUTCOME
// =============================================================================

/// Outcome of a successful [`delaunayize_by_flips`] call.
///
/// # Examples
///
/// ```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);
/// assert!(!outcome.used_fallback_rebuild);
/// ```
#[derive(Debug, Clone)]
#[non_exhaustive]
pub struct DelaunayizeOutcome<T, U, V, const D: usize>
where
    T: CoordinateScalar,
    U: DataType,
    V: DataType,
{
    /// Statistics from the PL-manifold topology repair pass.
    pub topology_repair: PlManifoldRepairStats<T, U, V, D>,
    /// Statistics from the flip-based Delaunay repair pass.
    pub delaunay_repair: DelaunayRepairStats,
    /// Whether the fallback vertex-set rebuild was used.
    pub used_fallback_rebuild: bool,
}

// =============================================================================
// ERRORS
// =============================================================================

/// Errors that can occur during the [`delaunayize_by_flips`] workflow.
///
/// There are two orthogonal failure modes:
/// - **Topology repair** failed (step 1).
/// - **Delaunay repair** failed (step 2), with optional context about a
///   fallback rebuild attempt.
///
/// # Examples
///
/// ```rust
/// use delaunay::triangulation::delaunayize::DelaunayizeError;
/// use delaunay::core::algorithms::flips::DelaunayRepairError;
/// use delaunay::core::triangulation::TopologyGuarantee;
///
/// let err = DelaunayizeError::DelaunayRepairFailed {
///     source: DelaunayRepairError::InvalidTopology {
///         required: TopologyGuarantee::PLManifold,
///         found: TopologyGuarantee::Pseudomanifold,
///         message: "requires manifold",
///     },
///     context: "fallback not enabled".to_string(),
/// };
/// assert!(err.to_string().contains("Delaunay repair failed"));
/// ```
#[derive(Clone, Debug, Error)]
#[non_exhaustive]
pub enum DelaunayizeError {
    /// PL-manifold topology repair failed.
    #[error("Topology repair failed: {0}")]
    TopologyRepairFailed(#[from] PlManifoldRepairError),

    /// Delaunay flip repair failed (and fallback was not enabled or also failed).
    #[error("Delaunay repair failed ({context}): {source}")]
    DelaunayRepairFailed {
        /// The underlying flip-repair error.
        #[source]
        source: DelaunayRepairError,
        /// Operational context (e.g. "fallback not enabled" or
        /// "fallback rebuild also failed: ...").
        context: String,
    },
}

impl From<DelaunayRepairError> for DelaunayizeError {
    fn from(err: DelaunayRepairError) -> Self {
        Self::DelaunayRepairFailed {
            source: err,
            context: "fallback not enabled".to_string(),
        }
    }
}

// =============================================================================
// PUBLIC API
// =============================================================================

/// Performs bounded topology repair followed by flip-based Delaunay repair.
///
/// This is the primary public entrypoint for the "repair then delaunayize"
/// workflow described in the [module documentation](self).
///
/// # Type Constraints
///
/// The kernel must implement [`ExactPredicates`] (required by the underlying
/// Delaunay flip-repair engine). The default [`AdaptiveKernel`](crate::geometry::kernel::AdaptiveKernel)
/// satisfies this requirement.
///
/// # Errors
///
/// Returns [`DelaunayizeError`] if:
/// - Topology repair fails ([`TopologyRepairFailed`](DelaunayizeError::TopologyRepairFailed)).
/// - Delaunay flip repair fails and fallback is disabled
///   ([`DelaunayRepairFailed`](DelaunayizeError::DelaunayRepairFailed)).
/// - Fallback rebuild also fails (reported via
///   [`DelaunayRepairFailed`](DelaunayizeError::DelaunayRepairFailed) with context).
///
/// # Examples
///
/// ```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);
/// ```
pub fn delaunayize_by_flips<K, U, V, const D: usize>(
    dt: &mut DelaunayTriangulation<K, U, V, D>,
    config: DelaunayizeConfig,
) -> Result<DelaunayizeOutcome<K::Scalar, U, V, D>, DelaunayizeError>
where
    K: Kernel<D> + ExactPredicates,
    K::Scalar: ScalarAccumulative,
    U: DataType,
    V: DataType,
{
    // Snapshot vertex set before repair so fallback rebuild can use it.
    let fallback_vertices: Option<Vec<Vertex<K::Scalar, U, D>>> = if config.fallback_rebuild {
        Some(
            dt.as_triangulation()
                .tds
                .vertices()
                .map(|(_, v)| Vertex::new_with_uuid(*v.point(), v.uuid(), v.data))
                .collect(),
        )
    } else {
        None
    };

    // Step 1: PL-manifold topology repair (facet over-sharing).
    let pl_config = PlManifoldRepairConfig {
        max_iterations: config.topology_max_iterations,
        max_cells_removed: config.topology_max_cells_removed,
    };
    let topology_stats = match repair_facet_oversharing(
        &mut dt.as_triangulation_mut().tds,
        &pl_config,
    ) {
        Ok(stats) => stats,
        Err(topo_err) => {
            if let Some(ref verts) = fallback_vertices {
                // Topology repair failed but fallback is enabled — try rebuilding.
                match DelaunayTriangulation::with_kernel(&dt.as_triangulation().kernel, verts) {
                    Ok(rebuilt) => {
                        *dt = rebuilt;
                        return Ok(DelaunayizeOutcome {
                            topology_repair: PlManifoldRepairStats::default(),
                            delaunay_repair: DelaunayRepairStats::default(),
                            used_fallback_rebuild: true,
                        });
                    }
                    Err(rebuild_err) => {
                        return Err(DelaunayizeError::TopologyRepairFailed(
                            PlManifoldRepairError::Tds(
                                crate::core::triangulation_data_structure::TdsError::InconsistentDataStructure {
                                    message: format!(
                                        "topology repair failed ({topo_err}); fallback rebuild also failed: {rebuild_err}"
                                    ),
                                },
                            ),
                        ));
                    }
                }
            }
            return Err(topo_err.into());
        }
    };

    // Step 2: Flip-based Delaunay repair.
    let delaunay_result = if let Some(max_flips) = config.delaunay_max_flips {
        dt.repair_delaunay_with_flips_advanced(DelaunayRepairHeuristicConfig {
            max_flips: Some(max_flips),
            ..DelaunayRepairHeuristicConfig::default()
        })
        .map(|outcome| outcome.stats)
    } else {
        dt.repair_delaunay_with_flips()
    };

    match delaunay_result {
        Ok(delaunay_stats) => Ok(DelaunayizeOutcome {
            topology_repair: topology_stats,
            delaunay_repair: delaunay_stats,
            used_fallback_rebuild: false,
        }),
        Err(repair_err) => {
            if config.fallback_rebuild {
                // Step 3 (optional): rebuild from vertex set.
                let vertices: Vec<Vertex<K::Scalar, U, D>> = dt
                    .as_triangulation()
                    .tds
                    .vertices()
                    .map(|(_, v)| Vertex::new_with_uuid(*v.point(), v.uuid(), v.data))
                    .collect();

                match DelaunayTriangulation::with_kernel(&dt.as_triangulation().kernel, &vertices) {
                    Ok(rebuilt) => {
                        *dt = rebuilt;
                        // The rebuild succeeded — return stats reflecting the fallback.
                        Ok(DelaunayizeOutcome {
                            topology_repair: topology_stats,
                            delaunay_repair: DelaunayRepairStats::default(),
                            used_fallback_rebuild: true,
                        })
                    }
                    Err(rebuild_err) => Err(DelaunayizeError::DelaunayRepairFailed {
                        source: repair_err,
                        context: format!("fallback rebuild also failed: {rebuild_err}"),
                    }),
                }
            } else {
                Err(DelaunayizeError::from(repair_err))
            }
        }
    }
}

// =============================================================================
// TESTS
// =============================================================================

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

    // =============================================================================
    // HELPER FUNCTIONS
    // =============================================================================

    fn init_tracing() {
        let _ = tracing_subscriber::fmt::try_init();
    }

    // =============================================================================
    // CONFIG DEFAULT TESTS
    // =============================================================================

    #[test]
    fn test_config_defaults() {
        init_tracing();
        let config = DelaunayizeConfig::default();
        assert_eq!(config.topology_max_iterations, 64);
        assert_eq!(config.topology_max_cells_removed, 10_000);
        assert!(!config.fallback_rebuild);
    }

    // =============================================================================
    // SUCCESS PATH TESTS
    // =============================================================================

    #[test]
    fn test_already_delaunay_3d() {
        init_tracing();
        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);
        assert!(!outcome.used_fallback_rebuild);
        assert!(dt.validate().is_ok());
    }

    #[test]
    fn test_already_delaunay_2d() {
        init_tracing();
        let vertices = vec![
            vertex!([0.0, 0.0]),
            vertex!([1.0, 0.0]),
            vertex!([0.0, 1.0]),
            vertex!([1.0, 1.0]),
        ];
        let mut dt: DelaunayTriangulation<_, (), (), 2> =
            DelaunayTriangulation::new(&vertices).unwrap();

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

    // =============================================================================
    // OUTCOME POPULATION TESTS
    // =============================================================================

    #[test]
    fn test_outcome_populated_on_success() {
        init_tracing();
        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.5]),
        ];
        let mut dt: DelaunayTriangulation<_, (), (), 3> =
            DelaunayTriangulation::new(&vertices).unwrap();

        let outcome = delaunayize_by_flips(&mut dt, DelaunayizeConfig::default()).unwrap();
        assert!(outcome.topology_repair.succeeded);
        assert_eq!(outcome.topology_repair.cells_removed, 0);
        assert!(!outcome.used_fallback_rebuild);
    }

    // =============================================================================
    // FALLBACK BEHAVIOR TESTS
    // =============================================================================

    #[test]
    fn test_fallback_disabled_by_default() {
        init_tracing();
        let config = DelaunayizeConfig::default();
        assert!(!config.fallback_rebuild);
    }

    #[test]
    fn test_fallback_enabled_on_valid_triangulation() {
        init_tracing();
        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();

        // Fallback should not be triggered on a valid triangulation.
        let config = DelaunayizeConfig {
            fallback_rebuild: true,
            ..DelaunayizeConfig::default()
        };
        let outcome = delaunayize_by_flips(&mut dt, config).unwrap();
        assert!(!outcome.used_fallback_rebuild);
    }

    // =============================================================================
    // DETERMINISM TESTS
    // =============================================================================

    #[test]
    fn test_deterministic_repeated_runs() {
        init_tracing();
        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.5]),
        ];

        let config = DelaunayizeConfig::default();

        let mut dt1: DelaunayTriangulation<_, (), (), 3> =
            DelaunayTriangulation::new(&vertices).unwrap();
        let outcome1 = delaunayize_by_flips(&mut dt1, config).unwrap();

        let mut dt2: DelaunayTriangulation<_, (), (), 3> =
            DelaunayTriangulation::new(&vertices).unwrap();
        let outcome2 = delaunayize_by_flips(&mut dt2, config).unwrap();

        assert_eq!(
            outcome1.topology_repair.cells_removed,
            outcome2.topology_repair.cells_removed
        );
        assert_eq!(
            outcome1.used_fallback_rebuild,
            outcome2.used_fallback_rebuild
        );
    }
}