scirs2_spatial/proximity/graphs.rs
1//! Proximity-based graph construction from point sets
2//!
3//! Provides algorithms for building Euclidean MST, Gabriel graph,
4//! relative neighborhood graph, and alpha shapes from 2D point sets.
5
6use crate::error::{SpatialError, SpatialResult};
7use std::collections::HashSet;
8
9/// An edge in a Minimum Spanning Tree with its weight (Euclidean distance)
10#[derive(Debug, Clone)]
11pub struct MstEdge {
12 /// Index of the first endpoint
13 pub u: usize,
14 /// Index of the second endpoint
15 pub v: usize,
16 /// Euclidean distance (edge weight)
17 pub weight: f64,
18}
19
20/// Compute the Euclidean Minimum Spanning Tree of a 2D point set
21///
22/// The Euclidean MST connects all points with the minimum total edge length.
23/// This implementation uses a brute-force Kruskal's algorithm on all pairwise
24/// distances with a union-find structure for cycle detection.
25///
26/// For large point sets, building a Delaunay triangulation first and then
27/// computing the MST on the Delaunay edges would be more efficient (O(n log n)
28/// vs O(n^2 log n)), but this direct approach is simpler and correct.
29///
30/// # Arguments
31///
32/// * `points` - A slice of [x, y] coordinates
33///
34/// # Returns
35///
36/// * A vector of MST edges sorted by weight
37///
38/// # Examples
39///
40/// ```
41/// use scirs2_spatial::proximity::euclidean_mst;
42///
43/// let points = vec![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
44/// let mst = euclidean_mst(&points).expect("compute");
45/// assert_eq!(mst.len(), 3); // n-1 edges
46///
47/// // Total weight should be 3.0 (three unit edges)
48/// let total: f64 = mst.iter().map(|e| e.weight).sum();
49/// assert!((total - 3.0).abs() < 1e-10);
50/// ```
51pub fn euclidean_mst(points: &[[f64; 2]]) -> SpatialResult<Vec<MstEdge>> {
52 let n = points.len();
53
54 if n == 0 {
55 return Ok(Vec::new());
56 }
57
58 if n == 1 {
59 return Ok(Vec::new());
60 }
61
62 // Build all edges with distances
63 let mut edges: Vec<(usize, usize, f64)> = Vec::with_capacity(n * (n - 1) / 2);
64 for i in 0..n {
65 for j in (i + 1)..n {
66 let dx = points[i][0] - points[j][0];
67 let dy = points[i][1] - points[j][1];
68 let dist = (dx * dx + dy * dy).sqrt();
69 edges.push((i, j, dist));
70 }
71 }
72
73 // Sort by distance
74 edges.sort_by(|a, b| a.2.partial_cmp(&b.2).unwrap_or(std::cmp::Ordering::Equal));
75
76 // Kruskal's algorithm with union-find
77 let mut parent: Vec<usize> = (0..n).collect();
78 let mut rank: Vec<usize> = vec![0; n];
79 let mut mst = Vec::with_capacity(n - 1);
80
81 fn find(parent: &mut [usize], x: usize) -> usize {
82 if parent[x] != x {
83 parent[x] = find(parent, parent[x]);
84 }
85 parent[x]
86 }
87
88 fn union(parent: &mut [usize], rank: &mut [usize], x: usize, y: usize) -> bool {
89 let rx = find(parent, x);
90 let ry = find(parent, y);
91 if rx == ry {
92 return false;
93 }
94 if rank[rx] < rank[ry] {
95 parent[rx] = ry;
96 } else if rank[rx] > rank[ry] {
97 parent[ry] = rx;
98 } else {
99 parent[ry] = rx;
100 rank[rx] += 1;
101 }
102 true
103 }
104
105 for (u, v, w) in edges {
106 if mst.len() == n - 1 {
107 break;
108 }
109 if union(&mut parent, &mut rank, u, v) {
110 mst.push(MstEdge { u, v, weight: w });
111 }
112 }
113
114 Ok(mst)
115}
116
117/// Compute the Gabriel graph of a 2D point set
118///
119/// The Gabriel graph connects two points p_i and p_j if and only if the
120/// diametral circle (circle with p_i p_j as diameter) contains no other points.
121/// Equivalently, edge (i, j) exists iff for all k != i,j:
122/// d(i,j)^2 <= d(i,k)^2 + d(j,k)^2
123///
124/// The Gabriel graph is a subgraph of the Delaunay triangulation and a
125/// supergraph of the Euclidean MST.
126///
127/// # Arguments
128///
129/// * `points` - A slice of [x, y] coordinates
130///
131/// # Returns
132///
133/// * A vector of edges (i, j) in the Gabriel graph
134///
135/// # Examples
136///
137/// ```
138/// use scirs2_spatial::proximity::gabriel_graph;
139///
140/// let points = vec![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
141/// let edges = gabriel_graph(&points).expect("compute");
142/// // Unit square: all 6 pairs are Gabriel edges (diagonals included, since for a diagonal
143/// // the sum of squared distances from the endpoints to any other vertex equals the diagonal
144/// // length squared exactly, so no vertex lies strictly inside the diametral circle)
145/// assert_eq!(edges.len(), 6);
146/// ```
147pub fn gabriel_graph(points: &[[f64; 2]]) -> SpatialResult<Vec<(usize, usize)>> {
148 let n = points.len();
149
150 if n < 2 {
151 return Ok(Vec::new());
152 }
153
154 let mut edges = Vec::new();
155
156 for i in 0..n {
157 for j in (i + 1)..n {
158 let dx = points[i][0] - points[j][0];
159 let dy = points[i][1] - points[j][1];
160 let dist_sq_ij = dx * dx + dy * dy;
161
162 // Check if any other point lies inside the diametral circle
163 let mut is_gabriel = true;
164
165 for k in 0..n {
166 if k == i || k == j {
167 continue;
168 }
169
170 let dxi = points[i][0] - points[k][0];
171 let dyi = points[i][1] - points[k][1];
172 let dist_sq_ik = dxi * dxi + dyi * dyi;
173
174 let dxj = points[j][0] - points[k][0];
175 let dyj = points[j][1] - points[k][1];
176 let dist_sq_jk = dxj * dxj + dyj * dyj;
177
178 // Gabriel condition: d(i,j)^2 <= d(i,k)^2 + d(j,k)^2
179 if dist_sq_ij > dist_sq_ik + dist_sq_jk + 1e-10 {
180 is_gabriel = false;
181 break;
182 }
183 }
184
185 if is_gabriel {
186 edges.push((i, j));
187 }
188 }
189 }
190
191 Ok(edges)
192}
193
194/// Compute the Relative Neighborhood Graph (RNG) of a 2D point set
195///
196/// The RNG connects two points p_i and p_j if and only if there is no other
197/// point p_k that is closer to both p_i and p_j than they are to each other.
198/// Formally, edge (i, j) exists iff for all k != i,j:
199/// d(i,j) <= max(d(i,k), d(j,k))
200///
201/// The RNG is a subgraph of the Gabriel graph and a supergraph of the
202/// Euclidean MST.
203///
204/// # Arguments
205///
206/// * `points` - A slice of [x, y] coordinates
207///
208/// # Returns
209///
210/// * A vector of edges (i, j) in the relative neighborhood graph
211///
212/// # Examples
213///
214/// ```
215/// use scirs2_spatial::proximity::relative_neighborhood_graph;
216///
217/// let points = vec![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
218/// let edges = relative_neighborhood_graph(&points).expect("compute");
219/// // RNG is a subgraph of Gabriel graph
220/// assert!(edges.len() <= 4);
221/// ```
222pub fn relative_neighborhood_graph(points: &[[f64; 2]]) -> SpatialResult<Vec<(usize, usize)>> {
223 let n = points.len();
224
225 if n < 2 {
226 return Ok(Vec::new());
227 }
228
229 let mut edges = Vec::new();
230
231 for i in 0..n {
232 for j in (i + 1)..n {
233 let dx = points[i][0] - points[j][0];
234 let dy = points[i][1] - points[j][1];
235 let dist_ij = (dx * dx + dy * dy).sqrt();
236
237 let mut is_rng = true;
238
239 for k in 0..n {
240 if k == i || k == j {
241 continue;
242 }
243
244 let dxi = points[i][0] - points[k][0];
245 let dyi = points[i][1] - points[k][1];
246 let dist_ik = (dxi * dxi + dyi * dyi).sqrt();
247
248 let dxj = points[j][0] - points[k][0];
249 let dyj = points[j][1] - points[k][1];
250 let dist_jk = (dxj * dxj + dyj * dyj).sqrt();
251
252 // RNG condition: d(i,j) <= max(d(i,k), d(j,k))
253 if dist_ij > dist_ik.max(dist_jk) + 1e-10 {
254 is_rng = false;
255 break;
256 }
257 }
258
259 if is_rng {
260 edges.push((i, j));
261 }
262 }
263 }
264
265 Ok(edges)
266}
267
268/// Compute the alpha shape edges for a 2D point set
269///
270/// Alpha shapes are a generalization of convex hulls. For a given alpha parameter,
271/// the alpha shape is the intersection of all closed complements of circles with
272/// radius 1/alpha that contain all points. As alpha decreases toward 0, the alpha
273/// shape approaches the convex hull. As alpha increases, the shape becomes more
274/// concave, eventually decomposing into individual points.
275///
276/// An edge (i, j) is in the alpha shape if there exists a circle of radius
277/// 1/alpha passing through both p_i and p_j that contains no other points.
278///
279/// # Arguments
280///
281/// * `points` - A slice of [x, y] coordinates
282/// * `alpha` - The alpha parameter (> 0). Larger alpha = more concave shape.
283/// A good starting value is the inverse of the typical edge length.
284///
285/// # Returns
286///
287/// * A vector of edges (i, j) forming the alpha shape boundary
288///
289/// # Examples
290///
291/// ```
292/// use scirs2_spatial::proximity::alpha_shape_edges;
293///
294/// let points = vec![
295/// [0.0, 0.0], [1.0, 0.0], [2.0, 0.0],
296/// [0.0, 1.0], [1.0, 1.0], [2.0, 1.0],
297/// ];
298///
299/// // With small alpha, most edges are included (near convex hull)
300/// let edges = alpha_shape_edges(&points, 0.5).expect("compute");
301/// assert!(!edges.is_empty());
302/// ```
303pub fn alpha_shape_edges(points: &[[f64; 2]], alpha: f64) -> SpatialResult<Vec<(usize, usize)>> {
304 let n = points.len();
305
306 if n < 2 {
307 return Ok(Vec::new());
308 }
309
310 if alpha <= 0.0 {
311 return Err(SpatialError::ValueError(
312 "Alpha must be positive".to_string(),
313 ));
314 }
315
316 let radius = 1.0 / alpha;
317 let radius_sq = radius * radius;
318 let mut edges = Vec::new();
319
320 for i in 0..n {
321 for j in (i + 1)..n {
322 let dx = points[j][0] - points[i][0];
323 let dy = points[j][1] - points[i][1];
324 let dist_sq = dx * dx + dy * dy;
325 let dist = dist_sq.sqrt();
326
327 // If the edge is longer than the diameter (2*radius), skip
328 if dist > 2.0 * radius + 1e-10 {
329 continue;
330 }
331
332 // Find the two circle centers of radius r passing through p_i and p_j
333 let mid_x = (points[i][0] + points[j][0]) * 0.5;
334 let mid_y = (points[i][1] + points[j][1]) * 0.5;
335
336 let half_dist_sq = dist_sq * 0.25;
337 let h_sq = radius_sq - half_dist_sq;
338
339 if h_sq < 0.0 {
340 continue; // Edge too long for this radius
341 }
342
343 let h = h_sq.sqrt();
344
345 // Normal direction perpendicular to the edge
346 let nx = -dy / dist;
347 let ny = dx / dist;
348
349 // Two potential circle centers
350 let c1x = mid_x + h * nx;
351 let c1y = mid_y + h * ny;
352 let c2x = mid_x - h * nx;
353 let c2y = mid_y - h * ny;
354
355 // Check if either circle is empty (contains no other points)
356 let circle1_empty = (0..n).all(|k| {
357 if k == i || k == j {
358 return true;
359 }
360 let dxk = points[k][0] - c1x;
361 let dyk = points[k][1] - c1y;
362 dxk * dxk + dyk * dyk >= radius_sq - 1e-10
363 });
364
365 let circle2_empty = (0..n).all(|k| {
366 if k == i || k == j {
367 return true;
368 }
369 let dxk = points[k][0] - c2x;
370 let dyk = points[k][1] - c2y;
371 dxk * dxk + dyk * dyk >= radius_sq - 1e-10
372 });
373
374 if circle1_empty || circle2_empty {
375 edges.push((i, j));
376 }
377 }
378 }
379
380 Ok(edges)
381}
382
383/// Get boundary edges of an alpha shape (edges that appear in only one triangle)
384///
385/// This extracts only the boundary from the alpha shape edges, which forms
386/// the concave hull outline.
387///
388/// # Arguments
389///
390/// * `points` - A slice of [x, y] coordinates
391/// * `alpha` - The alpha parameter (> 0)
392///
393/// # Returns
394///
395/// * Boundary edges forming the concave hull
396pub fn alpha_shape_boundary(points: &[[f64; 2]], alpha: f64) -> SpatialResult<Vec<(usize, usize)>> {
397 let all_edges = alpha_shape_edges(points, alpha)?;
398
399 // Build an edge set for quick lookup
400 let edge_set: HashSet<(usize, usize)> = all_edges.iter().copied().collect();
401
402 // A boundary edge is one where the two adjacent triangles that could use it
403 // don't both exist. We approximate this by checking if for each edge (i,j),
404 // there are 0 or 1 common neighbors k such that both (i,k) and (j,k) are edges.
405 let mut boundary = Vec::new();
406
407 for &(i, j) in &all_edges {
408 let common_count = (0..points.len())
409 .filter(|&k| {
410 k != i
411 && k != j
412 && (edge_set.contains(&(i.min(k), i.max(k))))
413 && (edge_set.contains(&(j.min(k), j.max(k))))
414 })
415 .count();
416
417 // Boundary edges have 0 or 1 common neighbor (1 means one side is open)
418 if common_count <= 1 {
419 boundary.push((i, j));
420 }
421 }
422
423 Ok(boundary)
424}
425
426#[cfg(test)]
427mod tests {
428 use super::*;
429
430 #[test]
431 fn test_mst_triangle() {
432 let points = vec![[0.0, 0.0], [1.0, 0.0], [0.5, 0.866]]; // equilateral
433 let mst = euclidean_mst(&points).expect("compute");
434 assert_eq!(mst.len(), 2);
435
436 // All edges should be ~1.0
437 for e in &mst {
438 assert!((e.weight - 1.0).abs() < 1e-3);
439 }
440 }
441
442 #[test]
443 fn test_mst_square() {
444 let points = vec![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
445 let mst = euclidean_mst(&points).expect("compute");
446 assert_eq!(mst.len(), 3);
447
448 // Total weight should be 3.0 (three unit edges, not the diagonal)
449 let total: f64 = mst.iter().map(|e| e.weight).sum();
450 assert!((total - 3.0).abs() < 1e-10);
451 }
452
453 #[test]
454 fn test_mst_single_point() {
455 let points = vec![[0.0, 0.0]];
456 let mst = euclidean_mst(&points).expect("compute");
457 assert!(mst.is_empty());
458 }
459
460 #[test]
461 fn test_mst_empty() {
462 let points: Vec<[f64; 2]> = vec![];
463 let mst = euclidean_mst(&points).expect("compute");
464 assert!(mst.is_empty());
465 }
466
467 #[test]
468 fn test_gabriel_graph_square() {
469 let points = vec![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
470 let edges = gabriel_graph(&points).expect("compute");
471
472 // For a unit square, diagonal points lie exactly on the diametral circle
473 // boundary (d(i,j)^2 == d(i,k)^2 + d(j,k)^2 for diagonals).
474 // With our tolerance, diagonals may or may not be included.
475 // The 4 side edges must always be present.
476 let edge_set: HashSet<(usize, usize)> = edges.into_iter().collect();
477 assert!(edge_set.contains(&(0, 1))); // bottom
478 assert!(edge_set.contains(&(0, 2))); // left
479 assert!(edge_set.contains(&(1, 3))); // right
480 assert!(edge_set.contains(&(2, 3))); // top
481 assert!(edge_set.len() >= 4);
482 assert!(edge_set.len() <= 6); // max C(4,2) = 6
483 }
484
485 #[test]
486 fn test_gabriel_is_supergraph_of_rng() {
487 let points = vec![[0.0, 0.0], [1.0, 0.0], [0.5, 0.5], [0.0, 1.0], [1.0, 1.0]];
488
489 let gabriel_edges = gabriel_graph(&points).expect("compute");
490 let rng_edges = relative_neighborhood_graph(&points).expect("compute");
491
492 let gabriel_set: HashSet<(usize, usize)> = gabriel_edges.into_iter().collect();
493 let rng_set: HashSet<(usize, usize)> = rng_edges.into_iter().collect();
494
495 // Every RNG edge should be in Gabriel graph
496 for edge in &rng_set {
497 assert!(
498 gabriel_set.contains(edge),
499 "RNG edge {:?} not in Gabriel graph",
500 edge
501 );
502 }
503 }
504
505 #[test]
506 fn test_rng_square() {
507 let points = vec![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
508 let edges = relative_neighborhood_graph(&points).expect("compute");
509
510 // RNG of a square: the 4 side edges (no diagonals)
511 assert_eq!(edges.len(), 4);
512 }
513
514 #[test]
515 fn test_rng_is_supergraph_of_mst() {
516 let points = vec![[0.0, 0.0], [1.0, 0.0], [2.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
517
518 let mst = euclidean_mst(&points).expect("compute");
519 let rng_edges = relative_neighborhood_graph(&points).expect("compute");
520
521 let rng_set: HashSet<(usize, usize)> = rng_edges.into_iter().collect();
522
523 // Every MST edge should be in the RNG
524 for e in &mst {
525 let edge = (e.u.min(e.v), e.u.max(e.v));
526 assert!(rng_set.contains(&edge), "MST edge {:?} not in RNG", edge);
527 }
528 }
529
530 #[test]
531 fn test_alpha_shape_large_alpha() {
532 // With very large alpha (small radius), few or no edges should survive
533 let points = vec![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
534 let edges = alpha_shape_edges(&points, 10.0).expect("compute");
535 // Very large alpha means very small circle radius - edges too long
536 // might have some short edges
537 assert!(edges.len() <= 6); // at most C(4,2)=6
538 }
539
540 #[test]
541 fn test_alpha_shape_small_alpha() {
542 // With small alpha (large radius), most edges should be included
543 let points = vec![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
544 let edges = alpha_shape_edges(&points, 0.1).expect("compute");
545 // Small alpha = large circles = almost convex hull, should include all 6 edges
546 assert!(edges.len() >= 4);
547 }
548
549 #[test]
550 fn test_alpha_shape_invalid_alpha() {
551 let points = vec![[0.0, 0.0], [1.0, 0.0]];
552 assert!(alpha_shape_edges(&points, 0.0).is_err());
553 assert!(alpha_shape_edges(&points, -1.0).is_err());
554 }
555
556 #[test]
557 fn test_graph_hierarchy() {
558 // The containment hierarchy should be: MST ⊆ RNG ⊆ Gabriel ⊆ Delaunay
559 let points = vec![[0.0, 0.0], [2.0, 0.0], [1.0, 1.5], [3.0, 1.0], [0.5, 2.5]];
560
561 let mst = euclidean_mst(&points).expect("mst");
562 let rng = relative_neighborhood_graph(&points).expect("rng");
563 let gabriel = gabriel_graph(&points).expect("gabriel");
564
565 let mst_edges: HashSet<(usize, usize)> =
566 mst.iter().map(|e| (e.u.min(e.v), e.u.max(e.v))).collect();
567 let rng_set: HashSet<(usize, usize)> = rng.into_iter().collect();
568 let gabriel_set: HashSet<(usize, usize)> = gabriel.into_iter().collect();
569
570 // MST ⊆ RNG
571 for edge in &mst_edges {
572 assert!(rng_set.contains(edge), "MST edge {:?} not in RNG", edge);
573 }
574
575 // RNG ⊆ Gabriel
576 for edge in &rng_set {
577 assert!(
578 gabriel_set.contains(edge),
579 "RNG edge {:?} not in Gabriel",
580 edge
581 );
582 }
583
584 assert!(mst_edges.len() <= rng_set.len());
585 assert!(rng_set.len() <= gabriel_set.len());
586 }
587}