delaunay 0.7.2

A d-dimensional Delaunay triangulation library with float coordinate support
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
//! Facet key utilities.

#![forbid(unsafe_code)]

use crate::core::facet::FacetError;
use crate::core::traits::data_type::DataType;
use crate::core::triangulation_data_structure::{CellKey, Tds, VertexKey};
use crate::core::util::hashing::stable_hash_u64_slice;
use crate::geometry::traits::coordinate::CoordinateScalar;
use slotmap::Key;
use thiserror::Error;

/// Derives a facet key directly from vertex keys.
///
/// Computes the canonical facet key for lookup in facet-to-cells mappings. This is useful
/// in hot paths like visibility checking in convex hull algorithms and boundary analysis.
///
/// If you have `Vertex` instances instead of `VertexKey`s, obtain the keys via the TDS
/// using the vertex's UUID (e.g., `tds.vertex_key_from_uuid(vertex.uuid())`).
///
/// # Arguments
///
/// * `facet_vertex_keys` - The vertex keys that make up the facet
///
/// # Returns
///
/// A `Result` containing the facet key or a `FacetError` if validation fails.
///
/// # Errors
///
/// Returns `FacetError::InsufficientVertices` if the vertex count doesn't equal `D`
///
/// # Examples
///
/// ```
/// use delaunay::prelude::query::*;
/// use delaunay::core::util::checked_facet_key_from_vertex_keys;
///
/// 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::new(&vertices).unwrap();
/// let tds = dt.tds();
///
/// // Get facet vertex keys from a cell - no need to materialize Vertex objects
/// if let Some(cell) = tds.cells().map(|(_, cell)| cell).next() {
///     let facet_vertex_keys: Vec<_> = cell.vertices().iter().skip(1).copied().collect(); // Skip 1 vertex to get D vertices
///     assert_eq!(facet_vertex_keys.len(), 3); // For 3D triangulation, facet has 3 vertices
///     let facet_key = checked_facet_key_from_vertex_keys::<3>(&facet_vertex_keys).unwrap();
///     println!("Facet key: {}", facet_key);
/// }
/// ```
///
/// # Performance
///
/// - Time Complexity: O(D log D) where D is the facet dimension (for sorting vertex keys)
/// - Space Complexity: O(D) for the temporary sorted buffer (stack-allocated via `SmallBuffer`)
/// - Improves cache locality by working only with compact `VertexKey` types
///
/// # See Also
///
/// - [`crate::core::facet::facet_key_from_vertices`] - Low-level function that computes the hash from keys
pub fn checked_facet_key_from_vertex_keys<const D: usize>(
    facet_vertex_keys: &[VertexKey],
) -> Result<u64, FacetError> {
    use crate::core::facet::facet_key_from_vertices;

    // Validate that the number of vertex keys matches the expected dimension
    // In a D-dimensional triangulation, a facet should have exactly D vertices
    if facet_vertex_keys.len() != D {
        return Err(FacetError::InsufficientVertices {
            expected: D,
            actual: facet_vertex_keys.len(),
            dimension: D,
        });
    }

    // Directly compute the facet key from vertex keys
    // facet_key_from_vertices handles the sorting internally
    Ok(facet_key_from_vertices(facet_vertex_keys))
}

/// Errors that can occur while deriving a periodic facet signature from lifted vertices.
#[derive(Clone, Debug, Error, PartialEq, Eq)]
pub(crate) enum PeriodicFacetKeyDerivationError {
    /// The lifted cell does not have the expected simplex arity (`D + 1` vertices).
    #[error("Invalid lifted cell arity: expected {expected} vertices, got {actual}")]
    InvalidLiftedCellArity {
        /// Expected number of lifted vertices (`D + 1`).
        expected: usize,
        /// Actual lifted vertex count provided by the caller.
        actual: usize,
    },

    /// The requested facet index exceeds the lifted-vertex count.
    #[error("Facet index {facet_index} out of bounds for lifted vertex count {vertex_count}")]
    FacetIndexOutOfBounds {
        /// Requested facet index.
        facet_index: usize,
        /// Number of lifted vertices available.
        vertex_count: usize,
    },

