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
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
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
//! 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).

#![forbid(unsafe_code)]

// Re-export outcome/error field types so users can name the public contract
// without reaching into lower-level modules.
pub use crate::core::algorithms::flips::{DelaunayRepairError, DelaunayRepairStats};
pub use crate::core::algorithms::pl_manifold_repair::{
    PlManifoldRepairError, PlManifoldRepairStats,
};
pub use crate::core::cell::CellValidationError;
pub use crate::triangulation::delaunay::DelaunayTriangulationConstructionError;

use crate::core::algorithms::pl_manifold_repair::{
    PlManifoldRepairConfig, repair_facet_oversharing,
};
use crate::core::cell::Cell;
use crate::core::collections::{Entry, FastHashMap, Uuid};
use crate::core::tds::{CellKey, Tds};
use crate::core::traits::data_type::DataType;
use crate::core::vertex::Vertex;
use crate::geometry::kernel::{ExactPredicates, Kernel};
use crate::geometry::traits::coordinate::CoordinateScalar;
use crate::triangulation::delaunay::{DelaunayRepairHeuristicConfig, DelaunayTriangulation};
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.
    ///
    /// Cell-level user data (`V`) is restored for rebuilt cells whose sorted
    /// vertex UUID set matches exactly one original cell. Cells that change
    /// during rebuild, have no original payload, or have ambiguous duplicate
    /// original signatures receive `None`.
    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> {
    /// 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 four orthogonal failure modes:
/// - **Topology repair** failed (step 1).
/// - **Delaunay repair** failed (step 2), with optional context about a
///   fallback rebuild attempt.
/// - **Fallback snapshot** failed before any repair could run.
/// - **Fallback cell-data recovery** failed while snapshotting or restoring
///   cell payloads after a repair failure.
///
/// # Orthogonality
///
/// The variants are mutually exclusive by failure mode:
/// - Topology repair, fallback not attempted -> [`TopologyRepairFailed`](Self::TopologyRepairFailed).
/// - Topology repair, fallback also failed   -> [`TopologyRepairFailedWithRebuild`](Self::TopologyRepairFailedWithRebuild).
/// - Topology repair, fallback rebuild succeeded but payload restore failed -> [`TopologyRepairFailedWithRebuildRestore`](Self::TopologyRepairFailedWithRebuildRestore).
/// - Delaunay repair, fallback not attempted -> [`DelaunayRepairFailed`](Self::DelaunayRepairFailed).
/// - Delaunay repair, fallback payload snapshot failed -> [`DelaunayRepairFailedWithCellDataSnapshot`](Self::DelaunayRepairFailedWithCellDataSnapshot).
/// - Delaunay repair, fallback also failed   -> [`DelaunayRepairFailedWithRebuild`](Self::DelaunayRepairFailedWithRebuild).
/// - Delaunay repair, fallback rebuild succeeded but payload restore failed -> [`DelaunayRepairFailedWithRebuildRestore`](Self::DelaunayRepairFailedWithRebuildRestore).
/// - Fallback was enabled, but the pre-repair payload snapshot failed -> [`FallbackCellDataSnapshotFailed`](Self::FallbackCellDataSnapshotFailed).
///
/// Variants with secondary fallback failures preserve **both** the primary
/// repair error and the secondary construction, snapshot, or restore error as
/// typed values (no stringification), so consumers can inspect both errors via
/// pattern matching; the primary repair error is exposed via
/// [`Error::source`](std::error::Error::source).
///
/// # Examples
///
/// ```rust
/// use delaunay::prelude::triangulation::delaunayize::*;
///
/// let err = DelaunayizeError::DelaunayRepairFailed {
///     source: DelaunayRepairError::PostconditionFailed {
///         message: "still non-Delaunay after repair".to_string(),
///     },
/// };
/// assert!(err.to_string().contains("Delaunay repair failed"));
/// ```
#[derive(Clone, Debug, Error, PartialEq, Eq)]
#[non_exhaustive]
pub enum DelaunayizeError {
    /// PL-manifold topology repair failed; no fallback rebuild was attempted
    /// (fallback disabled, or the caller's config did not request one).
    #[error("Topology repair failed: {source}")]
    TopologyRepairFailed {
        /// The underlying topology-repair error.
        #[from]
        #[source]
        source: PlManifoldRepairError,
    },

    /// PL-manifold topology repair failed **and** the fallback vertex-set
    /// rebuild also failed.  Both errors are preserved as typed values.
    #[error("Topology repair failed ({source}); fallback rebuild also failed: {rebuild_error}")]
    TopologyRepairFailedWithRebuild {
        /// The underlying topology-repair error that triggered the fallback.
        #[source]
        source: PlManifoldRepairError,
        /// The construction error from the subsequent vertex-set rebuild attempt.
        rebuild_error: DelaunayTriangulationConstructionError,
    },

    /// PL-manifold topology repair failed, the fallback vertex-set rebuild
    /// succeeded, but cell-payload restoration from the rebuilt topology failed.
    #[error(
        "Topology repair failed ({source}); fallback rebuild succeeded but cell-data restore failed: {restore_error}"
    )]
    TopologyRepairFailedWithRebuildRestore {
        /// The underlying topology-repair error that triggered the fallback.
        #[source]
        source: PlManifoldRepairError,
        /// The cell-data restoration error from the rebuilt triangulation.
        restore_error: CellValidationError,
    },

    /// Delaunay flip repair failed; no fallback rebuild was attempted
    /// (fallback disabled, or the caller's config did not request one).
    #[error("Delaunay repair failed: {source}")]
    DelaunayRepairFailed {
        /// The underlying flip-repair error.
        #[from]
        #[source]
        source: DelaunayRepairError,
    },

    /// Delaunay flip repair failed and the fallback payload snapshot could not
    /// be collected from the current triangulation.
    #[error(
        "Delaunay repair failed ({source}); fallback cell-data snapshot failed: {snapshot_error}"
    )]
    DelaunayRepairFailedWithCellDataSnapshot {
        /// The underlying flip-repair error that triggered the fallback.
        #[source]
        source: DelaunayRepairError,
        /// The cell-data snapshot error from the current triangulation.
        snapshot_error: CellValidationError,
    },

    /// Delaunay flip repair failed **and** the fallback vertex-set rebuild
    /// also failed.  Both errors are preserved as typed values.
    #[error("Delaunay repair failed ({source}); fallback rebuild also failed: {rebuild_error}")]
    DelaunayRepairFailedWithRebuild {
        /// The underlying flip-repair error that triggered the fallback.
        #[source]
        source: DelaunayRepairError,
        /// The construction error from the subsequent vertex-set rebuild attempt.
        rebuild_error: DelaunayTriangulationConstructionError,
    },

    /// Delaunay flip repair failed, the fallback vertex-set rebuild succeeded,
    /// but cell-payload restoration from the rebuilt topology failed.
    #[error(
        "Delaunay repair failed ({source}); fallback rebuild succeeded but cell-data restore failed: {restore_error}"
    )]
    DelaunayRepairFailedWithRebuildRestore {
        /// The underlying flip-repair error that triggered the fallback.
        #[source]
        source: DelaunayRepairError,
        /// The cell-data restoration error from the rebuilt triangulation.
        restore_error: CellValidationError,
    },

    /// Fallback rebuild was enabled, but the pre-repair cell-payload snapshot
    /// could not be collected from the input triangulation. No topology or
    /// Delaunay repair was attempted.
    #[error("Fallback cell-data snapshot failed before repair; no repair attempted: {source}")]
    FallbackCellDataSnapshotFailed {
        /// The cell-data snapshot error from the input triangulation.
        #[from]
        #[source]
        source: CellValidationError,
    },
}

