#![forbid(unsafe_code)]
use super::collections::{MAX_PRACTICAL_DIMENSION_SIZE, SmallBuffer};
use super::traits::data_type::DataType;
use super::util::{stable_hash_u64_slice, usize_to_u8};
use super::{
cell::Cell,
triangulation_data_structure::{CellKey, Tds, TdsValidationError, VertexKey},
vertex::Vertex,
};
use crate::geometry::traits::coordinate::CoordinateScalar;
use slotmap::Key;
use std::fmt::{self, Debug};
use std::sync::Arc;
use thiserror::Error;
#[derive(Clone, Debug, Error, PartialEq, Eq)]
#[non_exhaustive]
pub enum FacetError {
#[error("The cell does not contain the vertex!")]
CellDoesNotContainVertex,
#[error("Vertex UUID not found in mapping: {uuid}")]
VertexNotFound {
uuid: uuid::Uuid,
},
#[error(
"Facet must have exactly {expected} vertices for {dimension}D triangulation, got {actual}"
)]
InsufficientVertices {
expected: usize,
actual: usize,
dimension: usize,
},
#[error("Facet not found in triangulation")]
FacetNotFoundInTriangulation,
#[error(
"Facet key {facet_key:016x} not found in cache with {cache_size} entries - possible invariant violation or key derivation mismatch. Vertex UUIDs: {vertex_uuids:?}"
)]
FacetKeyNotFoundInCache {
facet_key: u64,
cache_size: usize,
vertex_uuids: Vec<uuid::Uuid>,
},
#[error("Expected exactly 1 adjacent cell for boundary facet, found {found}")]
InvalidAdjacentCellCount {
found: usize,
},
#[error("Adjacent cell not found")]
AdjacentCellNotFound,
#[error("Could not find inside vertex for boundary facet")]
InsideVertexNotFound,
#[error("Failed to compute orientation: {details}")]
OrientationComputationFailed {
details: String,
},
#[error("Invalid facet index {index} for cell with {facet_count} facets")]
InvalidFacetIndex {
index: u8,
facet_count: usize,
},
#[error(
"Invalid facet index {original_index} (too large for u8 conversion) for {facet_count} facets"
)]
InvalidFacetIndexOverflow {
original_index: usize,
facet_count: usize,
},
#[error("Cell not found in triangulation (potential data corruption)")]
CellNotFoundInTriangulation,
#[error("Vertex key not found in triangulation: {key:?}")]
VertexKeyNotFoundInTriangulation {
key: VertexKey,
},
#[error(
"Facet with key {facet_key:016x} has invalid multiplicity {found}, expected 1 (boundary) or 2 (internal)"
)]
InvalidFacetMultiplicity {
facet_key: u64,
found: usize,
},
#[error("Failed to retrieve boundary facets: {source}")]
BoundaryFacetRetrievalFailed {
#[source]
source: Arc<TdsValidationError>,
},
#[error("Cell operation failed: {source}")]
CellOperationFailed {
#[source]
source: Arc<TdsValidationError>,
},
}
#[derive(Copy, Clone, Debug, PartialEq, Eq, Hash)]
pub struct FacetHandle {
cell_key: CellKey,
facet_index: u8,
}
impl FacetHandle {
#[must_use]
pub const fn new(cell_key: CellKey, facet_index: u8) -> Self {
Self {
cell_key,
facet_index,
}
}
#[must_use]
pub const fn cell_key(&self) -> CellKey {
self.cell_key
}
#[must_use]
pub const fn facet_index(&self) -> u8 {
self.facet_index
}
}
pub struct FacetView<'tds, T, U, V, const D: usize>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
tds: &'tds Tds<T, U, V, D>,
cell_key: CellKey,
facet_index: u8,
}
impl<'tds, T, U, V, const D: usize> FacetView<'tds, T, U, V, D>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
#[inline]
#[must_use]
pub const fn cell_key(&self) -> CellKey {
self.cell_key
}
#[inline]
#[must_use]
pub const fn facet_index(&self) -> u8 {
self.facet_index
}
#[inline]
#[must_use]
pub const fn tds(&self) -> &'tds Tds<T, U, V, D> {
self.tds
}
pub fn new(
tds: &'tds Tds<T, U, V, D>,
cell_key: CellKey,
facet_index: u8,
) -> Result<Self, FacetError> {
let cell = tds
.get_cell(cell_key)
.ok_or(FacetError::CellNotFoundInTriangulation)?;
let vertex_count = cell.number_of_vertices();
if usize::from(facet_index) >= vertex_count {
return Err(FacetError::InvalidFacetIndex {
index: facet_index,
facet_count: vertex_count,
});
}
Ok(Self {
tds,
cell_key,
facet_index,
})
}
pub fn vertices(&self) -> Result<impl Iterator<Item = &'tds Vertex<T, U, D>>, FacetError> {
let cell = self
.tds
.get_cell(self.cell_key)
.ok_or(FacetError::CellNotFoundInTriangulation)?;
let facet_index = usize::from(self.facet_index);
let mut refs: SmallBuffer<&'tds Vertex<T, U, D>, MAX_PRACTICAL_DIMENSION_SIZE> =
SmallBuffer::with_capacity(cell.number_of_vertices().saturating_sub(1));
for (i, &vkey) in cell.vertices().iter().enumerate() {
if i == facet_index {
continue;
}
refs.push(
self.tds
.get_vertex_by_key(vkey)
.ok_or(FacetError::VertexKeyNotFoundInTriangulation { key: vkey })?,
);
}
Ok(refs.into_iter())
}
pub fn opposite_vertex(&self) -> Result<&'tds Vertex<T, U, D>, FacetError> {
let cell = self
.tds
.get_cell(self.cell_key)
.ok_or(FacetError::CellNotFoundInTriangulation)?;
let vertices = cell.vertices();
let facet_index = usize::from(self.facet_index);
let vkey = vertices
.get(facet_index)
.ok_or(FacetError::InvalidFacetIndex {
index: self.facet_index,
facet_count: vertices.len(),
})?;
self.tds
.get_vertex_by_key(*vkey)
.ok_or(FacetError::VertexKeyNotFoundInTriangulation { key: *vkey })
}
pub fn cell(&self) -> Result<&'tds Cell<T, U, V, D>, FacetError> {
self.tds
.get_cell(self.cell_key)
.ok_or(FacetError::CellNotFoundInTriangulation)
}
pub fn key(&self) -> Result<u64, FacetError> {
self.tds
.facet_key_for_cell_facet(self.cell_key, usize::from(self.facet_index))
.map_err(|e| FacetError::CellOperationFailed {
source: Arc::new(e),
})
}
}
impl<T, U, V, const D: usize> Debug for FacetView<'_, T, U, V, D>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
f.debug_struct("FacetView")
.field("cell_key", &self.cell_key)
.field("facet_index", &self.facet_index)
.field("dimension", &D)
.finish()
}
}
#[expect(clippy::expl_impl_clone_on_copy)]
#[expect(clippy::non_canonical_clone_impl)]
impl<T, U, V, const D: usize> Clone for FacetView<'_, T, U, V, D>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
fn clone(&self) -> Self {
Self {
tds: self.tds,
cell_key: self.cell_key,
facet_index: self.facet_index,
}
}
}
impl<T, U, V, const D: usize> Copy for FacetView<'_, T, U, V, D>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
}
impl<T, U, V, const D: usize> PartialEq for FacetView<'_, T, U, V, D>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
fn eq(&self, other: &Self) -> bool {
std::ptr::eq(self.tds, other.tds)
&& self.cell_key == other.cell_key
&& self.facet_index == other.facet_index
}
}
impl<T, U, V, const D: usize> Eq for FacetView<'_, T, U, V, D>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
}
pub fn all_facets_for_cell<T, U, V, const D: usize>(
tds: &Tds<T, U, V, D>,
cell_key: CellKey,
) -> Result<Vec<FacetView<'_, T, U, V, D>>, FacetError>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
let cell = tds
.get_cell(cell_key)
.ok_or(FacetError::CellNotFoundInTriangulation)?;
let vertex_count = cell.number_of_vertices();
let mut facet_views = Vec::with_capacity(vertex_count);
for facet_index in 0..vertex_count {
let idx = facet_index; let facet_view = FacetView::new(tds, cell_key, usize_to_u8(idx, vertex_count)?)?;
facet_views.push(facet_view);
}
Ok(facet_views)
}
#[derive(Clone)]
pub struct AllFacetsIter<'tds, T, U, V, const D: usize>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
tds: &'tds Tds<T, U, V, D>,
cell_keys: std::vec::IntoIter<CellKey>,
current_cell_key: Option<CellKey>,
current_facet_index: usize,
current_cell_facet_count: usize,
}
impl<'tds, T, U, V, const D: usize> AllFacetsIter<'tds, T, U, V, D>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
#[must_use]
pub fn new(tds: &'tds Tds<T, U, V, D>) -> Self {
assert!(
D <= 255,
"Dimension D={D} exceeds maximum of 255 for u8 facet indices"
);
#[expect(clippy::needless_collect)]
let cell_keys: Vec<CellKey> = tds.cell_keys().collect();
Self {
tds,
cell_keys: cell_keys.into_iter(),
current_cell_key: None,
current_facet_index: 0,
current_cell_facet_count: 0,
}
}
}
impl<'tds, T, U, V, const D: usize> Iterator for AllFacetsIter<'tds, T, U, V, D>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
type Item = FacetView<'tds, T, U, V, D>;
fn next(&mut self) -> Option<Self::Item> {
loop {
if let Some(cell_key) = self.current_cell_key
&& self.current_facet_index < self.current_cell_facet_count
{
let facet_index = self.current_facet_index;
self.current_facet_index += 1;
let Ok(facet_u8) = usize_to_u8(facet_index, self.current_cell_facet_count) else {
return None;
};
if let Ok(facet_view) = FacetView::new(self.tds, cell_key, facet_u8) {
return Some(facet_view);
}
}
if let Some(next_cell_key) = self.cell_keys.next() {
if let Some(cell) = self.tds.get_cell(next_cell_key) {
self.current_cell_key = Some(next_cell_key);
self.current_facet_index = 0;
self.current_cell_facet_count = cell.number_of_vertices();
} else {
}
} else {
return None;
}
}
}
}
#[derive(Clone)]
pub struct BoundaryFacetsIter<'tds, T, U, V, const D: usize>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
all_facets: AllFacetsIter<'tds, T, U, V, D>,
facet_to_cells_map: crate::core::collections::FacetToCellsMap,
}
impl<'tds, T, U, V, const D: usize> BoundaryFacetsIter<'tds, T, U, V, D>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
#[must_use]
pub fn new(
tds: &'tds Tds<T, U, V, D>,
facet_to_cells_map: crate::core::collections::FacetToCellsMap,
) -> Self {
Self {
all_facets: AllFacetsIter::new(tds),
facet_to_cells_map,
}
}
}
impl<'tds, T, U, V, const D: usize> Iterator for BoundaryFacetsIter<'tds, T, U, V, D>
where
T: CoordinateScalar,
U: DataType,
V: DataType,
{
type Item = FacetView<'tds, T, U, V, D>;
fn next(&mut self) -> Option<Self::Item> {
self.all_facets.find(|facet_view| {
if let Ok(facet_key) = facet_view.key()
&& let Some(cell_list) = self.facet_to_cells_map.get(&facet_key)
{
return cell_list.len() == 1;
}
false
})
}
}
#[must_use]
pub fn facet_key_from_vertices(vertices: &[VertexKey]) -> u64 {
if vertices.is_empty() {
return 0;
}
let mut key_values: SmallBuffer<u64, MAX_PRACTICAL_DIMENSION_SIZE> =
vertices.iter().map(|key| key.data().as_ffi()).collect();
key_values.sort_unstable();
stable_hash_u64_slice(&key_values)
}
#[cfg(test)]
mod tests {
use super::*;
use crate::core::delaunay_triangulation::DelaunayTriangulation;
use crate::core::triangulation_data_structure::VertexKey;
use crate::vertex;
#[test]
fn test_usize_to_u8_conversion() {
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));
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");
}
}
#[test]
fn test_facet_error_handling() {
let vertices = vec![vertex!([0.0]), vertex!([1.0])];
let dt = DelaunayTriangulation::new(&vertices).unwrap();
let cell_key = dt.cells().next().unwrap().0;
assert!(matches!(
FacetView::new(dt.tds(), cell_key, 99),
Err(FacetError::InvalidFacetIndex { .. })
));
}
#[test]
fn facet_new() {
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 cell_key = dt.cells().next().unwrap().0;
let facet = FacetView::new(dt.tds(), cell_key, 0).unwrap();
assert_eq!(facet.cell_key(), cell_key);
assert_eq!(facet.facet_index(), 0);
println!(
"FacetView: cell_key={:?}, facet_index={}",
facet.cell_key(),
facet.facet_index()
);
}
#[test]
fn test_facet_new_success_coverage() {
let vertices_2d = vec![
vertex!([0.0, 0.0]),
vertex!([1.0, 0.0]),
vertex!([0.5, 1.0]),
];
let dt_2d = DelaunayTriangulation::new(&vertices_2d).unwrap();
let cell_key_2d = dt_2d.cells().next().unwrap().0;
let result_2d = FacetView::new(dt_2d.tds(), cell_key_2d, 0);
assert!(result_2d.is_ok());
let facet_2d = result_2d.unwrap();
assert_eq!(facet_2d.vertices().unwrap().count(), 2);
let vertices_1d = vec![vertex!([0.0]), vertex!([1.0])];
let dt_1d = DelaunayTriangulation::new(&vertices_1d).unwrap();
let cell_key_1d = dt_1d.cells().next().unwrap().0;
let result_1d = FacetView::new(dt_1d.tds(), cell_key_1d, 0);
assert!(result_1d.is_ok());
let facet_1d = result_1d.unwrap();
assert_eq!(facet_1d.vertices().unwrap().count(), 1); }
#[test]
fn facet_new_with_incorrect_vertex() {
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 cell_key = dt.cells().next().unwrap().0;
assert!(FacetView::new(dt.tds(), cell_key, 4).is_err());
}
#[test]
fn facet_vertices() {
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 cell_key = dt.cells().next().unwrap().0;
let facet = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet_vertices: Vec<_> = facet.vertices().unwrap().collect();
assert_eq!(facet_vertices.len(), 3);
println!(
"FacetView: facet_index={}, vertex_count={}",
facet.facet_index(),
facet_vertices.len()
);
}
#[test]
fn facet_partial_eq() {
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 cell_key = dt.cells().next().unwrap().0;
let facet1 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet2 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet3 = FacetView::new(dt.tds(), cell_key, 1).unwrap();
assert_eq!(facet1, facet2);
assert_ne!(facet1, facet3);
}
#[test]
fn facet_clone() {
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 cell_key = dt.cells().next().unwrap().0;
let facet = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let cloned_facet = facet;
assert_eq!(facet, cloned_facet);
assert_eq!(facet.cell_key(), cloned_facet.cell_key());
assert_eq!(facet.facet_index(), cloned_facet.facet_index());
let cell1 = facet.cell().unwrap();
let cell2 = cloned_facet.cell().unwrap();
assert_eq!(cell1.uuid(), cell2.uuid());
let vertex1 = facet.opposite_vertex().unwrap();
let vertex2 = cloned_facet.opposite_vertex().unwrap();
assert_eq!(vertex1.uuid(), vertex2.uuid());
}
#[test]
fn facet_debug() {
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 cell_key = dt.cells().next().unwrap().0;
let facet = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let debug_str = format!("{facet:?}");
assert!(debug_str.contains("FacetView"));
assert!(debug_str.contains("cell_key"));
assert!(debug_str.contains("facet_index"));
assert!(debug_str.contains("dimension"));
}
#[test]
fn facet_with_typed_data() {
use crate::core::vertex::Vertex;
use crate::geometry::kernel::FastKernel;
let vertices: Vec<Vertex<f64, i32, 3>> = vec![
vertex!([0.0, 0.0, 0.0], 1),
vertex!([1.0, 0.0, 0.0], 2),
vertex!([0.0, 1.0, 0.0], 3),
vertex!([0.0, 0.0, 1.0], 4),
];
let dt: DelaunayTriangulation<FastKernel<f64>, i32, (), 3> =
DelaunayTriangulation::with_kernel(&FastKernel::new(), &vertices).unwrap();
let cell_key = dt.cells().next().unwrap().0;
let facet = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet_vertices: Vec<_> = facet.vertices().unwrap().collect();
assert_eq!(facet_vertices.len(), 3); assert!(facet_vertices.iter().any(|v| v.data == Some(2)));
assert!(facet_vertices.iter().any(|v| v.data == Some(3)));
assert!(facet_vertices.iter().any(|v| v.data == Some(4)));
}
macro_rules! test_facet_dimensions {
($(
$test_name:ident => $dim:expr => $desc:expr => $expected_facet_vertices:expr => $vertices:expr
),+ $(,)?) => {
$(
#[test]
fn $test_name() {
let vertices = $vertices;
let dt = DelaunayTriangulation::new(&vertices).unwrap();
let cell_key = dt.cells().next().unwrap().0;
let facet = FacetView::new(dt.tds(), cell_key, 0).unwrap();
assert_eq!(facet.vertices().unwrap().count(), $expected_facet_vertices,
"Facet of {}D {} should have {} vertices", $dim, $desc, $expected_facet_vertices);
}
pastey::paste! {
#[test]
fn [<$test_name _key_consistency>]() {
let vertices = $vertices;
let dt = DelaunayTriangulation::new(&vertices).unwrap();
let cell_key = dt.cells().next().unwrap().0;
let facet1 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet2 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
assert_eq!(facet1.key().unwrap(), facet2.key().unwrap(),
"Same facet should produce same key");
let facet3 = FacetView::new(dt.tds(), cell_key, 1).unwrap();
assert_ne!(facet1.key().unwrap(), facet3.key().unwrap(),
"Different facets should produce different keys");
}
#[test]
fn [<$test_name _equality>]() {
let vertices = $vertices;
let dt = DelaunayTriangulation::new(&vertices).unwrap();
let cell_key = dt.cells().next().unwrap().0;
let facet1 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet2 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet3 = FacetView::new(dt.tds(), cell_key, 1).unwrap();
assert_eq!(facet1, facet2, "Same facet should be equal");
assert_ne!(facet1, facet3, "Different facets should not be equal");
}
#[test]
fn [<$test_name _all_facets>]() {
let vertices = $vertices;
let dt = DelaunayTriangulation::new(&vertices).unwrap();
let cell_key = dt.cells().next().unwrap().0;
let expected_facets = $dim + 1;
let mut facet_keys = std::collections::HashSet::new();
for i in 0..expected_facets {
let facet = FacetView::new(dt.tds(), cell_key, u8::try_from(i).unwrap()).unwrap();
facet_keys.insert(facet.key().unwrap());
}
assert_eq!(facet_keys.len(), expected_facets,
"{}D cell should have {} unique facets", $dim, expected_facets);
}
}
)+
};
}
test_facet_dimensions! {
facet_2d_triangle => 2 => "triangle" => 2 => vec![
vertex!([0.0, 0.0]),
vertex!([1.0, 0.0]),
vertex!([0.5, 1.0]),
],
facet_3d_tetrahedron => 3 => "tetrahedron" => 3 => 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]),
],
facet_4d_simplex => 4 => "4-simplex" => 4 => vec![
vertex!([0.0, 0.0, 0.0, 0.0]),
vertex!([1.0, 0.0, 0.0, 0.0]),
vertex!([0.0, 1.0, 0.0, 0.0]),
vertex!([0.0, 0.0, 1.0, 0.0]),
vertex!([0.0, 0.0, 0.0, 1.0]),
],
facet_5d_simplex => 5 => "5-simplex" => 5 => vec![
vertex!([0.0, 0.0, 0.0, 0.0, 0.0]),
vertex!([1.0, 0.0, 0.0, 0.0, 0.0]),
vertex!([0.0, 1.0, 0.0, 0.0, 0.0]),
vertex!([0.0, 0.0, 1.0, 0.0, 0.0]),
vertex!([0.0, 0.0, 0.0, 1.0, 0.0]),
vertex!([0.0, 0.0, 0.0, 0.0, 1.0]),
],
}
#[test]
fn facet_1d_edge() {
let vertices = vec![vertex!([0.0]), vertex!([1.0])];
let dt = DelaunayTriangulation::new(&vertices).unwrap();
let cell_key = dt.cells().next().unwrap().0;
let facet = FacetView::new(dt.tds(), cell_key, 0).unwrap();
assert_eq!(facet.vertices().unwrap().count(), 1);
}
#[test]
fn facet_error_display() {
let cell_error = FacetError::CellDoesNotContainVertex;
assert_eq!(
cell_error.to_string(),
"The cell does not contain the vertex!"
);
}
#[test]
fn facet_error_debug() {
let cell_error = FacetError::CellDoesNotContainVertex;
let cell_debug = format!("{cell_error:?}");
assert!(cell_debug.contains("CellDoesNotContainVertex"));
}
#[test]
fn test_facet_key_consistency() {
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 cell_key = dt.cells().next().unwrap().0;
let facet1 = FacetView::new(dt.tds(), cell_key, 0).unwrap(); let facet2 = FacetView::new(dt.tds(), cell_key, 0).unwrap(); let facet3 = FacetView::new(dt.tds(), cell_key, 1).unwrap();
assert_eq!(
facet1.key().unwrap(),
facet2.key().unwrap(),
"Keys should be consistent for the same facet"
);
assert_ne!(
facet1.key().unwrap(),
facet3.key().unwrap(),
"Keys should be different for facets with different vertices"
);
}
#[test]
fn facet_vertices_empty_cell() {
let vertices = vec![vertex!([0.0]), vertex!([1.0])];
let dt = DelaunayTriangulation::new(&vertices).unwrap();
let cell_key = dt.cells().next().unwrap().0;
let facet = FacetView::new(dt.tds(), cell_key, 0).unwrap();
assert_eq!(facet.vertices().unwrap().count(), 1);
let other_facet = FacetView::new(dt.tds(), cell_key, 1).unwrap();
assert_eq!(other_facet.vertices().unwrap().count(), 1);
}
#[test]
fn facet_vertices_ordering() {
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 cell_key = dt.cells().next().unwrap().0;
let facet = FacetView::new(dt.tds(), cell_key, 2).unwrap();
assert_eq!(facet.vertices().unwrap().count(), 3);
}
#[test]
fn facet_eq_different_vertices() {
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 cell_key = dt.cells().next().unwrap().0;
let facet1 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet2 = FacetView::new(dt.tds(), cell_key, 1).unwrap();
let facet3 = FacetView::new(dt.tds(), cell_key, 2).unwrap();
let facet4 = FacetView::new(dt.tds(), cell_key, 3).unwrap();
assert_ne!(facet1, facet2);
assert_ne!(facet1, facet3);
assert_ne!(facet1, facet4);
assert_ne!(facet2, facet3);
assert_ne!(facet2, facet4);
assert_ne!(facet3, facet4);
}
#[test]
fn facet_key_hash() {
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 cell_key = dt.cells().next().unwrap().0;
let facet1 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet2 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet3 = FacetView::new(dt.tds(), cell_key, 1).unwrap();
assert_eq!(facet1.key().unwrap(), facet2.key().unwrap());
assert_ne!(facet1.key().unwrap(), facet3.key().unwrap());
}
#[test]
fn test_facet_key_from_vertices() {
use slotmap::SlotMap;
let mut temp_vertices: SlotMap<VertexKey, ()> = SlotMap::with_key();
let vertices = vec![
temp_vertices.insert(()),
temp_vertices.insert(()),
temp_vertices.insert(()),
];
let key1 = facet_key_from_vertices(&vertices);
let mut reversed_keys = vertices;
reversed_keys.reverse();
let key2 = facet_key_from_vertices(&reversed_keys);
assert_eq!(
key1, key2,
"Facet keys should be identical for the same vertices in different order"
);
let different_keys = vec![
temp_vertices.insert(()),
temp_vertices.insert(()),
temp_vertices.insert(()),
];
let key3 = facet_key_from_vertices(&different_keys);
assert_ne!(
key1, key3,
"Different vertices should produce different keys"
);
let empty_keys: Vec<VertexKey> = vec![];
let key_empty = facet_key_from_vertices(&empty_keys);
assert_eq!(key_empty, 0, "Empty vertex keys should produce key 0");
}
#[test]
fn test_facet_view_creation() {
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 cell_key = dt.cells().next().unwrap().0;
let facet_view = FacetView::new(dt.tds(), cell_key, 0).unwrap();
assert_eq!(facet_view.cell_key(), cell_key);
assert_eq!(facet_view.facet_index(), 0);
let result = FacetView::new(dt.tds(), cell_key, 10);
assert!(matches!(result, Err(FacetError::InvalidFacetIndex { .. })));
}
#[test]
fn test_facet_view_vertices_iteration() {
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 cell_key = dt.cells().next().unwrap().0;
let facet_view = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet_vertices: Vec<_> = facet_view.vertices().unwrap().collect();
assert_eq!(facet_vertices.len(), 3);
let original_vertices = vertices;
let opposite_vertex = &original_vertices[0];
assert!(
!facet_vertices
.iter()
.any(|v| v.uuid() == opposite_vertex.uuid())
);
}
#[test]
fn test_facet_view_opposite_vertex() {
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 cell_key = dt.cells().next().unwrap().0;
let facet_view = FacetView::new(dt.tds(), cell_key, 1).unwrap();
let opposite = facet_view.opposite_vertex().unwrap();
let cell = dt.tds().get_cell(cell_key).expect("cell exists");
let cell_vertex_keys = cell.vertices();
let expected_vertex = dt
.tds()
.get_vertex_by_key(cell_vertex_keys[1])
.expect("vertex exists");
assert_eq!(opposite.uuid(), expected_vertex.uuid());
}
#[test]
fn test_facet_view_key_computation() {
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 cell_key = dt.cells().next().unwrap().0;
let facet_view = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let key = facet_view.key().unwrap();
assert_ne!(key, 0);
}
#[test]
fn test_all_facets_for_cell() {
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 cell_key = dt.cells().next().unwrap().0;
let facet_views = all_facets_for_cell(dt.tds(), cell_key).unwrap();
assert_eq!(facet_views.len(), 4);
for (i, facet_view) in facet_views.iter().enumerate() {
assert_eq!(
facet_view.facet_index(),
usize_to_u8(i, facet_views.len()).unwrap()
);
assert_eq!(facet_view.cell_key(), cell_key);
}
}
#[test]
fn test_facet_view_equality() {
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 cell_key = dt.cells().next().unwrap().0;
let facet_view1 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet_view2 = FacetView::new(dt.tds(), cell_key, 0).unwrap();
let facet_view3 = FacetView::new(dt.tds(), cell_key, 1).unwrap();
assert_eq!(facet_view1, facet_view2);
assert_ne!(facet_view1, facet_view3);
}
#[test]
fn test_facet_view_debug() {
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 cell_key = dt.cells().next().unwrap().0;
let facet_view = FacetView::new(dt.tds(), cell_key, 1).unwrap();
let debug_str = format!("{facet_view:?}");
assert!(debug_str.contains("FacetView"));
assert!(debug_str.contains("cell_key"));
assert!(debug_str.contains("facet_index"));
assert!(debug_str.contains("dimension"));
}
#[test]
fn test_facet_view_memory_efficiency() {
use std::mem;
let lightweight_size = mem::size_of::<FacetView<f64, (), (), 3>>();
println!("Lightweight FacetView size: {lightweight_size} bytes");
assert!(lightweight_size <= 24);
}
}