#![forbid(unsafe_code)]
use super::point_generation::{generate_random_points, generate_random_points_seeded};
use crate::core::delaunay_triangulation::{
ConstructionOptions, DelaunayTriangulation, DelaunayTriangulationConstructionError,
InsertionOrderStrategy, RetryPolicy,
};
use crate::core::traits::data_type::DataType;
use crate::core::triangulation::{TopologyGuarantee, TriangulationConstructionError};
use crate::core::vertex::{Vertex, VertexBuilder};
use crate::geometry::kernel::RobustKernel;
use crate::geometry::point::Point;
use crate::geometry::robust_predicates::config_presets;
use crate::geometry::traits::coordinate::{CoordinateScalar, ScalarAccumulative};
use rand::SeedableRng;
use rand::distr::uniform::SampleUniform;
use rand::rngs::StdRng;
use rand::seq::SliceRandom;
const RANDOM_TRIANGULATION_MAX_SHUFFLE_ATTEMPTS: usize = 6;
const RANDOM_TRIANGULATION_MAX_POINTSET_ATTEMPTS: usize = 6;
const RANDOM_TRIANGULATION_POINTSET_SEED_MIX: u64 = 0x9E37_79B9_7F4A_7C15;
fn validate_random_triangulation<K, U, V, const D: usize>(
dt: DelaunayTriangulation<K, U, V, D>,
) -> Result<DelaunayTriangulation<K, U, V, D>, DelaunayTriangulationConstructionError>
where
K: crate::geometry::kernel::Kernel<D>,
K::Scalar: ScalarAccumulative,
U: DataType,
V: DataType,
{
dt.as_triangulation().validate().map_err(|e| {
TriangulationConstructionError::GeometricDegeneracy {
message: format!("Random triangulation failed topology/Euler validation: {e}"),
}
})?;
Ok(dt)
}
fn random_triangulation_is_acceptable<K, U, V, const D: usize>(
dt: &DelaunayTriangulation<K, U, V, D>,
min_vertices: usize,
) -> bool
where
K: crate::geometry::kernel::Kernel<D>,
K::Scalar: ScalarAccumulative,
U: DataType,
V: DataType,
{
dt.number_of_vertices() >= min_vertices
}
fn random_triangulation_try_build<K, T, U, V, const D: usize>(
kernel: &K,
vertices: &[Vertex<T, U, D>],
min_vertices: usize,
topology_guarantee: TopologyGuarantee,
) -> Option<DelaunayTriangulation<K, U, V, D>>
where
K: crate::geometry::kernel::Kernel<D, Scalar = T>,
T: ScalarAccumulative,
U: DataType,
V: DataType,
{
let options = ConstructionOptions::default()
.with_insertion_order(InsertionOrderStrategy::Input)
.with_retry_policy(RetryPolicy::Disabled);
let dt = DelaunayTriangulation::with_topology_guarantee_and_options(
kernel,
vertices,
topology_guarantee,
options,
)
.ok()?;
let dt = validate_random_triangulation(dt).ok()?;
random_triangulation_is_acceptable(&dt, min_vertices).then_some(dt)
}
fn random_triangulation_build_vertices<T, U, const D: usize>(
points: Vec<Point<T, D>>,
vertex_data: Option<U>,
) -> Vec<Vertex<T, U, D>>
where
T: CoordinateScalar,
U: DataType,
{
points
.into_iter()
.map(|point| {
vertex_data.map_or_else(
|| {
VertexBuilder::default()
.point(point)
.build()
.expect("Failed to build vertex without data")
},
|data| {
VertexBuilder::default()
.point(point)
.data(data)
.build()
.expect("Failed to build vertex with data")
},
)
})
.collect()
}
fn make_robust_kernel<T: ScalarAccumulative>() -> RobustKernel<T> {
RobustKernel::with_config(config_presets::degenerate_robust::<T>())
}
fn random_triangulation_try_with_vertices<T, U, V, const D: usize>(
vertices: &[Vertex<T, U, D>],
min_vertices: usize,
shuffle_seed: Option<u64>,
topology_guarantee: TopologyGuarantee,
) -> Option<DelaunayTriangulation<RobustKernel<T>, U, V, D>>
where
T: ScalarAccumulative,
U: DataType,
V: DataType,
{
let robust_kernel = make_robust_kernel::<T>();
if let Some(dt) =
random_triangulation_try_build(&robust_kernel, vertices, min_vertices, topology_guarantee)
{
return Some(dt);
}
for attempt in 0..RANDOM_TRIANGULATION_MAX_SHUFFLE_ATTEMPTS {
let mut shuffled = vertices.to_vec();
if let Some(seed_value) = shuffle_seed {
let mix = seed_value.wrapping_add(attempt as u64 + 1);
let mut rng = StdRng::seed_from_u64(mix);
shuffled.shuffle(&mut rng);
} else {
let mut rng = rand::rng();
shuffled.shuffle(&mut rng);
}
if let Some(dt) = random_triangulation_try_build(
&robust_kernel,
&shuffled,
min_vertices,
topology_guarantee,
) {
return Some(dt);
}
}
None
}
pub fn generate_random_triangulation<T, U, V, const D: usize>(
n_points: usize,
bounds: (T, T),
vertex_data: Option<U>,
seed: Option<u64>,
) -> Result<DelaunayTriangulation<RobustKernel<T>, U, V, D>, DelaunayTriangulationConstructionError>
where
T: ScalarAccumulative + SampleUniform,
U: DataType,
V: DataType,
{
#[cfg(debug_assertions)]
if std::env::var_os("DELAUNAY_DEBUG_UNUSED_IMPORTS").is_some() {
eprintln!(
"triangulation_generation::generate_random_triangulation called (n_points={n_points}, D={D}, seed={seed:?})"
);
}
generate_random_triangulation_with_topology_guarantee(
n_points,
bounds,
vertex_data,
seed,
TopologyGuarantee::DEFAULT,
)
}
#[allow(clippy::too_many_lines)]
pub fn generate_random_triangulation_with_topology_guarantee<T, U, V, const D: usize>(
n_points: usize,
bounds: (T, T),
vertex_data: Option<U>,
seed: Option<u64>,
topology_guarantee: TopologyGuarantee,
) -> Result<DelaunayTriangulation<RobustKernel<T>, U, V, D>, DelaunayTriangulationConstructionError>
where
T: ScalarAccumulative + SampleUniform,
U: DataType,
V: DataType,
{
if n_points == 0 {
return Ok(
DelaunayTriangulation::with_empty_kernel_and_topology_guarantee(
make_robust_kernel::<T>(),
topology_guarantee,
),
);
}
if n_points < D + 1 {
return Err(TriangulationConstructionError::InsufficientVertices {
dimension: D,
source: crate::core::cell::CellValidationError::InsufficientVertices {
actual: n_points,
expected: D + 1,
dimension: D,
},
}
.into());
}
let points: Vec<Point<T, D>> =
match seed {
Some(seed_value) => generate_random_points_seeded(n_points, bounds, seed_value)
.map_err(|e| TriangulationConstructionError::GeometricDegeneracy {
message: format!("Random point generation failed: {e}"),
})?,
None => generate_random_points(n_points, bounds).map_err(|e| {
TriangulationConstructionError::GeometricDegeneracy {
message: format!("Random point generation failed: {e}"),
}
})?,
};
let min_vertices = (n_points / 6).max(D + 1);
let mut initial_points = Some(points);
for attempt in 0..RANDOM_TRIANGULATION_MAX_POINTSET_ATTEMPTS {
#[cfg(debug_assertions)]
if std::env::var_os("DELAUNAY_DEBUG_RANDOM_POINTSET_RETRIES").is_some() {
eprintln!(
"random_triangulation: pointset attempt {attempt} of {RANDOM_TRIANGULATION_MAX_POINTSET_ATTEMPTS} (0-based)"
);
}
let point_seed = seed.map(|base| {
if attempt > 0 {
base ^ RANDOM_TRIANGULATION_POINTSET_SEED_MIX.wrapping_mul(attempt as u64)
} else {
base
}
});
let points = if attempt == 0 {
initial_points
.take()
.expect("initial points already consumed")
} else {
match point_seed {
Some(seed_value) => generate_random_points_seeded(n_points, bounds, seed_value)
.map_err(|e| TriangulationConstructionError::GeometricDegeneracy {
message: format!("Random point generation failed: {e}"),
})?,
None => generate_random_points(n_points, bounds).map_err(|e| {
TriangulationConstructionError::GeometricDegeneracy {
message: format!("Random point generation failed: {e}"),
}
})?,
}
};
let vertices = random_triangulation_build_vertices(points, vertex_data);
if let Some(dt) = random_triangulation_try_with_vertices(
&vertices,
min_vertices,
point_seed,
topology_guarantee,
) {
return Ok(dt);
}
}
Err(TriangulationConstructionError::GeometricDegeneracy {
message: "Random triangulation failed validation after robust fallback".to_string(),
}
.into())
}
pub struct RandomTriangulationBuilder<T> {
n_points: usize,
bounds: (T, T),
seed: Option<u64>,
topology_guarantee: TopologyGuarantee,
construction_options: ConstructionOptions,
}
impl<T> RandomTriangulationBuilder<T>
where
T: ScalarAccumulative + SampleUniform,
{
#[must_use]
pub fn new(n_points: usize, bounds: (T, T)) -> Self {
Self {
n_points,
bounds,
seed: None,
topology_guarantee: TopologyGuarantee::DEFAULT,
construction_options: ConstructionOptions::default(),
}
}
#[must_use]
pub const fn seed(mut self, seed: u64) -> Self {
self.seed = Some(seed);
self
}
#[must_use]
pub const fn topology_guarantee(mut self, topology_guarantee: TopologyGuarantee) -> Self {
self.topology_guarantee = topology_guarantee;
self
}
#[must_use]
pub const fn insertion_order(mut self, strategy: InsertionOrderStrategy) -> Self {
self.construction_options = self.construction_options.with_insertion_order(strategy);
self
}
#[must_use]
pub const fn construction_options(mut self, options: ConstructionOptions) -> Self {
self.construction_options = options;
self
}
pub fn build<U, V, const D: usize>(
self,
) -> Result<
DelaunayTriangulation<RobustKernel<T>, U, V, D>,
DelaunayTriangulationConstructionError,
>
where
U: DataType,
V: DataType,
{
self.build_with_vertex_data(None)
}
pub fn build_with_vertex_data<U, V, const D: usize>(
self,
vertex_data: Option<U>,
) -> Result<
DelaunayTriangulation<RobustKernel<T>, U, V, D>,
DelaunayTriangulationConstructionError,
>
where
U: DataType,
V: DataType,
{
if self.n_points == 0 {
return Ok(
DelaunayTriangulation::with_empty_kernel_and_topology_guarantee(
make_robust_kernel::<T>(),
self.topology_guarantee,
),
);
}
if self.n_points < D + 1 {
return Err(TriangulationConstructionError::InsufficientVertices {
dimension: D,
source: crate::core::cell::CellValidationError::InsufficientVertices {
actual: self.n_points,
expected: D + 1,
dimension: D,
},
}
.into());
}
let points: Vec<Point<T, D>> = match self.seed {
Some(seed_value) => {
generate_random_points_seeded(self.n_points, self.bounds, seed_value).map_err(
|e| TriangulationConstructionError::GeometricDegeneracy {
message: format!("Random point generation failed: {e}"),
},
)?
}
None => generate_random_points(self.n_points, self.bounds).map_err(|e| {
TriangulationConstructionError::GeometricDegeneracy {
message: format!("Random point generation failed: {e}"),
}
})?,
};
let vertices = random_triangulation_build_vertices(points, vertex_data);
#[cfg(debug_assertions)]
if std::env::var_os("DELAUNAY_DEBUG_RANDOM_BUILDER").is_some() {
eprintln!(
"random_triangulation_builder: single call to with_topology_guarantee_and_options with n_points={}, topology_guarantee={:?}, insertion_order={:?}, dedup_policy={:?}, retry_policy={:?}",
self.n_points,
self.topology_guarantee,
self.construction_options.insertion_order(),
self.construction_options.dedup_policy(),
self.construction_options.retry_policy(),
);
}
let dt = DelaunayTriangulation::with_topology_guarantee_and_options(
&make_robust_kernel::<T>(),
&vertices,
self.topology_guarantee,
self.construction_options,
)?;
validate_random_triangulation(dt)
}
}
#[cfg(test)]
mod tests {
use super::*;
use crate::vertex;
#[test]
#[ignore = "Flaky: unseeded random generation occasionally fails with cavity filling errors - needs investigation"]
fn test_generate_random_triangulation_basic() {
let triangulation_2d =
generate_random_triangulation::<f64, (), (), 2>(10, (-5.0, 5.0), None, Some(42))
.unwrap();
assert!(
triangulation_2d.number_of_vertices() >= 3,
"Expected at least 3 vertices in 2D triangulation, got {}",
triangulation_2d.number_of_vertices()
);
assert_eq!(triangulation_2d.dim(), 2);
let valid_2d = triangulation_2d.is_valid();
if let Err(e) = &valid_2d {
println!("test_generate_random_triangulation_basic (2D): TDS invalid: {e}");
}
assert!(valid_2d.is_ok());
let triangulation_3d =
generate_random_triangulation::<f64, i32, (), 3>(8, (0.0, 1.0), Some(123), Some(456))
.unwrap();
assert!(
triangulation_3d.number_of_vertices() >= 4,
"Expected at least 4 vertices in 3D triangulation, got {}",
triangulation_3d.number_of_vertices()
);
assert_eq!(triangulation_3d.dim(), 3);
let valid_3d = triangulation_3d.is_valid();
if let Err(e) = &valid_3d {
println!("test_generate_random_triangulation_basic (3D): TDS invalid: {e}");
}
assert!(valid_3d.is_ok());
let triangulation_seeded =
generate_random_triangulation::<f64, (), (), 2>(5, (-1.0, 1.0), None, Some(789))
.unwrap();
let triangulation_unseeded =
generate_random_triangulation::<f64, (), (), 2>(5, (-1.0, 1.0), None, None).unwrap();
let valid_seeded = triangulation_seeded.is_valid();
if let Err(e) = &valid_seeded {
println!("test_generate_random_triangulation_basic (seeded 2D): TDS invalid: {e}");
}
assert!(valid_seeded.is_ok());
let valid_unseeded = triangulation_unseeded.is_valid();
if let Err(e) = &valid_unseeded {
println!("test_generate_random_triangulation_basic (unseeded 2D): TDS invalid: {e}");
}
assert!(valid_unseeded.is_ok());
assert!(
triangulation_seeded.number_of_vertices() >= 3,
"Expected at least 3 vertices in seeded 2D triangulation, got {}",
triangulation_seeded.number_of_vertices()
);
assert!(
triangulation_unseeded.number_of_vertices() >= 3,
"Expected at least 3 vertices in unseeded 2D triangulation, got {}",
triangulation_unseeded.number_of_vertices()
);
}
#[test]
fn test_generate_random_triangulation_error_cases() {
let result = generate_random_triangulation::<f64, (), (), 2>(
10,
(5.0, 1.0), None,
Some(42),
);
assert!(result.is_err());
let result =
generate_random_triangulation::<f64, (), (), 2>(0, (-1.0, 1.0), None, Some(42));
assert!(result.is_ok());
let triangulation = result.unwrap();
assert_eq!(triangulation.number_of_vertices(), 0);
assert_eq!(triangulation.dim(), -1);
}
#[test]
fn test_generate_random_triangulation_reproducibility() {
let triangulation1 =
generate_random_triangulation::<f64, (), (), 3>(6, (-2.0, 2.0), None, Some(12345))
.unwrap();
let triangulation2 =
generate_random_triangulation::<f64, (), (), 3>(6, (-2.0, 2.0), None, Some(12345))
.unwrap();
assert_eq!(
triangulation1.number_of_vertices(),
triangulation2.number_of_vertices()
);
assert_eq!(
triangulation1.number_of_cells(),
triangulation2.number_of_cells()
);
assert_eq!(triangulation1.dim(), triangulation2.dim());
}
#[test]
fn test_random_triangulation_try_with_vertices_exercises_fallbacks() {
let vertices: Vec<Vertex<f64, (), 2>> = vec![
vertex!([0.0, 0.0]),
vertex!([1.0, 0.0]),
vertex!([0.0, 1.0]),
];
let result = random_triangulation_try_with_vertices::<f64, (), (), 2>(
&vertices,
vertices.len() + 1,
Some(7),
TopologyGuarantee::PLManifold,
);
assert!(result.is_none());
}
#[test]
#[ignore = "High-dimensional random triangulations occasionally hit cavity filling errors in CI - needs robust predicate investigation"]
fn test_generate_random_triangulation_dimensions() {
let tri_2d =
generate_random_triangulation::<f64, (), (), 2>(15, (0.0, 10.0), None, Some(555))
.unwrap();
assert_eq!(tri_2d.dim(), 2);
assert!(tri_2d.number_of_cells() > 0);
let tri_3d =
generate_random_triangulation::<f64, (), (), 3>(20, (-3.0, 3.0), None, Some(666))
.unwrap();
assert_eq!(tri_3d.dim(), 3);
assert!(tri_3d.number_of_cells() > 0);
let tri_4d =
generate_random_triangulation::<f64, (), (), 4>(12, (-1.0, 1.0), None, Some(777))
.unwrap();
assert_eq!(tri_4d.dim(), 4);
assert!(tri_4d.number_of_cells() > 0);
let tri_5d =
generate_random_triangulation::<f64, (), (), 5>(10, (0.0, 5.0), None, Some(888))
.unwrap();
assert_eq!(tri_5d.dim(), 5);
assert!(tri_5d.number_of_cells() > 0);
}
#[test]
fn test_generate_random_triangulation_with_data() {
let tri_with_char_array = generate_random_triangulation::<f64, [char; 8], (), 2>(
6,
(-2.0, 2.0),
Some(['v', 'e', 'r', 't', 'e', 'x', '_', 'd']),
Some(888),
)
.unwrap();
assert!(
tri_with_char_array.number_of_vertices() >= 3,
"Expected at least 3 vertices in 2D triangulation with data, got {}",
tri_with_char_array.number_of_vertices()
);
let valid_char = tri_with_char_array.tds().is_valid();
if let Err(e) = &valid_char {
println!(
"test_generate_random_triangulation_with_data (2D char data): TDS invalid: {e}"
);
}
assert!(valid_char.is_ok());
let char_array_data = ['v', 'e', 'r', 't', 'e', 'x', '_', 'd'];
let string_representation: String = char_array_data.iter().collect();
assert_eq!(string_representation, "vertex_d");
let seeds = [999_u64, 123, 456, 789, 2024];
let mut tri_with_int_data: Option<DelaunayTriangulation<RobustKernel<f64>, u32, (), 3>> =
None;
let mut last_err: Option<String> = None;
for seed in seeds {
match generate_random_triangulation::<f64, u32, (), 3>(
8,
(0.0, 5.0),
Some(42u32),
Some(seed),
) {
Ok(tri) => {
tri_with_int_data = Some(tri);
break;
}
Err(e) => {
last_err = Some(format!("{e}"));
}
}
}
let tri_with_int_data = tri_with_int_data.unwrap_or_else(|| {
panic!("All seeds failed to generate 3D triangulation with int data: {last_err:?}")
});
assert!(
tri_with_int_data.number_of_vertices() >= 4,
"Expected at least 4 vertices in 3D triangulation with data, got {}",
tri_with_int_data.number_of_vertices()
);
let valid_int = tri_with_int_data.tds().is_valid();
if let Err(e) = &valid_int {
println!(
"test_generate_random_triangulation_with_data (3D int data): TDS invalid: {e}"
);
}
assert!(valid_int.is_ok());
let tri_no_data =
generate_random_triangulation::<f64, (), (), 2>(5, (-1.0, 1.0), None, Some(111))
.unwrap();
assert!(
tri_no_data.number_of_vertices() >= 3,
"Expected at least 3 vertices in 2D triangulation without data, got {}",
tri_no_data.number_of_vertices()
);
let valid_no_data = tri_no_data.tds().is_valid();
if let Err(e) = &valid_no_data {
println!("test_generate_random_triangulation_with_data (2D no data): TDS invalid: {e}");
}
assert!(valid_no_data.is_ok());
}
}