// =============================================================================
// HELPERS
// =============================================================================

#[derive(Clone, Debug, PartialEq, Eq)]
enum CellDataMatch<V> {
    Unique(Option<V>),
    Ambiguous,
}

type CellDataByVertexUuids<V> = FastHashMap<Vec<Uuid>, CellDataMatch<V>>;
type FallbackRebuildState<T, U, V, const D: usize> =
    (Vec<Vertex<T, U, D>>, CellDataByVertexUuids<V>);

/// Captures the fallback rebuild inputs from the current TDS, including typed
/// failure if any cell cannot resolve its vertex UUID identity.
fn snapshot_rebuild_state<T, U, V, const D: usize>(
    tds: &Tds<T, U, V, D>,
) -> Result<FallbackRebuildState<T, U, V, D>, CellValidationError>
where
    T: CoordinateScalar,
    U: DataType,
    V: DataType,
{
    let vertices = tds
        .vertices()
        .map(|(_, v)| Vertex::new_with_uuid(*v.point(), v.uuid(), v.data))
        .collect::<Vec<_>>();
    let cell_data = collect_cell_data(tds)?;
    Ok((vertices, cell_data))
}

/// Hashes cell payloads by sorted vertex UUIDs so fallback rebuilds can
/// recover payloads for cells whose vertex set survives unchanged.
fn collect_cell_data<T, U, V, const D: usize>(
    tds: &Tds<T, U, V, D>,
) -> Result<CellDataByVertexUuids<V>, CellValidationError>
where
    U: DataType,
    V: DataType,
{
    let mut cell_data = FastHashMap::default();
    for (_, cell) in tds.cells() {
        let vertex_uuids = cell_vertex_uuids(tds, cell)?;
        match cell_data.entry(vertex_uuids) {
            Entry::Vacant(entry) => {
                entry.insert(CellDataMatch::Unique(cell.data().copied()));
            }
            Entry::Occupied(mut entry) => {
                entry.insert(CellDataMatch::Ambiguous);
            }
        }
    }
    Ok(cell_data)
}

