1use crate::error::{SpatialError, SpatialResult};
41use qhull::Qh;
42use scirs2_core::ndarray::Array2;
43use std::collections::{HashMap, HashSet};
44use std::fmt::Debug;
45
46pub struct Delaunay {
55 points: Array2<f64>,
57
58 ndim: usize,
60
61 npoints: usize,
63
64 simplices: Vec<Vec<usize>>,
67
68 neighbors: Vec<Vec<i64>>,
72
73 #[allow(dead_code)]
75 qh: Option<Qh<'static>>,
76
77 constraints: Vec<(usize, usize)>,
80}
81
82impl Debug for Delaunay {
83 fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
84 f.debug_struct("Delaunay")
85 .field("points", &self.points.shape())
86 .field("ndim", &self.ndim)
87 .field("npoints", &self.npoints)
88 .field("simplices", &self.simplices.len())
89 .field("neighbors", &self.neighbors.len())
90 .field("constraints", &self.constraints.len())
91 .finish()
92 }
93}
94
95impl Clone for Delaunay {
96 fn clone(&self) -> Self {
97 Self {
98 points: self.points.clone(),
99 ndim: self.ndim,
100 npoints: self.npoints,
101 simplices: self.simplices.clone(),
102 neighbors: self.neighbors.clone(),
103 qh: None, constraints: self.constraints.clone(),
105 }
106 }
107}
108
109impl Delaunay {
110 pub fn new(points: &Array2<f64>) -> SpatialResult<Self> {
138 let npoints = points.nrows();
139 let ndim = points.ncols();
140
141 if npoints <= ndim {
143 return Err(SpatialError::ValueError(format!(
144 "Need at least {ndim_plus_1} _points in {ndim} dimensions for triangulation",
145 ndim_plus_1 = ndim + 1
146 )));
147 }
148
149 if ndim == 2 && npoints == 3 {
151 let simplex = vec![0, 1, 2];
152 let simplices = vec![simplex];
153 let neighbors = vec![vec![-1, -1, -1]]; return Ok(Delaunay {
156 points: points.clone(),
157 ndim,
158 npoints,
159 simplices,
160 neighbors,
161 qh: None,
162 constraints: Vec::new(),
163 });
164 }
165
166 let _points_vec: Vec<Vec<f64>> = (0..npoints).map(|i| points.row(i).to_vec()).collect();
168
169 let qh_result = Qh::new_delaunay(_points_vec.clone());
171
172 let qh = match qh_result {
173 Ok(qh) => qh,
174 Err(e) => {
175 if ndim == 2 && npoints == 4 {
177 let simplex1 = vec![0, 1, 2];
179 let simplex2 = vec![1, 2, 3];
180 let simplices = vec![simplex1, simplex2];
181 let neighbors = vec![vec![-1, 1, -1], vec![-1, -1, 0]];
182
183 return Ok(Delaunay {
184 points: points.clone(),
185 ndim,
186 npoints,
187 simplices,
188 neighbors,
189 qh: None,
190 constraints: Vec::new(),
191 });
192 }
193
194 let mut perturbed_points = vec![];
196 use scirs2_core::random::Rng;
197 let mut rng = rand::rng();
198
199 for i in 0..npoints {
200 let mut pt = points.row(i).to_vec();
201 for val in pt.iter_mut().take(ndim) {
202 *val += rng.gen_range(-0.0001..0.0001);
203 }
204 perturbed_points.push(pt);
205 }
206
207 match Qh::new_delaunay(perturbed_points) {
209 Ok(qh2) => qh2,
210 Err(_) => {
211 return Err(SpatialError::ComputationError(format!(
212 "Qhull error (even with perturbation): {e}"
213 )));
214 }
215 }
216 }
217 };
218
219 let simplices = Self::extract_simplices(&qh, ndim);
221
222 let neighbors = Self::calculate_neighbors(&simplices, ndim + 1);
224
225 Ok(Delaunay {
226 points: points.clone(),
227 ndim,
228 npoints,
229 simplices,
230 neighbors,
231 qh: Some(qh),
232 constraints: Vec::new(),
233 })
234 }
235
236 pub fn new_constrained(
274 points: &Array2<f64>,
275 constraints: Vec<(usize, usize)>,
276 ) -> SpatialResult<Self> {
277 let ndim = points.ncols();
278
279 if ndim != 2 && ndim != 3 {
282 return Err(SpatialError::NotImplementedError(
283 "Constrained Delaunay triangulation only supports 2D and 3D points".to_string(),
284 ));
285 }
286
287 let npoints = points.nrows();
289 for &(i, j) in &constraints {
290 if i >= npoints || j >= npoints {
291 return Err(SpatialError::ValueError(format!(
292 "Constraint edge ({i}, {j}) contains invalid point indices"
293 )));
294 }
295 if i == j {
296 return Err(SpatialError::ValueError(format!(
297 "Constraint edge ({i}, {j}) connects a point to itself"
298 )));
299 }
300 }
301
302 let mut delaunay = Self::new(points)?;
304 delaunay.constraints = constraints.clone();
305
306 delaunay.insert_constraints()?;
308
309 Ok(delaunay)
310 }
311
312 fn insert_constraints(&mut self) -> SpatialResult<()> {
314 for &(i, j) in &self.constraints.clone() {
315 self.insert_constraint_edge(i, j)?;
316 }
317 Ok(())
318 }
319
320 fn insert_constraint_edge(&mut self, start: usize, end: usize) -> SpatialResult<()> {
322 if self.edge_exists(start, end) {
324 return Ok(()); }
326
327 let intersecting_edges = self.find_intersecting_edges(start, end)?;
329
330 if intersecting_edges.is_empty() {
331 return Err(SpatialError::ComputationError(
333 "Constraint edge has no intersections but doesn't exist in triangulation"
334 .to_string(),
335 ));
336 }
337
338 let affected_triangles = self.find_triangles_with_edges(&intersecting_edges);
340 self.remove_triangles(&affected_triangles);
341
342 self.retriangulate_with_constraint(start, end, &affected_triangles)?;
344
345 Ok(())
346 }
347
348 fn edge_exists(&self, start: usize, end: usize) -> bool {
350 for simplex in &self.simplices {
351 let simplex_size = simplex.len();
352 for i in 0..simplex_size {
354 for j in (i + 1)..simplex_size {
355 let v1 = simplex[i];
356 let v2 = simplex[j];
357 if (v1 == start && v2 == end) || (v1 == end && v2 == start) {
358 return true;
359 }
360 }
361 }
362 }
363 false
364 }
365
366 fn find_intersecting_edges(
368 &self,
369 start: usize,
370 end: usize,
371 ) -> SpatialResult<Vec<(usize, usize)>> {
372 let mut intersecting = Vec::new();
373
374 let p1: Vec<f64> = self.points.row(start).to_vec();
376 let p2: Vec<f64> = self.points.row(end).to_vec();
377
378 let mut checked_edges = HashSet::new();
380
381 for simplex in &self.simplices {
382 let simplex_size = simplex.len();
383
384 for i in 0..simplex_size {
386 for j in (i + 1)..simplex_size {
387 let v1 = simplex[i];
388 let v2 = simplex[j];
389
390 let edge = if v1 < v2 { (v1, v2) } else { (v2, v1) };
392 if checked_edges.contains(&edge) {
393 continue;
394 }
395 checked_edges.insert(edge);
396
397 if v1 == start || v1 == end || v2 == start || v2 == end {
399 continue;
400 }
401
402 let q1: Vec<f64> = self.points.row(v1).to_vec();
403 let q2: Vec<f64> = self.points.row(v2).to_vec();
404
405 if self.ndim == 2 {
406 let p1_2d = [p1[0], p1[1]];
408 let p2_2d = [p2[0], p2[1]];
409 let q1_2d = [q1[0], q1[1]];
410 let q2_2d = [q2[0], q2[1]];
411
412 if Self::segments_intersect(p1_2d, p2_2d, q1_2d, q2_2d) {
413 intersecting.push((v1, v2));
414 }
415 } else if self.ndim == 3 {
416 if Self::edges_interfere_3d(&p1, &p2, &q1, &q2) {
419 intersecting.push((v1, v2));
420 }
421 }
422 }
423 }
424 }
425
426 Ok(intersecting)
427 }
428
429 fn segments_intersect(p1: [f64; 2], p2: [f64; 2], q1: [f64; 2], q2: [f64; 2]) -> bool {
431 fn orientation(p: [f64; 2], q: [f64; 2], r: [f64; 2]) -> i32 {
432 let val = (q[1] - p[1]) * (r[0] - q[0]) - (q[0] - p[0]) * (r[1] - q[1]);
433 if val.abs() < 1e-10 {
434 0
435 }
436 else if val > 0.0 {
438 1
439 }
440 else {
442 2
443 } }
445
446 fn on_segment(p: [f64; 2], q: [f64; 2], r: [f64; 2]) -> bool {
447 q[0] <= p[0].max(r[0])
448 && q[0] >= p[0].min(r[0])
449 && q[1] <= p[1].max(r[1])
450 && q[1] >= p[1].min(r[1])
451 }
452
453 let o1 = orientation(p1, p2, q1);
454 let o2 = orientation(p1, p2, q2);
455 let o3 = orientation(q1, q2, p1);
456 let o4 = orientation(q1, q2, p2);
457
458 if o1 != o2 && o3 != o4 {
460 return true;
461 }
462
463 if o1 == 0 && on_segment(p1, q1, p2) {
465 return true;
466 }
467 if o2 == 0 && on_segment(p1, q2, p2) {
468 return true;
469 }
470 if o3 == 0 && on_segment(q1, p1, q2) {
471 return true;
472 }
473 if o4 == 0 && on_segment(q1, p2, q2) {
474 return true;
475 }
476
477 false
478 }
479
480 fn edges_interfere_3d(p1: &[f64], p2: &[f64], q1: &[f64], q2: &[f64]) -> bool {
483 let eps = 1e-6; let u = [p2[0] - p1[0], p2[1] - p1[1], p2[2] - p1[2]];
488 let v = [q2[0] - q1[0], q2[1] - q1[1], q2[2] - q1[2]];
490 let w = [q1[0] - p1[0], q1[1] - p1[1], q1[2] - p1[2]];
492
493 let u_dot_u = u[0] * u[0] + u[1] * u[1] + u[2] * u[2];
494 let v_dot_v = v[0] * v[0] + v[1] * v[1] + v[2] * v[2];
495 let u_dot_v = u[0] * v[0] + u[1] * v[1] + u[2] * v[2];
496 let u_dot_w = u[0] * w[0] + u[1] * w[1] + u[2] * w[2];
497 let v_dot_w = v[0] * w[0] + v[1] * w[1] + v[2] * w[2];
498
499 let denom = u_dot_u * v_dot_v - u_dot_v * u_dot_v;
500
501 if denom.abs() < eps {
503 let cross_u_w = [
505 u[1] * w[2] - u[2] * w[1],
506 u[2] * w[0] - u[0] * w[2],
507 u[0] * w[1] - u[1] * w[0],
508 ];
509 let dist_sq = (cross_u_w[0] * cross_u_w[0]
510 + cross_u_w[1] * cross_u_w[1]
511 + cross_u_w[2] * cross_u_w[2])
512 / u_dot_u;
513 return dist_sq < eps * eps;
514 }
515
516 let s = (u_dot_v * v_dot_w - v_dot_v * u_dot_w) / denom;
518 let t = (u_dot_u * v_dot_w - u_dot_v * u_dot_w) / denom;
519
520 let s_clamped = s.clamp(0.0, 1.0);
522 let t_clamped = t.clamp(0.0, 1.0);
523
524 let closest_p = [
526 p1[0] + s_clamped * u[0],
527 p1[1] + s_clamped * u[1],
528 p1[2] + s_clamped * u[2],
529 ];
530 let closest_q = [
531 q1[0] + t_clamped * v[0],
532 q1[1] + t_clamped * v[1],
533 q1[2] + t_clamped * v[2],
534 ];
535
536 let dist_sq = (closest_p[0] - closest_q[0]) * (closest_p[0] - closest_q[0])
538 + (closest_p[1] - closest_q[1]) * (closest_p[1] - closest_q[1])
539 + (closest_p[2] - closest_q[2]) * (closest_p[2] - closest_q[2]);
540
541 dist_sq < eps * eps
542 }
543
544 fn find_triangles_with_edges(&self, edges: &[(usize, usize)]) -> Vec<usize> {
546 let mut triangles = HashSet::new();
547
548 for (i, simplex) in self.simplices.iter().enumerate() {
549 for &(e1, e2) in edges {
550 if self.triangle_contains_edge(simplex, e1, e2) {
551 triangles.insert(i);
552 }
553 }
554 }
555
556 triangles.into_iter().collect()
557 }
558
559 fn triangle_contains_edge(&self, triangle: &[usize], v1: usize, v2: usize) -> bool {
561 for i in 0..3 {
562 let j = (i + 1) % 3;
563 let t1 = triangle[i];
564 let t2 = triangle[j];
565 if (t1 == v1 && t2 == v2) || (t1 == v2 && t2 == v1) {
566 return true;
567 }
568 }
569 false
570 }
571
572 fn remove_triangles(&mut self, _triangleindices: &[usize]) {
574 let mut sorted_indices = _triangleindices.to_vec();
576 sorted_indices.sort_by(|a, b| b.cmp(a));
577
578 for &idx in &sorted_indices {
579 if idx < self.simplices.len() {
580 self.simplices.remove(idx);
581 self.neighbors.remove(idx);
582 }
583 }
584 }
585
586 fn retriangulate_with_constraint(
588 &mut self,
589 start: usize,
590 end: usize,
591 affected_triangles: &[usize],
592 ) -> SpatialResult<()> {
593 if affected_triangles.is_empty() {
594 return Ok(());
595 }
596
597 let cavity_vertices = self.extract_cavity_vertices(affected_triangles);
599
600 let boundary_edges = self.find_cavity_boundary(affected_triangles, start, end)?;
602
603 let new_triangles =
605 self.fan_triangulate_cavity(&cavity_vertices, &boundary_edges, start, end)?;
606
607 for triangle in new_triangles {
609 self.simplices.push(triangle);
610 }
611
612 self.compute_neighbors();
614
615 Ok(())
616 }
617
618 fn extract_cavity_vertices(&self, _affectedtriangles: &[usize]) -> Vec<usize> {
620 let mut vertices = HashSet::new();
621
622 for &triangle_idx in _affectedtriangles {
623 if triangle_idx < self.simplices.len() {
624 for &vertex in &self.simplices[triangle_idx] {
625 vertices.insert(vertex);
626 }
627 }
628 }
629
630 vertices.into_iter().collect()
631 }
632
633 fn find_cavity_boundary(
635 &self,
636 affected_triangles: &[usize],
637 start: usize,
638 end: usize,
639 ) -> SpatialResult<Vec<(usize, usize)>> {
640 let affected_set: HashSet<usize> = affected_triangles.iter().cloned().collect();
641 let mut boundary_edges = Vec::new();
642
643 for &triangle_idx in affected_triangles {
645 if triangle_idx >= self.simplices.len() {
646 continue;
647 }
648
649 let simplex = &self.simplices[triangle_idx];
650 if simplex.len() < 3 {
651 continue;
652 }
653
654 for i in 0..simplex.len() {
656 let v1 = simplex[i];
657 let v2 = simplex[(i + 1) % simplex.len()];
658
659 if (v1 == start && v2 == end) || (v1 == end && v2 == start) {
661 continue;
662 }
663
664 if self.is_boundary_edge(v1, v2, &affected_set, triangle_idx) {
666 boundary_edges.push((v1, v2));
667 }
668 }
669 }
670
671 Ok(boundary_edges)
672 }
673
674 fn is_boundary_edge(
676 &self,
677 v1: usize,
678 v2: usize,
679 affected_set: &HashSet<usize>,
680 current_triangle: usize,
681 ) -> bool {
682 for (tri_idx, simplex) in self.simplices.iter().enumerate() {
684 if tri_idx == current_triangle || affected_set.contains(&tri_idx) {
685 continue;
686 }
687
688 if self.triangle_contains_edge(simplex, v1, v2) {
690 return false; }
692 }
693
694 true }
696
697 fn fan_triangulate_cavity(
699 &self,
700 cavity_vertices: &[usize],
701 boundary_edges: &[(usize, usize)],
702 start: usize,
703 end: usize,
704 ) -> SpatialResult<Vec<Vec<usize>>> {
705 let mut new_triangles = Vec::new();
706
707 let mut interior_vertices = Vec::new();
709 for &vertex in cavity_vertices {
710 if vertex != start && vertex != end {
711 interior_vertices.push(vertex);
712 }
713 }
714
715 if !interior_vertices.is_empty() {
717 for i in 0..interior_vertices.len() {
719 for j in (i + 1)..interior_vertices.len() {
720 let v1 = interior_vertices[i];
721 let v2 = interior_vertices[j];
722
723 if self.is_valid_triangle_in_cavity(start, v1, v2, boundary_edges) {
725 new_triangles.push(vec![start, v1, v2]);
726 }
727
728 if self.is_valid_triangle_in_cavity(end, v1, v2, boundary_edges) {
729 new_triangles.push(vec![end, v1, v2]);
730 }
731 }
732 }
733 }
734
735 if new_triangles.is_empty() && !interior_vertices.is_empty() {
737 let v = interior_vertices[0];
738 new_triangles.push(vec![start, end, v]);
739 }
740
741 for &(v1, v2) in boundary_edges {
743 if v1 != start && v1 != end && v2 != start && v2 != end {
744 if self.points_form_valid_triangle(start, v1, v2) {
746 new_triangles.push(vec![start, v1, v2]);
747 }
748 if self.points_form_valid_triangle(end, v1, v2) {
749 new_triangles.push(vec![end, v1, v2]);
750 }
751 }
752 }
753
754 Ok(new_triangles)
755 }
756
757 fn points_form_valid_triangle(&self, v1: usize, v2: usize, v3: usize) -> bool {
759 if v1 >= self.npoints || v2 >= self.npoints || v3 >= self.npoints {
760 return false;
761 }
762
763 let p1 = self.points.row(v1);
764 let p2 = self.points.row(v2);
765 let p3 = self.points.row(v3);
766
767 let dx1 = p2[0] - p1[0];
769 let dy1 = p2[1] - p1[1];
770 let dx2 = p3[0] - p1[0];
771 let dy2 = p3[1] - p1[1];
772
773 let cross = dx1 * dy2 - dy1 * dx2;
774 cross.abs() > 1e-10 }
776
777 fn is_valid_triangle_in_cavity(
779 &self,
780 v1: usize,
781 v2: usize,
782 v3: usize,
783 _boundary_edges: &[(usize, usize)],
784 ) -> bool {
785 self.points_form_valid_triangle(v1, v2, v3)
787 }
788
789 fn compute_neighbors(&mut self) {
791 self.neighbors = Self::calculate_neighbors(&self.simplices, self.ndim + 1);
792 }
793
794 pub fn constraints(&self) -> &[(usize, usize)] {
800 &self.constraints
801 }
802
803 fn extract_simplices(qh: &Qh, ndim: usize) -> Vec<Vec<usize>> {
814 qh.simplices()
816 .filter(|f| !f.upper_delaunay())
817 .filter_map(|f| {
818 let vertices = match f.vertices() {
819 Some(v) => v,
820 None => return None,
821 };
822 let indices: Vec<usize> = vertices.iter().filter_map(|v| v.index(qh)).collect();
824
825 if indices.len() == ndim + 1 {
827 Some(indices)
828 } else {
829 None
830 }
831 })
832 .collect()
833 }
834
835 fn calculate_neighbors(simplices: &[Vec<usize>], n: usize) -> Vec<Vec<i64>> {
846 let nsimplex = simplices.len();
847 let mut neighbors = vec![vec![-1; n]; nsimplex];
848
849 let mut face_to_simplex: HashMap<Vec<usize>, Vec<(usize, usize)>> = HashMap::new();
852
853 for (i, simplex) in simplices.iter().enumerate() {
854 for j in 0..n {
855 let mut face: Vec<usize> = simplex
857 .iter()
858 .enumerate()
859 .filter(|&(k_, _)| k_ != j)
860 .map(|(_, &v)| v)
861 .collect();
862
863 face.sort();
865
866 face_to_simplex.entry(face).or_default().push((i, j));
868 }
869 }
870
871 for (_, simplex_info) in face_to_simplex.iter() {
873 if simplex_info.len() == 2 {
874 let (i1, j1) = simplex_info[0];
875 let (i2, j2) = simplex_info[1];
876
877 neighbors[i1][j1] = i2 as i64;
878 neighbors[i2][j2] = i1 as i64;
879 }
880 }
881
882 neighbors
883 }
884
885 pub fn npoints(&self) -> usize {
891 self.npoints
892 }
893
894 pub fn ndim(&self) -> usize {
900 self.ndim
901 }
902
903 pub fn points(&self) -> &Array2<f64> {
909 &self.points
910 }
911
912 pub fn simplices(&self) -> &[Vec<usize>] {
918 &self.simplices
919 }
920
921 pub fn neighbors(&self) -> &[Vec<i64>] {
927 &self.neighbors
928 }
929
930 pub fn find_simplex(&self, point: &[f64]) -> Option<usize> {
960 if point.len() != self.ndim {
961 return None;
962 }
963
964 if self.simplices.is_empty() {
965 return None;
966 }
967
968 for (i, simplex) in self.simplices.iter().enumerate() {
972 if self.point_in_simplex(point, simplex) {
973 return Some(i);
974 }
975 }
976
977 None
978 }
979
980 fn point_in_simplex(&self, point: &[f64], simplex: &[usize]) -> bool {
991 if self.ndim == 2 {
992 let a = self.points.row(simplex[0]).to_vec();
994 let b = self.points.row(simplex[1]).to_vec();
995 let c = self.points.row(simplex[2]).to_vec();
996
997 let v0x = b[0] - a[0];
998 let v0y = b[1] - a[1];
999 let v1x = c[0] - a[0];
1000 let v1y = c[1] - a[1];
1001 let v2x = point[0] - a[0];
1002 let v2y = point[1] - a[1];
1003
1004 let d00 = v0x * v0x + v0y * v0y;
1005 let d01 = v0x * v1x + v0y * v1y;
1006 let d11 = v1x * v1x + v1y * v1y;
1007 let d20 = v2x * v0x + v2y * v0y;
1008 let d21 = v2x * v1x + v2y * v1y;
1009
1010 let denom = d00 * d11 - d01 * d01;
1011 if denom.abs() < 1e-10 {
1012 return false; }
1014
1015 let v = (d11 * d20 - d01 * d21) / denom;
1016 let w = (d00 * d21 - d01 * d20) / denom;
1017 let u = 1.0 - v - w;
1018
1019 let eps = 1e-10;
1022 return u >= -eps && v >= -eps && w >= -eps;
1023 } else if self.ndim == 3 {
1024 let a = self.points.row(simplex[0]).to_vec();
1026 let b = self.points.row(simplex[1]).to_vec();
1027 let c = self.points.row(simplex[2]).to_vec();
1028 let d = self.points.row(simplex[3]).to_vec();
1029
1030 let mut bary = [0.0; 4];
1032
1033 let v0 = [b[0] - a[0], b[1] - a[1], b[2] - a[2]];
1035 let v1 = [c[0] - a[0], c[1] - a[1], c[2] - a[2]];
1036 let v2 = [d[0] - a[0], d[1] - a[1], d[2] - a[2]];
1037
1038 let vol = v0[0] * (v1[1] * v2[2] - v1[2] * v2[1])
1040 - v0[1] * (v1[0] * v2[2] - v1[2] * v2[0])
1041 + v0[2] * (v1[0] * v2[1] - v1[1] * v2[0]);
1042
1043 if vol.abs() < 1e-10 {
1044 return false; }
1046
1047 let _vp = [point[0] - a[0], point[1] - a[1], point[2] - a[2]];
1049
1050 let v3 = [c[0] - b[0], c[1] - b[1], c[2] - b[2]];
1051 let v4 = [d[0] - b[0], d[1] - b[1], d[2] - b[2]];
1052 let v5 = [point[0] - b[0], point[1] - b[1], point[2] - b[2]];
1053
1054 bary[0] = (v3[0] * (v4[1] * v5[2] - v4[2] * v5[1])
1055 - v3[1] * (v4[0] * v5[2] - v4[2] * v5[0])
1056 + v3[2] * (v4[0] * v5[1] - v4[1] * v5[0]))
1057 / vol;
1058
1059 let v3 = [b[0] - a[0], b[1] - a[1], b[2] - a[2]];
1060 let v4 = [d[0] - a[0], d[1] - a[1], d[2] - a[2]];
1061 let v5 = [point[0] - a[0], point[1] - a[1], point[2] - a[2]];
1062
1063 bary[1] = (v3[0] * (v4[1] * v5[2] - v4[2] * v5[1])
1064 - v3[1] * (v4[0] * v5[2] - v4[2] * v5[0])
1065 + v3[2] * (v4[0] * v5[1] - v4[1] * v5[0]))
1066 / vol;
1067
1068 let v3 = [b[0] - a[0], b[1] - a[1], b[2] - a[2]];
1069 let v4 = [c[0] - a[0], c[1] - a[1], c[2] - a[2]];
1070 let v5 = [point[0] - a[0], point[1] - a[1], point[2] - a[2]];
1071
1072 bary[2] = (v3[0] * (v4[1] * v5[2] - v4[2] * v5[1])
1073 - v3[1] * (v4[0] * v5[2] - v4[2] * v5[0])
1074 + v3[2] * (v4[0] * v5[1] - v4[1] * v5[0]))
1075 / vol;
1076
1077 bary[3] = 1.0 - bary[0] - bary[1] - bary[2];
1078
1079 let eps = 1e-10;
1081 return bary.iter().all(|&b| b >= -eps);
1082 }
1083
1084 false
1086 }
1087
1088 pub fn convex_hull(&self) -> Vec<usize> {
1114 let mut hull = HashSet::new();
1115
1116 for (i, neighbors) in self.neighbors.iter().enumerate() {
1118 for (j, &neighbor) in neighbors.iter().enumerate() {
1119 if neighbor == -1 {
1120 for k in 0..self.ndim + 1 {
1123 if k != j {
1124 hull.insert(self.simplices[i][k]);
1125 }
1126 }
1127 }
1128 }
1129 }
1130
1131 let mut hull_vec: Vec<usize> = hull.into_iter().collect();
1133 hull_vec.sort();
1134
1135 hull_vec
1136 }
1137}
1138
1139#[cfg(test)]
1140mod tests {
1141 use super::*;
1142 use scirs2_core::ndarray::arr2;
1143 use scirs2_core::random::Rng;
1144 #[test]
1147 fn test_delaunay_simple() {
1148 let points = arr2(&[[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]]);
1149
1150 let tri = Delaunay::new(&points).unwrap();
1151
1152 assert_eq!(tri.simplices().len(), 2);
1154
1155 for simplex in tri.simplices() {
1157 assert_eq!(simplex.len(), 3);
1158
1159 for &idx in simplex {
1161 assert!(idx < points.nrows());
1162 }
1163 }
1164
1165 let hull = tri.convex_hull();
1167 assert_eq!(hull.len(), 4); }
1169
1170 #[test]
1171 fn test_delaunay_with_interior_point() {
1172 let points = arr2(&[
1173 [0.0, 0.0],
1174 [1.0, 0.0],
1175 [0.0, 1.0],
1176 [0.5, 0.5], ]);
1178
1179 let tri = Delaunay::new(&points).unwrap();
1180
1181 assert_eq!(tri.simplices().len(), 3);
1183
1184 let hull = tri.convex_hull();
1186 assert_eq!(hull.len(), 3); assert!(!hull.contains(&3));
1190 }
1191
1192 #[test]
1193 fn test_delaunay_3d() {
1194 let points = arr2(&[
1195 [0.0, 0.0, 0.0],
1196 [1.0, 0.0, 0.0],
1197 [0.0, 1.0, 0.0],
1198 [0.0, 0.0, 1.0],
1199 [1.0, 1.0, 1.0],
1200 ]);
1201
1202 let tri = Delaunay::new(&points).unwrap();
1203
1204 for simplex in tri.simplices() {
1206 assert_eq!(simplex.len(), 4);
1207 }
1208 }
1209
1210 #[test]
1211 fn test_find_simplex() {
1212 let points = arr2(&[[0.0, 0.0], [1.0, 0.0], [0.0, 1.0]]);
1213
1214 let tri = Delaunay::new(&points).unwrap();
1215
1216 let inside_point = [0.3, 0.3];
1218 assert!(tri.find_simplex(&inside_point).is_some());
1219
1220 let outside_point = [1.5, 1.5];
1222 assert!(tri.find_simplex(&outside_point).is_none());
1223 }
1224
1225 #[test]
1226 fn test_random_points_2d() {
1227 let mut rng = rand::rng();
1229
1230 let n = 20;
1231 let mut points_data = Vec::with_capacity(n * 2);
1232
1233 for _ in 0..n {
1234 points_data.push(rng.gen_range(0.0..1.0));
1235 points_data.push(rng.gen_range(0.0..1.0));
1236 }
1237
1238 let points = Array2::from_shape_vec((n, 2), points_data).unwrap();
1239
1240 let tri = Delaunay::new(&points).unwrap();
1241
1242 assert_eq!(tri.ndim(), 2);
1244 assert_eq!(tri.npoints(), n);
1245
1246 for simplex in tri.simplices() {
1248 assert_eq!(simplex.len(), 3);
1249 for &idx in simplex {
1250 assert!(idx < n);
1251 }
1252 }
1253 }
1254
1255 #[test]
1256 fn test_constrained_delaunay_basic() {
1257 let points = arr2(&[
1258 [0.0, 0.0],
1259 [1.0, 0.0],
1260 [1.0, 1.0],
1261 [0.0, 1.0],
1262 [0.5, 0.5], ]);
1264
1265 let constraints = vec![(0, 1), (1, 2), (2, 3), (3, 0)];
1267
1268 let tri = Delaunay::new_constrained(&points, constraints.clone()).unwrap();
1269
1270 assert_eq!(tri.constraints().len(), 4);
1272 for &constraint in &constraints {
1273 assert!(tri.constraints().contains(&constraint));
1274 }
1275
1276 assert!(tri.simplices().len() >= 2); }
1279
1280 #[test]
1281 fn test_constrained_delaunay_invalid_constraints() {
1282 let points = arr2(&[[0.0, 0.0], [1.0, 0.0], [1.0, 1.0]]);
1283
1284 let invalid_constraints = vec![(0, 5)];
1286 let result = Delaunay::new_constrained(&points, invalid_constraints);
1287 assert!(result.is_err());
1288
1289 let self_constraint = vec![(0, 0)];
1291 let result = Delaunay::new_constrained(&points, self_constraint);
1292 assert!(result.is_err());
1293 }
1294
1295 #[test]
1296 fn test_constrained_delaunay_3d_error() {
1297 let points_3d = arr2(&[
1298 [0.0, 0.0, 0.0],
1299 [1.0, 0.0, 0.0],
1300 [0.0, 1.0, 0.0],
1301 [0.0, 0.0, 1.0],
1302 ]);
1303
1304 let constraints = vec![(0, 1)];
1305 let result = Delaunay::new_constrained(&points_3d, constraints);
1306 assert!(result.is_err()); }
1308
1309 #[test]
1310 fn test_edge_exists() {
1311 let points = arr2(&[[0.0, 0.0], [1.0, 0.0], [0.0, 1.0]]);
1312 let tri = Delaunay::new(&points).unwrap();
1313
1314 assert!(tri.edge_exists(0, 1) || tri.edge_exists(1, 0));
1316 assert!(tri.edge_exists(1, 2) || tri.edge_exists(2, 1));
1317 assert!(tri.edge_exists(0, 2) || tri.edge_exists(2, 0));
1318 }
1319
1320 #[test]
1321 fn test_segments_intersect() {
1322 let p1 = [0.0, 0.0];
1324 let p2 = [1.0, 1.0];
1325 let q1 = [0.0, 1.0];
1326 let q2 = [1.0, 0.0];
1327 assert!(Delaunay::segments_intersect(p1, p2, q1, q2));
1328
1329 let p1 = [0.0, 0.0];
1331 let p2 = [1.0, 0.0];
1332 let q1 = [0.0, 1.0];
1333 let q2 = [1.0, 1.0];
1334 assert!(!Delaunay::segments_intersect(p1, p2, q1, q2));
1335
1336 let p1 = [0.0, 0.0];
1338 let p2 = [2.0, 0.0];
1339 let q1 = [1.0, 0.0];
1340 let q2 = [3.0, 0.0];
1341 assert!(Delaunay::segments_intersect(p1, p2, q1, q2));
1342 }
1343}