    /// Relative periodic offset component is outside the encodable byte-delta range.
    #[error(
        "Periodic offset component {component} (axis {axis}) is out of encodable range 0..=255"
    )]
    RelativeOffsetOutOfRange {
        /// Axis whose shifted component failed.
        axis: usize,
        /// Shifted component value.
        component: i16,
    },
}

/// Computes a translation-invariant periodic facet key from lifted cell vertices.
///
/// `lifted_vertices` must represent one full lifted cell as `(vertex_key, lattice_offset)` pairs.
/// The facet opposite `facet_index` is selected, vertex keys are sorted for permutation
/// invariance, and offsets are normalized against the first facet vertex so equivalent
/// periodic facets hash identically.
pub(crate) fn periodic_facet_key_from_lifted_vertices<const D: usize>(
    lifted_vertices: &[(VertexKey, [i8; D])],
    facet_index: usize,
) -> Result<u64, PeriodicFacetKeyDerivationError> {
    let expected_arity = D + 1;
    if lifted_vertices.len() != expected_arity {
        return Err(PeriodicFacetKeyDerivationError::InvalidLiftedCellArity {
            expected: expected_arity,
            actual: lifted_vertices.len(),
        });
    }

    if facet_index >= lifted_vertices.len() {
        return Err(PeriodicFacetKeyDerivationError::FacetIndexOutOfBounds {
            facet_index,
            vertex_count: lifted_vertices.len(),
        });
    }

    let mut lifted_facet: Vec<(u64, [i8; D])> =
        Vec::with_capacity(lifted_vertices.len().saturating_sub(1));
    for (idx, (vertex_key, offset)) in lifted_vertices.iter().enumerate() {
        if idx != facet_index {
            lifted_facet.push((vertex_key.data().as_ffi(), *offset));
        }
    }
    // Sort by (vertex_key_value, offset) to make ordering deterministic when keys are equal
    lifted_facet.sort_unstable_by(|(key_a, offset_a), (key_b, offset_b)| {
        key_a.cmp(key_b).then_with(|| offset_a.cmp(offset_b))
    });
    let facet_anchor_offset = lifted_facet
        .first()
        .map_or([0_i8; D], |(_, offset)| *offset);

    let mut packed_signature: Vec<u64> = Vec::with_capacity(lifted_facet.len() * (D + 1));
    for (vertex_key_value, offset) in lifted_facet {
        packed_signature.push(vertex_key_value);
        for axis in 0..D {
            let shifted = i16::from(offset[axis]) - i16::from(facet_anchor_offset[axis]) + 128;
            if !(0..=255).contains(&shifted) {
                return Err(PeriodicFacetKeyDerivationError::RelativeOffsetOutOfRange {
                    axis,
                    component: shifted,
                });
            }
            packed_signature
                .push(u64::try_from(shifted).expect("validated shifted offset is in 0..=255"));
        }
    }

    Ok(stable_hash_u64_slice(&packed_signature))
}