/// Builds the order-independent cell identity used to match original and
/// rebuilt cells across fallback reconstruction.
fn cell_vertex_uuids<T, U, V, const D: usize>(
    tds: &Tds<T, U, V, D>,
    cell: &Cell<T, U, V, D>,
) -> Result<Vec<Uuid>, CellValidationError>
where
    U: DataType,
    V: DataType,
{
    let mut vertex_uuids = cell
        .vertex_uuid_iter(tds)
        .collect::<Result<Vec<_>, CellValidationError>>()?;
    vertex_uuids.sort_unstable();
    Ok(vertex_uuids)
}

/// Reattaches original cell payloads to rebuilt cells that retain the same
/// vertex UUID set after fallback reconstruction.
fn restore_cell_data<K, U, V, const D: usize>(
    rebuilt: &mut DelaunayTriangulation<K, U, V, D>,
    original_cell_data: &CellDataByVertexUuids<V>,
) -> Result<(), CellValidationError>
where
    K: Kernel<D>,
    U: DataType,
    V: DataType,
{
    let mut assignments: Vec<(CellKey, V)> = Vec::new();
    for (cell_key, cell) in rebuilt.cells() {
        let vertex_uuids = cell_vertex_uuids(rebuilt.tds(), cell)?;
        let Some(CellDataMatch::Unique(Some(data))) = original_cell_data.get(&vertex_uuids) else {
            continue;
        };
        assignments.push((cell_key, *data));
    }

    for (cell_key, data) in assignments {
        rebuilt.set_cell_data(cell_key, Some(data));
    }

    Ok(())
}

#[derive(Clone, Debug, Error, PartialEq, Eq)]
enum FallbackRebuildError {
    #[error("fallback rebuild failed: {source}")]
    Construction {
        #[from]
        #[source]
        source: DelaunayTriangulationConstructionError,
    },
    #[error("fallback cell-data restore failed: {source}")]
    Restore {
        #[from]
        #[source]
        source: CellValidationError,
    },
}

/// Maps a fallback rebuild failure while handling a topology-repair failure
/// without erasing either typed source.
fn topology_rebuild_error(
    source: PlManifoldRepairError,
    fallback_error: FallbackRebuildError,
) -> DelaunayizeError {
    match fallback_error {
        FallbackRebuildError::Construction {
            source: rebuild_error,
        } => DelaunayizeError::TopologyRepairFailedWithRebuild {
            source,
            rebuild_error,
        },
        FallbackRebuildError::Restore {
            source: restore_error,
        } => DelaunayizeError::TopologyRepairFailedWithRebuildRestore {
            source,
            restore_error,
        },
    }
}

/// Maps a fallback rebuild failure while handling a Delaunay-repair failure
/// without erasing either typed source.
fn delaunay_rebuild_error(
    source: DelaunayRepairError,
    fallback_error: FallbackRebuildError,
) -> DelaunayizeError {
    match fallback_error {
        FallbackRebuildError::Construction {
            source: rebuild_error,
        } => DelaunayizeError::DelaunayRepairFailedWithRebuild {
            source,
            rebuild_error,
        },
        FallbackRebuildError::Restore {
            source: restore_error,
        } => DelaunayizeError::DelaunayRepairFailedWithRebuildRestore {
            source,
            restore_error,
        },
    }
}

/// Rebuilds a triangulation from preserved vertices while restoring any
/// cell payloads whose vertex UUID signatures survive the rebuild unchanged.
fn rebuild_preserving_data<K, U, V, const D: usize>(
    kernel: &K,
    vertices: &[Vertex<K::Scalar, U, D>],
    original_cell_data: &CellDataByVertexUuids<V>,
) -> Result<DelaunayTriangulation<K, U, V, D>, FallbackRebuildError>
where
    K: Kernel<D> + ExactPredicates,
    U: DataType,
    V: DataType,
{
    let mut rebuilt = DelaunayTriangulation::with_kernel(kernel, vertices)?;
    restore_cell_data(&mut rebuilt, original_cell_data)?;
    Ok(rebuilt)
}

// =============================================================================
// 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 and no fallback rebuild was attempted
///   ([`TopologyRepairFailed`](DelaunayizeError::TopologyRepairFailed)).
/// - Topology repair fails **and** the fallback vertex-set rebuild also
///   fails
///   ([`TopologyRepairFailedWithRebuild`](DelaunayizeError::TopologyRepairFailedWithRebuild)).
/// - Delaunay flip repair fails and no fallback rebuild was attempted
///   ([`DelaunayRepairFailed`](DelaunayizeError::DelaunayRepairFailed)).
/// - Delaunay flip repair fails **and** the fallback vertex-set rebuild also
///   fails
///   ([`DelaunayRepairFailedWithRebuild`](DelaunayizeError::DelaunayRepairFailedWithRebuild)).
///
/// The `*WithRebuild` variants preserve both errors as typed fields so
/// consumers can inspect both typed errors;
/// [`Error::source`](std::error::Error::source) exposes the primary repair error.
///
/// # 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);
/// ```
#[expect(
    clippy::result_large_err,
    reason = "DelaunayizeError preserves typed source and rebuild_error values on the *WithRebuild variants (no boxing) so callers can pattern-match both errors while Error::source exposes the primary repair error; this is a cold error path."
)]
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,
    U: DataType,
    V: DataType,
{
    // Snapshot vertex and cell data before repair so fallback rebuild can use it.
    let fallback_state: Option<FallbackRebuildState<K::Scalar, U, V, D>> =
        if config.fallback_rebuild {
            let tds = &dt.as_triangulation().tds;
            Some(snapshot_rebuild_state(tds)?)
        } 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,
            // Topology repair failed but fallback is enabled — try rebuilding.
            Err(topo_err) if let Some((ref verts, ref cell_data)) = fallback_state => {
                match rebuild_preserving_data(&dt.as_triangulation().kernel, verts, cell_data) {
                    Ok(rebuilt) => {
                        *dt = rebuilt;
                        return Ok(DelaunayizeOutcome {
                            topology_repair: PlManifoldRepairStats {
                                succeeded: true,
                                ..PlManifoldRepairStats::default()
                            },
                            delaunay_repair: DelaunayRepairStats::default(),
                            used_fallback_rebuild: true,
                        });
                    }
                    Err(fallback_error) => {
                        return Err(topology_rebuild_error(topo_err, fallback_error));
                    }
                }
            }
            Err(topo_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 tds = &dt.as_triangulation().tds;
                let (vertices, cell_data) = match snapshot_rebuild_state(tds) {
                    Ok(state) => state,
                    Err(snapshot_error) => {
                        return Err(DelaunayizeError::DelaunayRepairFailedWithCellDataSnapshot {
                            source: repair_err,
                            snapshot_error,
                        });
                    }
                };

                match rebuild_preserving_data(&dt.as_triangulation().kernel, &vertices, &cell_data)
                {
                    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(fallback_error) => Err(delaunay_rebuild_error(repair_err, fallback_error)),
                }
            } else {
                Err(DelaunayizeError::from(repair_err))
            }
        }
    }
}

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