/// Verifies facet index consistency between two neighboring cells.
///
/// This function checks that a shared facet computed from both cells' perspectives
/// produces the same facet key, which is critical for catching subtle neighbor
/// assignment errors in triangulation algorithms.
///
/// # Arguments
///
/// * `tds` - The triangulation data structure
/// * `cell1_key` - Key of the first cell
/// * `cell2_key` - Key of the second (neighboring) cell
/// * `facet_idx` - Index of the facet in cell1 that should match a facet in cell2
///
/// # Returns
///
/// `Ok(true)` if a matching facet is found in cell2 with the same facet key.
/// `Ok(false)` if no matching facet is found.
/// `Err(FacetError)` if there's an error accessing cell or facet data.
///
/// # Errors
///
/// Returns `FacetError` if:
/// - Either cell cannot be found in the TDS
/// - Facet views cannot be created from the cells
/// - Facet vertices cannot be accessed
/// - The facet index is out of bounds
///
/// # Examples
///
/// ```rust,no_run
/// use delaunay::core::util::verify_facet_index_consistency;
/// use delaunay::core::triangulation_data_structure::Tds;
///
/// fn validate_neighbor_consistency(
///     tds: &Tds<f64, (), (), 3>,
/// ) -> Result<(), String> {
///     // Get two neighboring cell keys
///     let cell_keys: Vec<_> = tds.cell_keys().take(2).collect();
///     if cell_keys.len() >= 2 {
///         // Check if facet 0 of cell1 matches a facet in cell2
///         let consistent = verify_facet_index_consistency(
///             tds,
///             cell_keys[0],
///             cell_keys[1],
///             0,
///         )
///         .map_err(|e| format!("Facet error: {}", e))?;
///
///         if consistent {
///             println!("Facet indices are consistent");
///             Ok(())
///         } else {
///             Err("No matching facet found in neighbor".to_string())
///         }
///     } else {
///         Ok(())
///     }
/// }
/// ```
///
/// # Performance
///
/// - Time Complexity: O(D²) where D is the dimension (iterates over facets and vertices)
/// - Space Complexity: O(D) for temporary vertex buffers
pub fn verify_facet_index_consistency<const D: usize>(
    tds: &Tds<impl CoordinateScalar, impl DataType, impl DataType, D>,
    cell1_key: CellKey,
    cell2_key: CellKey,
    facet_idx: usize,
) -> Result<bool, FacetError> {
    // Get facet views from both cells (validates cells exist)
    let cell1_facets = crate::core::cell::Cell::facet_views_from_tds(tds, cell1_key)?;
    let cell2_facets = crate::core::cell::Cell::facet_views_from_tds(tds, cell2_key)?;

    // Check facet index bounds
    if facet_idx >= cell1_facets.len() {
        // Saturate to u8::MAX for error reporting if index overflows u8
        let idx_u8 = u8::try_from(facet_idx).unwrap_or(u8::MAX);
        return Err(FacetError::InvalidFacetIndex {
            index: idx_u8,
            facet_count: cell1_facets.len(),
        });
    }

    // Get the facet from cell1 and compute its key
    let cell1_facet = &cell1_facets[facet_idx];
    let cell1_key_value = cell1_facet.key()?;

    // Find matching facet in cell2
    for cell2_facet in &cell2_facets {
        if cell1_key_value == cell2_facet.key()? {
            return Ok(true);
        }
    }

    Ok(false) // No matching facet found
}

/// Helper function to safely convert usize to u8 for facet indices.
///
/// This function provides a centralized, safe conversion from `usize` to `u8`
/// that is commonly needed throughout the triangulation codebase for facet indexing.
/// It handles the conversion error gracefully by returning appropriate `FacetError`
/// variants with detailed error information.
///
/// # Arguments
///
/// * `idx` - The usize index to convert
/// * `facet_count` - The number of facets (for error reporting)
///
/// # Returns
///
/// A `Result` containing the converted u8 index or a `FacetError`.
///
/// # Errors
///
/// Returns `FacetError::InvalidFacetIndexOverflow` if the index cannot fit in a u8.
///
/// # Examples
///
/// ```
/// use delaunay::core::util::usize_to_u8;
///
/// // Successful conversion
/// assert_eq!(usize_to_u8(0, 4), Ok(0));
/// assert_eq!(usize_to_u8(255, 256), Ok(255));
///
/// // Failed conversion
/// let result = usize_to_u8(256, 10);
/// assert!(result.is_err());
/// ```
pub fn usize_to_u8(idx: usize, facet_count: usize) -> Result<u8, FacetError> {
    u8::try_from(idx).map_err(|_| FacetError::InvalidFacetIndexOverflow {
        original_index: idx,
        facet_count,
    })
}

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

    use crate::core::util::measure_with_result;
    use crate::vertex;

    use std::thread;
    use std::time::Instant;

    #[test]
    fn periodic_facet_key_rejects_invalid_lifted_cell_arity() {
        let lifted_vertices = vec![
            (VertexKey::null(), [0_i8; 2]),
            (VertexKey::null(), [0_i8; 2]),
        ];
        let err = periodic_facet_key_from_lifted_vertices::<2>(&lifted_vertices, 0).unwrap_err();
        assert!(matches!(
            err,
            PeriodicFacetKeyDerivationError::InvalidLiftedCellArity {
                expected: 3,
                actual: 2
            }
        ));
    }

    #[test]
    fn periodic_facet_key_rejects_relative_offset_out_of_range() {
        let lifted_vertices = vec![
            (VertexKey::null(), [-128_i8, 0_i8]),
            (VertexKey::null(), [0_i8, 0_i8]),
            (VertexKey::null(), [127_i8, 0_i8]),
        ];
        let err = periodic_facet_key_from_lifted_vertices::<2>(&lifted_vertices, 1).unwrap_err();
        assert!(matches!(
            err,
            PeriodicFacetKeyDerivationError::RelativeOffsetOutOfRange { axis: 0, .. }
        ));
    }
    #[test]
    fn periodic_facet_key_happy_path_translation_invariance_and_distinctness() {
        let lifted_base = vec![
            (VertexKey::null(), [0_i8, 0_i8]),
            (VertexKey::null(), [1_i8, 0_i8]),
            (VertexKey::null(), [0_i8, 1_i8]),
        ];
        let lifted_translated = vec![
            (VertexKey::null(), [7_i8, -3_i8]),
            (VertexKey::null(), [8_i8, -3_i8]),
            (VertexKey::null(), [7_i8, -2_i8]),
        ];
        let lifted_distinct = vec![
            (VertexKey::null(), [0_i8, 0_i8]),
            (VertexKey::null(), [2_i8, 0_i8]),
            (VertexKey::null(), [0_i8, 1_i8]),
        ];

        let base_key: Result<u64, PeriodicFacetKeyDerivationError> =
            periodic_facet_key_from_lifted_vertices::<2>(&lifted_base, 0);
        let translated_key: Result<u64, PeriodicFacetKeyDerivationError> =
            periodic_facet_key_from_lifted_vertices::<2>(&lifted_translated, 0);
        let distinct_key: Result<u64, PeriodicFacetKeyDerivationError> =
            periodic_facet_key_from_lifted_vertices::<2>(&lifted_distinct, 0);

        assert!(base_key.is_ok());
        assert!(translated_key.is_ok());
        assert!(distinct_key.is_ok());
        let base_key = base_key.unwrap();
        let translated_key = translated_key.unwrap();
        let distinct_key = distinct_key.unwrap();
        assert_eq!(
            base_key, translated_key,
            "uniform global translation should preserve periodic facet key",
        );
        assert_ne!(
            base_key, distinct_key,
            "different relative offsets should produce a distinct periodic facet key",
        );
    }

    #[test]
    fn periodic_facet_key_rejects_out_of_bounds_facet_index() {
        let lifted_vertices = vec![
            (VertexKey::null(), [0_i8; 2]),
            (VertexKey::null(), [0_i8; 2]),
            (VertexKey::null(), [0_i8; 2]),
        ];
        let err = periodic_facet_key_from_lifted_vertices::<2>(&lifted_vertices, 3).unwrap_err();
        assert!(matches!(
            err,
            PeriodicFacetKeyDerivationError::FacetIndexOutOfBounds {
                facet_index: 3,
                vertex_count: 3
            }
        ));
    }

    #[test]
    #[expect(clippy::too_many_lines)]
    fn test_checked_facet_key_from_vertex_keys_comprehensive() {
        println!("Testing checked_facet_key_from_vertex_keys comprehensively");

        // Create a triangulation
        let vertices = vec![
            vertex!([0.0, 0.0, 0.0]),
            vertex!([1.0, 0.0, 0.0]),
            vertex!([0.0, 1.0, 0.0]),
            vertex!([0.0, 0.0, 1.0]),
        ];
        let dt =
            crate::core::delaunay_triangulation::DelaunayTriangulation::new(&vertices).unwrap();
        let tds = &dt.as_triangulation().tds;

        // Test 1: Basic functionality - successful key derivation
        println!("  Testing basic functionality...");
        let cell = tds.cells().map(|(_, cell)| cell).next().unwrap();
        let facet_vertex_keys: Vec<_> = cell.vertices().iter().skip(1).copied().collect();

        let result = checked_facet_key_from_vertex_keys::<3>(&facet_vertex_keys);
        assert!(
            result.is_ok(),
            "Facet key derivation should succeed for valid vertex keys"
        );

        let facet_key = result.unwrap();
        println!("    Derived facet key: {facet_key}");

        // Test deterministic behavior - same vertex keys produce same key
        let result2 = checked_facet_key_from_vertex_keys::<3>(&facet_vertex_keys);
        assert!(result2.is_ok(), "Second derivation should also succeed");
        assert_eq!(
            facet_key,
            result2.unwrap(),
            "Same vertex keys should produce same facet key"
        );

        // Test different vertex keys produce different keys
        let all_vertex_keys = cell.vertices();
        let different_facet_vertex_keys: Vec<_> = all_vertex_keys.iter().take(3).copied().collect();
        if different_facet_vertex_keys.len() == 3
            && different_facet_vertex_keys != facet_vertex_keys
        {
            let result3 = checked_facet_key_from_vertex_keys::<3>(&different_facet_vertex_keys);
            assert!(
                result3.is_ok(),
                "Different facet key derivation should succeed"
            );
            let different_facet_key = result3.unwrap();
            assert_ne!(
                facet_key, different_facet_key,
                "Different vertex keys should produce different facet keys"
            );
            println!("    Different facet key: {different_facet_key}");
        }

        // Test 2: Error cases
        println!("  Testing error handling...");

        // Wrong vertex key count
        let single_key: Vec<VertexKey> = vec![facet_vertex_keys[0]];
        let result_count = checked_facet_key_from_vertex_keys::<3>(&single_key);
        assert!(
            result_count.is_err(),
            "Should return error for wrong vertex key count"
        );
        if let Err(error) = result_count {
            match error {
                FacetError::InsufficientVertices {
                    expected,
                    actual,
                    dimension,
                } => {
                    assert_eq!(expected, 3, "Expected 3 vertex keys for 3D");
                    assert_eq!(actual, 1, "Got 1 vertex key");
                    assert_eq!(dimension, 3, "Dimension should be 3");
                }
                _ => panic!("Expected InsufficientVertices error, got: {error:?}"),
            }
        }

        // Empty vertex keys
        let empty_keys: Vec<VertexKey> = vec![];
        let result_empty = checked_facet_key_from_vertex_keys::<3>(&empty_keys);
        assert!(
            result_empty.is_err(),
            "Empty vertex keys should fail validation"
        );
        if let Err(error) = result_empty {
            match error {
                FacetError::InsufficientVertices {
                    expected,
                    actual,
                    dimension,
                } => {
                    assert_eq!(expected, 3, "Expected 3 vertex keys for 3D");
                    assert_eq!(actual, 0, "Got 0 vertex keys");
                    assert_eq!(dimension, 3, "Dimension should be 3");
                }
                _ => {
                    panic!(
                        "Expected InsufficientVertices error for empty vertex keys, got: {error:?}"
                    )
                }
            }
        }

        // Test 3: Consistency with TDS cache
        println!("  Testing consistency with TDS...");
        let cache = tds
            .build_facet_to_cells_map()
            .expect("Should build facet map in test");
        let mut keys_found = 0;
        let mut keys_tested = 0;

        for cell in tds.cells().map(|(_, cell)| cell) {
            let cell_vertex_keys = cell.vertices();
            for skip_vertex_idx in 0..cell_vertex_keys.len() {
                let facet_vertex_keys: Vec<_> = cell_vertex_keys
                    .iter()
                    .enumerate()
                    .filter(|(i, _)| *i != skip_vertex_idx)
                    .map(|(_, &vk)| vk)
                    .collect();

                if !facet_vertex_keys.is_empty() {
                    let key_result = checked_facet_key_from_vertex_keys::<3>(&facet_vertex_keys);
                    if let Ok(derived_key) = key_result {
                        keys_tested += 1;
                        if cache.contains_key(&derived_key) {
                            keys_found += 1;
                        }
                    }
                }
            }
        }

        println!("    Found {keys_found}/{keys_tested} derived keys in TDS cache");
        assert!(keys_tested > 0, "Should have tested some keys");
        println!("  ✓ All facet key derivation tests passed");
    }

    #[test]
    fn test_verify_facet_index_consistency_true_false_and_error_cases() {
        // True case: comparing a cell to itself.
        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 =
            crate::core::delaunay_triangulation::DelaunayTriangulation::new(&vertices).unwrap();
        let tds = &dt.as_triangulation().tds;
        let cell_key = tds.cell_keys().next().unwrap();
        assert!(verify_facet_index_consistency(tds, cell_key, cell_key, 0).unwrap());

        // Error case: facet index out of bounds.
        let err = verify_facet_index_consistency(tds, cell_key, cell_key, 99).unwrap_err();
        assert!(matches!(err, FacetError::InvalidFacetIndex { .. }));

        // Logging: demonstrate behavior for large out-of-bounds facet index
        let err_large = verify_facet_index_consistency(tds, cell_key, cell_key, 300).unwrap_err();
        println!("    Large facet_idx=300 error: {err_large:?}");
        assert!(matches!(err_large, FacetError::InvalidFacetIndex { .. }));

        // False case: two disjoint triangles in the same TDS share no facet keys.
        let mut tds2: Tds<f64, (), (), 2> = Tds::empty();
        let v_a = tds2
            .insert_vertex_with_mapping(vertex!([0.0, 0.0]))
            .unwrap();
        let v_b = tds2
            .insert_vertex_with_mapping(vertex!([1.0, 0.0]))
            .unwrap();
        let v_c = tds2
            .insert_vertex_with_mapping(vertex!([0.0, 1.0]))
            .unwrap();
        let v_d = tds2
            .insert_vertex_with_mapping(vertex!([10.0, 10.0]))
            .unwrap();
        let v_e = tds2
            .insert_vertex_with_mapping(vertex!([11.0, 10.0]))
            .unwrap();
        let v_f = tds2
            .insert_vertex_with_mapping(vertex!([10.0, 11.0]))
            .unwrap();

        let c1 = tds2
            .insert_cell_with_mapping(
                crate::core::cell::Cell::new(vec![v_a, v_b, v_c], None).unwrap(),
            )
            .unwrap();
        let c2 = tds2
            .insert_cell_with_mapping(
                crate::core::cell::Cell::new(vec![v_d, v_e, v_f], None).unwrap(),
            )
            .unwrap();

        assert!(!verify_facet_index_consistency(&tds2, c1, c2, 0).unwrap());
    }

    #[test]
    fn test_usize_to_u8_conversion_comprehensive() {
        // Sub-test: Successful conversions
        assert_eq!(usize_to_u8(0, 4), Ok(0));
        assert_eq!(usize_to_u8(1, 4), Ok(1));
        assert_eq!(usize_to_u8(255, 256), Ok(255));
        assert_eq!(usize_to_u8(u8::MAX as usize, 256), Ok(u8::MAX));

        // Sub-test: All valid u8 values (0..=255)
        for i in 0u8..=255 {
            let result = usize_to_u8(i as usize, 256);
            assert_eq!(result, Ok(i), "Failed to convert {i}");
        }

        // Sub-test: Edge cases around u8::MAX boundary
        assert_eq!(usize_to_u8(254, 300), Ok(254));
        assert_eq!(usize_to_u8(255, 300), Ok(255));
        assert!(usize_to_u8(256, 300).is_err());
        assert!(usize_to_u8(257, 300).is_err());

        // Sub-test: Different facet_count values
        let facet_counts = [1, 10, 100, 255, 256, 1000, usize::MAX];
        for &count in &facet_counts {
            assert_eq!(
                usize_to_u8(0, count),
                Ok(0),
                "Valid conversion should succeed"
            );
            let result_invalid = usize_to_u8(256, count);
            assert!(result_invalid.is_err(), "Invalid conversion should fail");
            if let Err(FacetError::InvalidFacetIndexOverflow { facet_count, .. }) = result_invalid {
                assert_eq!(facet_count, count);
            }
        }

        // Sub-test: Deterministic behavior
        for i in 0..10 {
            assert_eq!(
                usize_to_u8(i, 20),
                usize_to_u8(i, 20),
                "Should be deterministic"
            );
        }
        for i in [256, 1000, usize::MAX] {
            assert_eq!(
                usize_to_u8(i, 100),
                usize_to_u8(i, 100),
                "Error cases should be deterministic"
            );
        }
    }

    #[test]
    fn test_usize_to_u8_error_handling() {
        // Sub-test: Basic error cases
        let result = usize_to_u8(256, 10);
        assert!(result.is_err());
        if let Err(FacetError::InvalidFacetIndexOverflow {
            original_index,
            facet_count,
        }) = result
        {
            assert_eq!(original_index, 256);
            assert_eq!(facet_count, 10);
        } else {
            panic!("Expected InvalidFacetIndexOverflow error");
        }

        let result = usize_to_u8(usize::MAX, 5);
        assert!(result.is_err());
        if let Err(FacetError::InvalidFacetIndexOverflow {
            original_index,
            facet_count,
        }) = result
        {
            assert_eq!(original_index, usize::MAX);
            assert_eq!(facet_count, 5);
        } else {
            panic!("Expected InvalidFacetIndexOverflow error");
        }

        // Sub-test: Error consistency across various inputs
        let test_cases = [
            (256, 100),
            (300, 500),
            (1000, 1500),
            (usize::MAX, 42),
            (65536, 10),
        ];
        for &(idx, count) in &test_cases {
            let result = usize_to_u8(idx, count);
            assert!(result.is_err(), "Should fail for index {idx}");
            match result {
                Err(FacetError::InvalidFacetIndexOverflow {
                    original_index,
                    facet_count,
                }) => {
                    assert_eq!(original_index, idx, "Should preserve original index");
                    assert_eq!(facet_count, count, "Should preserve facet_count");
                }
                _ => panic!("Expected InvalidFacetIndex error for index {idx}"),
            }
        }

        // Sub-test: Large values
        let large_values = [257, 1000, 10000, 65536, usize::MAX];
        for &val in &large_values {
            let result = usize_to_u8(val, val);
            assert!(result.is_err(), "Should fail for value {val}");
            if let Err(FacetError::InvalidFacetIndexOverflow {
                original_index,
                facet_count,
            }) = result
            {
                assert_eq!(original_index, val);
                assert_eq!(facet_count, val);
            }
        }

        // Sub-test: Error message quality
        let result = usize_to_u8(300, 42);
        assert!(result.is_err());
        if let Err(error) = result {
            let error_string = format!("{error}");
            assert!(
                error_string.contains("InvalidFacetIndex") || error_string.contains("index"),
                "Error message should indicate invalid index: {error_string}"
            );
        }

        let result = usize_to_u8(usize::MAX, 7);
        if let Err(FacetError::InvalidFacetIndexOverflow {
            original_index,
            facet_count,
        }) = result
        {
            assert_eq!(original_index, usize::MAX);
            assert_eq!(facet_count, 7);
        } else {
            panic!("Expected InvalidFacetIndexOverflow error");
        }
    }

    #[test]
    fn test_usize_to_u8_performance_and_threading() {
        // Sub-test: Performance characteristics (informational only)
        let start = Instant::now();
        for i in 0..1000 {
            let _ = usize_to_u8(i % 256, 300);
        }
        let duration = start.elapsed();
        eprintln!("usize_to_u8 valid conversions: 1000 iters in {duration:?}");

        let start = Instant::now();
        for i in 256..1256 {
            let _ = usize_to_u8(i, 100);
        }
        let duration = start.elapsed();
        eprintln!("usize_to_u8 error conversions: 1000 iters in {duration:?}");

        // Sub-test: Memory efficiency (stack allocation only)
        let (result, _alloc_info) = measure_with_result(|| {
            let mut results = Vec::new();
            for i in 0..100 {
                results.push(usize_to_u8(i, 200));
            }
            results
        });
        for (i, result) in result.iter().enumerate() {
            assert_eq!(
                *result,
                Ok(u8::try_from(i).unwrap()),
                "Result should be correct for {i}"
            );
        }

        // Sub-test: Thread safety
        let handles: Vec<_> = (0..4)
            .map(|thread_id| {
                thread::spawn(move || {
                    let mut results = Vec::new();
                    for i in 0..100 {
                        let val = (thread_id * 50 + i) % 256;
                        results.push(usize_to_u8(val, 300));
                    }
                    results
                })
            })
            .collect();
        for handle in handles {
            let results = handle.join().expect("Thread should complete successfully");
            for result in results {
                assert!(result.is_ok(), "All results should be successful");
            }
        }
    }
}