#[cfg(test)]
mod tests {
    use super::*;
    use crate::core::{tds::VertexKey, triangulation::TriangulationConstructionError};
    use crate::geometry::kernel::AdaptiveKernel;
    use crate::triangulation::builder::DelaunayTriangulationBuilder;
    use crate::vertex;
    use slotmap::KeyData;
    use std::error::Error as StdError;

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

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

    fn construction_error() -> DelaunayTriangulationConstructionError {
        DelaunayTriangulationConstructionError::Triangulation(
            TriangulationConstructionError::FailedToCreateCell {
                message: "synthetic rebuild degeneracy".to_string(),
            },
        )
    }

    // =============================================================================
    // 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);
        assert!(config.delaunay_max_flips.is_none());
    }

    // =============================================================================
    // 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);
    }

    // =============================================================================
    // ERROR PATH TESTS
    // =============================================================================

    #[test]
    fn test_collect_cell_data_missing_vertex() {
        let mut tds: Tds<f64, (), i32, 2> = Tds::empty();
        let vertex_keys: Vec<_> = [
            vertex!([0.0, 0.0]),
            vertex!([1.0, 0.0]),
            vertex!([0.0, 1.0]),
        ]
        .iter()
        .map(|vertex| tds.insert_vertex_with_mapping(*vertex).unwrap())
        .collect();
        let missing = vertex_keys[0];
        let cell = Cell::new(vertex_keys, Some(7)).unwrap();
        let dangling_cell = cell.clone();
        let cell_key = tds.insert_cell_with_mapping(cell).unwrap();
        assert_eq!(tds.remove_cells_by_keys(&[cell_key]), 1);
        tds.remove_isolated_vertex(missing);
        tds.cells_mut().insert(dangling_cell);

        let err = collect_cell_data(&tds).unwrap_err();

        assert_eq!(err, CellValidationError::VertexKeyNotFound { key: missing });
    }

    #[test]
    fn test_snapshot_error_source() {
        let source = CellValidationError::VertexKeyNotFound {
            key: VertexKey::from(KeyData::from_ffi(0xBAD)),
        };
        let err = DelaunayizeError::FallbackCellDataSnapshotFailed {
            source: source.clone(),
        };

        assert_eq!(
            err,
            DelaunayizeError::FallbackCellDataSnapshotFailed {
                source: source.clone()
            }
        );
        assert!(
            err.to_string()
                .contains("Fallback cell-data snapshot failed")
        );
        assert!(err.to_string().contains("no repair attempted"));
        let error_source = StdError::source(&err).unwrap();
        assert_eq!(error_source.to_string(), source.to_string());
    }

    #[test]
    fn test_repair_snapshot_error_source() {
        let source = DelaunayRepairError::PostconditionFailed {
            message: "synthetic postcondition".to_string(),
        };
        let snapshot_error = CellValidationError::VertexKeyNotFound {
            key: VertexKey::from(KeyData::from_ffi(0xBAD)),
        };
        let err = DelaunayizeError::DelaunayRepairFailedWithCellDataSnapshot {
            source: source.clone(),
            snapshot_error: snapshot_error.clone(),
        };

        assert_eq!(
            err,
            DelaunayizeError::DelaunayRepairFailedWithCellDataSnapshot {
                source: source.clone(),
                snapshot_error,
            }
        );
        assert!(
            err.to_string()
                .contains("fallback cell-data snapshot failed")
        );
        let error_source = StdError::source(&err).unwrap();
        assert_eq!(error_source.to_string(), source.to_string());
    }

    #[test]
    fn test_restore_error_sources() {
        let topology_source = PlManifoldRepairError::NoProgress {
            over_shared_facets: 2,
            iterations: 3,
            cells_removed: 4,
        };
        let delaunay_source = DelaunayRepairError::PostconditionFailed {
            message: "synthetic postcondition".to_string(),
        };
        let restore_error = CellValidationError::VertexKeyNotFound {
            key: VertexKey::from(KeyData::from_ffi(0xBAD)),
        };

        let topology_err = DelaunayizeError::TopologyRepairFailedWithRebuildRestore {
            source: topology_source.clone(),
            restore_error: restore_error.clone(),
        };
        let delaunay_err = DelaunayizeError::DelaunayRepairFailedWithRebuildRestore {
            source: delaunay_source.clone(),
            restore_error,
        };

        assert!(
            topology_err
                .to_string()
                .contains("cell-data restore failed")
        );
        assert_eq!(
            StdError::source(&topology_err).unwrap().to_string(),
            topology_source.to_string()
        );
        assert!(
            delaunay_err
                .to_string()
                .contains("cell-data restore failed")
        );
        assert_eq!(
            StdError::source(&delaunay_err).unwrap().to_string(),
            delaunay_source.to_string()
        );
    }

    #[test]
    fn test_topology_rebuild_error_mapping() {
        let source = PlManifoldRepairError::NoProgress {
            over_shared_facets: 2,
            iterations: 3,
            cells_removed: 4,
        };
        let rebuild_error = construction_error();
        let restore_error = CellValidationError::VertexKeyNotFound {
            key: VertexKey::from(KeyData::from_ffi(0xBAD)),
        };

        let rebuild_err = topology_rebuild_error(
            source.clone(),
            FallbackRebuildError::Construction {
                source: rebuild_error.clone(),
            },
        );
        assert_eq!(
            rebuild_err,
            DelaunayizeError::TopologyRepairFailedWithRebuild {
                source: source.clone(),
                rebuild_error,
            }
        );
        assert!(
            rebuild_err
                .to_string()
                .contains("fallback rebuild also failed")
        );

        let restore_err = topology_rebuild_error(
            source.clone(),
            FallbackRebuildError::Restore {
                source: restore_error.clone(),
            },
        );
        assert_eq!(
            restore_err,
            DelaunayizeError::TopologyRepairFailedWithRebuildRestore {
                source,
                restore_error,
            }
        );
        assert!(
            restore_err
                .to_string()
                .contains("fallback rebuild succeeded but cell-data restore failed")
        );
    }

    #[test]
    fn test_delaunay_rebuild_error_mapping() {
        let source = DelaunayRepairError::PostconditionFailed {
            message: "synthetic postcondition".to_string(),
        };
        let rebuild_error = construction_error();
        let restore_error = CellValidationError::VertexKeyNotFound {
            key: VertexKey::from(KeyData::from_ffi(0xBAD)),
        };

        let rebuild_err = delaunay_rebuild_error(
            source.clone(),
            FallbackRebuildError::Construction {
                source: rebuild_error.clone(),
            },
        );
        assert_eq!(
            rebuild_err,
            DelaunayizeError::DelaunayRepairFailedWithRebuild {
                source: source.clone(),
                rebuild_error,
            }
        );
        assert!(
            rebuild_err
                .to_string()
                .contains("fallback rebuild also failed")
        );

        let restore_err = delaunay_rebuild_error(
            source.clone(),
            FallbackRebuildError::Restore {
                source: restore_error.clone(),
            },
        );
        assert_eq!(
            restore_err,
            DelaunayizeError::DelaunayRepairFailedWithRebuildRestore {
                source,
                restore_error,
            }
        );
        assert!(
            restore_err
                .to_string()
                .contains("fallback rebuild succeeded but cell-data restore failed")
        );
    }

    // =============================================================================
    // 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);
    }

    #[test]
    fn test_rebuild_restores_cell_data() {
        init_tracing();
        let vertices = [
            vertex!([0.0, 0.0]),
            vertex!([1.0, 0.0]),
            vertex!([0.0, 1.0]),
        ];
        let mut dt = DelaunayTriangulationBuilder::new(&vertices)
            .build::<i32>()
            .unwrap();
        let original_cell_key = dt.cells().next().unwrap().0;
        dt.set_cell_data(original_cell_key, Some(42));

        let tds = &dt.as_triangulation().tds;
        let vertices: Vec<_> = tds
            .vertices()
            .map(|(_, v)| Vertex::new_with_uuid(*v.point(), v.uuid(), v.data))
            .collect();
        let cell_data = collect_cell_data(tds).unwrap();

        let rebuilt =
            rebuild_preserving_data(&dt.as_triangulation().kernel, &vertices, &cell_data).unwrap();

        let (_, rebuilt_cell) = rebuilt.cells().next().unwrap();
        assert_eq!(rebuilt_cell.data(), Some(&42));
        assert!(rebuilt.validate().is_ok());
    }

    #[test]
    fn test_rebuild_drops_ambiguous_data() {
        init_tracing();
        let vertices = [
            vertex!([0.0, 0.0]),
            vertex!([1.0, 0.0]),
            vertex!([0.0, 1.0]),
        ];
        let mut tds: Tds<f64, (), i32, 2> = Tds::empty();
        let vertex_keys: Vec<_> = vertices
            .iter()
            .map(|vertex| tds.insert_vertex_with_mapping(*vertex).unwrap())
            .collect();

        let duplicate_a = Cell::new(vertex_keys.clone(), Some(42)).unwrap();
        let duplicate_b = Cell::new(vertex_keys, Some(42)).unwrap();
        tds.insert_cell_with_mapping(duplicate_a).unwrap();
        tds.insert_cell_with_mapping(duplicate_b).unwrap();

        let rebuild_vertices: Vec<_> = tds
            .vertices()
            .map(|(_, v)| Vertex::new_with_uuid(*v.point(), v.uuid(), v.data))
            .collect();
        let cell_data = collect_cell_data(&tds).unwrap();
        let kernel = AdaptiveKernel::new();
        let mut rebuilt: DelaunayTriangulation<_, (), i32, 2> =
            DelaunayTriangulation::with_kernel(&kernel, &rebuild_vertices).unwrap();

        restore_cell_data(&mut rebuilt, &cell_data).unwrap();

        let (_, rebuilt_cell) = rebuilt.cells().next().unwrap();
        assert_eq!(rebuilt_cell.data(), None);
        assert!(rebuilt.validate().is_ok());
    }

    // =============================================================================
    // 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
        );
    }
}