scirs2-spatial 0.4.2

Spatial algorithms module for SciRS2 (scirs2-spatial)
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
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
//! KD-Tree for efficient nearest neighbor searches
//!
//! This module provides a KD-Tree implementation for efficient
//! nearest neighbor and range searches in multidimensional spaces.
//!
//! The KD-Tree (k-dimensional tree) is a space-partitioning data structure
//! that organizes points in a k-dimensional space. It enables efficient range searches
//! and nearest neighbor searches.
//!
//! # Features
//!
//! * Fast nearest neighbor queries with customizable `k`
//! * Range queries with distance threshold
//! * Support for different distance metrics (Euclidean, Manhattan, Chebyshev, etc.)
//! * Parallel query processing (when using the `parallel` feature)
//! * Customizable leaf size for performance tuning
//!
//! # Examples
//!
//! ```
//! use scirs2_spatial::KDTree;
//! use scirs2_core::ndarray::array;
//!
//! // Create a KD-Tree with points in 2D space
//! let points = array![[1.0, 2.0], [3.0, 4.0], [5.0, 6.0], [7.0, 8.0]];
//! let kdtree = KDTree::new(&points).expect("Operation failed");
//!
//! // Find the nearest neighbor to [4.0, 5.0]
//! let (idx, dist) = kdtree.query(&[4.0, 5.0], 1).expect("Operation failed");
//! assert_eq!(idx.len(), 1); // Should return exactly one neighbor
//!
//! // Find all points within radius 3.0 of [4.0, 5.0]
//! let (indices, distances) = kdtree.query_radius(&[4.0, 5.0], 3.0).expect("Operation failed");
//! ```
//!
//! # Advanced Usage
//!
//! Using custom distance metrics:
//!
//! ```
//! use scirs2_spatial::KDTree;
//! use scirs2_spatial::distance::ManhattanDistance;
//! use scirs2_core::ndarray::array;
//!
//! // Create a KD-Tree with Manhattan distance metric
//! let points = array![[1.0, 2.0], [3.0, 4.0], [5.0, 6.0], [7.0, 8.0]];
//! let metric = ManhattanDistance::new();
//! let kdtree = KDTree::with_metric(&points, metric).expect("Operation failed");
//!
//! // Find the nearest neighbor to [4.0, 5.0] using Manhattan distance
//! let (idx, dist) = kdtree.query(&[4.0, 5.0], 1).expect("Operation failed");
//! ```
//!
//! Using custom leaf size for performance tuning:
//!
//! ```
//! use scirs2_spatial::KDTree;
//! use scirs2_core::ndarray::array;
//!
//! // Create a KD-Tree with custom leaf size
//! let points = array![[1.0, 2.0], [3.0, 4.0], [5.0, 6.0], [7.0, 8.0],
//!                      [9.0, 10.0], [11.0, 12.0], [13.0, 14.0], [15.0, 16.0]];
//! let leafsize = 2; // Default is 16
//! let kdtree = KDTree::with_leaf_size(&points, leafsize).expect("Operation failed");
//! ```

use crate::distance::{Distance, EuclideanDistance};
use crate::error::{SpatialError, SpatialResult};
use crate::safe_conversions::*;
use scirs2_core::ndarray::Array2;
use scirs2_core::numeric::Float;
use std::cmp::Ordering;

// Rayon parallel processing currently not used in this module
/// A rectangle representing a hyperrectangle in k-dimensional space
///
/// Used for efficient nearest-neighbor and range queries in KD-trees.
#[derive(Clone, Debug)]
pub struct Rectangle<T: Float> {
    /// Minimum coordinates for each dimension
    mins: Vec<T>,
    /// Maximum coordinates for each dimension
    maxes: Vec<T>,
}

impl<T: Float> Rectangle<T> {
    /// Create a new hyperrectangle
    ///
    /// # Arguments
    ///
    /// * `mins` - Minimum coordinates for each dimension
    /// * `maxes` - Maximum coordinates for each dimension
    ///
    /// # Returns
    ///
    /// * A new rectangle
    ///
    /// # Panics
    ///
    /// * If mins and maxes have different lengths
    /// * If any min value is greater than the corresponding max value
    pub fn new(mins: Vec<T>, maxes: Vec<T>) -> Self {
        assert_eq!(
            mins.len(),
            maxes.len(),
            "mins and maxes must have the same length"
        );

        for i in 0..mins.len() {
            assert!(
                mins[i] <= maxes[i],
                "min value must be less than or equal to max value"
            );
        }

        Rectangle { mins, maxes }
    }

    /// Get the minimum coordinates of the rectangle
    ///
    /// # Returns
    ///
    /// * A slice containing the minimum coordinate values for each dimension
    ///
    /// # Examples
    ///
    /// ```
    /// use scirs2_spatial::kdtree::Rectangle;
    ///
    /// let rect = Rectangle::new(vec![0.0, 0.0], vec![1.0, 1.0]);
    /// let mins = rect.mins();
    /// assert_eq!(mins, &[0.0, 0.0]);
    /// ```
    pub fn mins(&self) -> &[T] {
        &self.mins
    }

    /// Get the maximum coordinates of the rectangle
    ///
    /// # Returns
    ///
    /// * A slice containing the maximum coordinate values for each dimension
    ///
    /// # Examples
    ///
    /// ```
    /// use scirs2_spatial::kdtree::Rectangle;
    ///
    /// let rect = Rectangle::new(vec![0.0, 0.0], vec![1.0, 1.0]);
    /// let maxes = rect.maxes();
    /// assert_eq!(maxes, &[1.0, 1.0]);
    /// ```
    pub fn maxes(&self) -> &[T] {
        &self.maxes
    }

    /// Split the rectangle along a given dimension at a given value
    ///
    /// # Arguments
    ///
    /// * `dim` - The dimension to split on
    /// * `value` - The value to split at
    ///
    /// # Returns
    ///
    /// * A tuple of (left, right) rectangles
    pub fn split(&self, dim: usize, value: T) -> (Self, Self) {
        let mut left_maxes = self.maxes.clone();
        left_maxes[dim] = value;

        let mut right_mins = self.mins.clone();
        right_mins[dim] = value;

        let left = Rectangle::new(self.mins.clone(), left_maxes);
        let right = Rectangle::new(right_mins, self.maxes.clone());

        (left, right)
    }

    /// Check if the rectangle contains a point
    ///
    /// # Arguments
    ///
    /// * `point` - The point to check
    ///
    /// # Returns
    ///
    /// * true if the rectangle contains the point, false otherwise
    pub fn contains(&self, point: &[T]) -> bool {
        assert_eq!(
            point.len(),
            self.mins.len(),
            "point must have the same dimension as the rectangle"
        );

        for (i, &p) in point.iter().enumerate() {
            if p < self.mins[i] || p > self.maxes[i] {
                return false;
            }
        }

        true
    }

    /// Calculate the minimum distance from a point to the rectangle
    ///
    /// # Arguments
    ///
    /// * `point` - The point to calculate distance to
    /// * `metric` - The distance metric to use
    ///
    /// # Returns
    ///
    /// * The minimum distance from the point to any point in the rectangle
    pub fn min_distance<D: Distance<T>>(&self, point: &[T], metric: &D) -> T {
        metric.min_distance_point_rectangle(point, &self.mins, &self.maxes)
    }
}

/// A node in the KD-Tree
#[derive(Debug, Clone)]
struct KDNode<T: Float> {
    /// Index of the point in the original data array
    idx: usize,
    /// The value of the point along the splitting dimension
    value: T,
    /// The dimension used for splitting
    axis: usize,
    /// Left child node (values < median along splitting axis)
    left: Option<usize>,
    /// Right child node (values >= median along splitting axis)
    right: Option<usize>,
}

/// A KD-Tree for efficient nearest neighbor searches
///
/// # Type Parameters
///
/// * `T` - The floating point type for coordinates
/// * `D` - The distance metric type
#[derive(Debug, Clone)]
pub struct KDTree<T: Float + Send + Sync + 'static, D: Distance<T> + 'static> {
    /// The points stored in the KD-Tree
    points: Array2<T>,
    /// The nodes of the KD-Tree
    nodes: Vec<KDNode<T>>,
    /// The dimensionality of the points
    ndim: usize,
    /// The root node index
    root: Option<usize>,
    /// The distance metric
    metric: D,
    /// The leaf size (maximum number of points in a leaf node)
    leafsize: usize,
    /// Minimum bounding rectangle of the entire dataset
    bounds: Rectangle<T>,
}

impl<T: Float + Send + Sync + 'static> KDTree<T, EuclideanDistance<T>> {
    /// Create a new KD-Tree with default Euclidean distance metric
    ///
    /// # Arguments
    ///
    /// * `points` - Array of points, each row is a point
    ///
    /// # Returns
    ///
    /// * A new KD-Tree
    pub fn new(points: &Array2<T>) -> SpatialResult<Self> {
        let metric = EuclideanDistance::new();
        Self::with_metric(points, metric)
    }

    /// Create a new KD-Tree with custom leaf size (using Euclidean distance)
    ///
    /// # Arguments
    ///
    /// * `points` - Array of points, each row is a point
    /// * `leafsize` - The maximum number of points in a leaf node
    ///
    /// # Returns
    ///
    /// * A new KD-Tree
    pub fn with_leaf_size(points: &Array2<T>, leafsize: usize) -> SpatialResult<Self> {
        let metric = EuclideanDistance::new();
        Self::with_options(points, metric, leafsize)
    }
}

impl<T: Float + Send + Sync + 'static, D: Distance<T> + 'static> KDTree<T, D> {
    /// Create a new KD-Tree with custom distance metric
    ///
    /// # Arguments
    ///
    /// * `points` - Array of points, each row is a point
    /// * `metric` - The distance metric to use
    ///
    /// # Returns
    ///
    /// * A new KD-Tree
    pub fn with_metric(points: &Array2<T>, metric: D) -> SpatialResult<Self> {
        Self::with_options(points, metric, 16) // Default leaf size is 16
    }

    /// Create a new KD-Tree with custom distance metric and leaf size
    ///
    /// # Arguments
    ///
    /// * `points` - Array of points, each row is a point
    /// * `metric` - The distance metric to use
    /// * `leafsize` - The maximum number of points in a leaf node
    ///
    /// # Returns
    ///
    /// * A new KD-Tree
    pub fn with_options(points: &Array2<T>, metric: D, leafsize: usize) -> SpatialResult<Self> {
        let n = points.nrows();
        let ndim = points.ncols();

        if n == 0 {
            return Err(SpatialError::ValueError("Empty point set".to_string()));
        }

        if leafsize == 0 {
            return Err(SpatialError::ValueError(
                "Leaf _size must be greater than 0".to_string(),
            ));
        }

        // Calculate the bounds of the dataset
        let mut mins = vec![T::max_value(); ndim];
        let mut maxes = vec![T::min_value(); ndim];

        for i in 0..n {
            for j in 0..ndim {
                let val = points[[i, j]];
                if val < mins[j] {
                    mins[j] = val;
                }
                if val > maxes[j] {
                    maxes[j] = val;
                }
            }
        }

        let bounds = Rectangle::new(mins, maxes);

        let mut tree = KDTree {
            points: points.clone(),
            nodes: Vec::with_capacity(n),
            ndim,
            root: None,
            metric,
            leafsize,
            bounds,
        };

        // Create indices for the points
        let mut indices: Vec<usize> = (0..n).collect();

        // Build the tree recursively
        if n > 0 {
            let root = tree.build_tree(&mut indices, 0, 0, n)?;
            tree.root = Some(root);
        }

        Ok(tree)
    }

    /// Build the KD-Tree recursively
    ///
    /// # Arguments
    ///
    /// * `indices` - Indices of the points to consider
    /// * `depth` - Current depth in the tree
    /// * `start` - Start index in the indices array
    /// * `end` - End index in the indices array
    ///
    /// # Returns
    ///
    /// * Index of the root node of the subtree
    fn build_tree(
        &mut self,
        indices: &mut [usize],
        depth: usize,
        start: usize,
        end: usize,
    ) -> SpatialResult<usize> {
        let n = end - start;

        if n == 0 {
            return Err(SpatialError::ValueError(
                "Empty point set in build_tree".to_string(),
            ));
        }

        // Choose axis based on depth (cycle through axes)
        let axis = depth % self.ndim;

        // If we have only one point, create a leaf node
        let node_idx;
        if n == 1 {
            let idx = indices[start];
            let value = self.points[[idx, axis]];

            node_idx = self.nodes.len();
            self.nodes.push(KDNode {
                idx,
                value,
                axis,
                left: None,
                right: None,
            });

            return Ok(node_idx);
        }

        // Sort indices based on the axis
        indices[start..end].sort_by(|&i, &j| {
            let a = self.points[[i, axis]];
            let b = self.points[[j, axis]];
            a.partial_cmp(&b).unwrap_or(Ordering::Equal)
        });

        // Get the median index
        let mid = start + n / 2;
        let idx = indices[mid];
        let value = self.points[[idx, axis]];

        // Create node
        node_idx = self.nodes.len();
        self.nodes.push(KDNode {
            idx,
            value,
            axis,
            left: None,
            right: None,
        });

        // Build left and right subtrees
        if mid > start {
            let left_idx = self.build_tree(indices, depth + 1, start, mid)?;
            self.nodes[node_idx].left = Some(left_idx);
        }

        if mid + 1 < end {
            let right_idx = self.build_tree(indices, depth + 1, mid + 1, end)?;
            self.nodes[node_idx].right = Some(right_idx);
        }

        Ok(node_idx)
    }

    /// Find the k nearest neighbors to a query point
    ///
    /// # Arguments
    ///
    /// * `point` - Query point
    /// * `k` - Number of nearest neighbors to find
    ///
    /// # Returns
    ///
    /// * (indices, distances) of the k nearest neighbors
    ///
    /// # Examples
    ///
    /// ```
    /// use scirs2_spatial::KDTree;
    /// use scirs2_core::ndarray::array;
    ///
    /// // Create points for the KDTree - use the exact same points from test_kdtree_with_custom_leaf_size
    /// let points = array![[2.0, 3.0], [5.0, 4.0], [9.0, 6.0], [4.0, 7.0], [8.0, 1.0], [7.0, 2.0]];
    /// let kdtree = KDTree::new(&points).expect("Operation failed");
    ///
    /// // Find the 2 nearest neighbors to [0.5, 0.5]
    /// let (indices, distances) = kdtree.query(&[0.5, 0.5], 2).expect("Operation failed");
    /// assert_eq!(indices.len(), 2);
    /// assert_eq!(distances.len(), 2);
    /// ```
    pub fn query(&self, point: &[T], k: usize) -> SpatialResult<(Vec<usize>, Vec<T>)> {
        if point.len() != self.ndim {
            return Err(SpatialError::DimensionError(format!(
                "Query point dimension ({}) does not match tree dimension ({})",
                point.len(),
                self.ndim
            )));
        }

        if k == 0 {
            return Ok((vec![], vec![]));
        }

        if self.points.nrows() == 0 {
            return Ok((vec![], vec![]));
        }

        // Initialize priority queue for k nearest neighbors
        // We use a max-heap so we can efficiently replace the furthest point when we find a closer one
        let mut neighbors: Vec<(T, usize)> = Vec::with_capacity(k + 1);

        // Keep track of the maximum distance in the heap, for early termination
        let mut max_dist = T::infinity();

        if let Some(root) = self.root {
            // Search recursively
            self.query_recursive(root, point, k, &mut neighbors, &mut max_dist);

            // Sort by distance (ascending), with index as tiebreaker
            neighbors.sort_by(|a, b| {
                match safe_partial_cmp(&a.0, &b.0, "kdtree sort neighbors") {
                    Ok(std::cmp::Ordering::Equal) => a.1.cmp(&b.1), // Use index as tiebreaker
                    Ok(ord) => ord,
                    Err(_) => std::cmp::Ordering::Equal,
                }
            });

            // Trim to k elements if needed
            if neighbors.len() > k {
                neighbors.truncate(k);
            }

            // Convert to sorted lists of indices and distances
            let mut indices = Vec::with_capacity(neighbors.len());
            let mut distances = Vec::with_capacity(neighbors.len());

            for (dist, idx) in neighbors {
                indices.push(idx);
                distances.push(dist);
            }

            Ok((indices, distances))
        } else {
            Err(SpatialError::ValueError("Empty tree".to_string()))
        }
    }

    /// Recursive helper for query
    fn query_recursive(
        &self,
        node_idx: usize,
        point: &[T],
        k: usize,
        neighbors: &mut Vec<(T, usize)>,
        max_dist: &mut T,
    ) {
        let node = &self.nodes[node_idx];
        let idx = node.idx;
        let axis = node.axis;

        // Calculate distance to current point
        let node_point = self.points.row(idx).to_vec();
        let _dist = self.metric.distance(&node_point, point);

        // Update neighbors if needed
        if neighbors.len() < k {
            neighbors.push((_dist, idx));

            // Sort if we just filled to capacity to establish max-heap
            if neighbors.len() == k {
                neighbors.sort_by(|a, b| {
                    match safe_partial_cmp(&b.0, &a.0, "kdtree sort max-heap") {
                        Ok(std::cmp::Ordering::Equal) => b.1.cmp(&a.1), // Use index as tiebreaker
                        Ok(ord) => ord,
                        Err(_) => std::cmp::Ordering::Equal,
                    }
                });
                *max_dist = neighbors[0].0;
            }
        } else if &_dist < max_dist {
            // Replace the worst neighbor with this one
            neighbors[0] = (_dist, idx);

            // Re-sort to maintain max-heap property
            neighbors.sort_by(|a, b| {
                match safe_partial_cmp(&b.0, &a.0, "kdtree re-sort max-heap") {
                    Ok(std::cmp::Ordering::Equal) => b.1.cmp(&a.1), // Use index as tiebreaker
                    Ok(ord) => ord,
                    Err(_) => std::cmp::Ordering::Equal,
                }
            });
            *max_dist = neighbors[0].0;
        }

        // Determine which subtree to search first
        let diff = point[axis] - node.value;
        let (first, second) = if diff < T::zero() {
            (node.left, node.right)
        } else {
            (node.right, node.left)
        };

        // Search the near subtree
        if let Some(first_idx) = first {
            self.query_recursive(first_idx, point, k, neighbors, max_dist);
        }

        // Only search the far subtree if it could contain closer points
        let axis_dist = if diff < T::zero() {
            // Point is to the left of the splitting hyperplane
            T::zero() // No need to calculate distance if we're considering the left subtree next
        } else {
            // Point is to the right of the splitting hyperplane
            diff
        };

        if let Some(second_idx) = second {
            // Only search the second subtree if necessary
            if neighbors.len() < k || axis_dist < *max_dist {
                self.query_recursive(second_idx, point, k, neighbors, max_dist);
            }
        }
    }

    /// Find all points within a radius of a query point
    ///
    /// # Arguments
    ///
    /// * `point` - Query point
    /// * `radius` - Search radius
    ///
    /// # Returns
    ///
    /// * (indices, distances) of points within the radius
    ///
    /// # Examples
    ///
    /// ```
    /// use scirs2_spatial::KDTree;
    /// use scirs2_core::ndarray::array;
    ///
    /// # fn example() -> Result<(), Box<dyn std::error::Error>> {
    /// let points = array![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
    /// let kdtree = KDTree::new(&points)?;
    ///
    /// // Find all points within radius 0.7 of [0.5, 0.5]
    /// let (indices, distances) = kdtree.query_radius(&[0.5, 0.5], 0.7)?;
    /// assert_eq!(indices.len(), 4); // All points are within 0.7 units of [0.5, 0.5]
    /// # Ok(())
    /// # }
    /// ```
    pub fn query_radius(&self, point: &[T], radius: T) -> SpatialResult<(Vec<usize>, Vec<T>)> {
        if point.len() != self.ndim {
            return Err(SpatialError::DimensionError(format!(
                "Query point dimension ({}) does not match tree dimension ({})",
                point.len(),
                self.ndim
            )));
        }

        if radius < T::zero() {
            return Err(SpatialError::ValueError(
                "Radius must be non-negative".to_string(),
            ));
        }

        let mut indices = Vec::new();
        let mut distances = Vec::new();

        if let Some(root) = self.root {
            // If the radius is outside the bounds of the entire dataset, just return an empty result
            let bounds_dist = self.bounds.min_distance(point, &self.metric);
            if bounds_dist > radius {
                return Ok((indices, distances));
            }

            // Search recursively
            self.query_radius_recursive(root, point, radius, &mut indices, &mut distances);

            // Sort by distance
            if !indices.is_empty() {
                let mut idx_dist: Vec<(usize, T)> = indices.into_iter().zip(distances).collect();
                idx_dist.sort_by(|a, b| {
                    safe_partial_cmp(&a.1, &b.1, "kdtree sort radius results")
                        .unwrap_or(std::cmp::Ordering::Equal)
                });

                indices = idx_dist.iter().map(|(idx_, _)| *idx_).collect();
                distances = idx_dist.iter().map(|(_, dist)| *dist).collect();
            }
        }

        Ok((indices, distances))
    }

    /// Recursive helper for query_radius
    fn query_radius_recursive(
        &self,
        node_idx: usize,
        point: &[T],
        radius: T,
        indices: &mut Vec<usize>,
        distances: &mut Vec<T>,
    ) {
        let node = &self.nodes[node_idx];
        let idx = node.idx;
        let axis = node.axis;

        // Calculate distance to current point
        let node_point = self.points.row(idx).to_vec();
        let dist = self.metric.distance(&node_point, point);

        // If point is within radius, add it to results
        if dist <= radius {
            indices.push(idx);
            distances.push(dist);
        }

        // Determine which subtrees need to be searched
        let diff = point[axis] - node.value;

        // Always search the near subtree
        let (near, far) = if diff < T::zero() {
            (node.left, node.right)
        } else {
            (node.right, node.left)
        };

        if let Some(near_idx) = near {
            self.query_radius_recursive(near_idx, point, radius, indices, distances);
        }

        // Only search the far subtree if it could contain points within radius
        if diff.abs() <= radius {
            if let Some(far_idx) = far {
                self.query_radius_recursive(far_idx, point, radius, indices, distances);
            }
        }
    }

    /// Count the number of points within a radius of a query point
    ///
    /// This method is more efficient than query_radius when only the count is needed.
    ///
    /// # Arguments
    ///
    /// * `point` - Query point
    /// * `radius` - Search radius
    ///
    /// # Returns
    ///
    /// * Number of points within the radius
    ///
    /// # Examples
    ///
    /// ```
    /// use scirs2_spatial::KDTree;
    /// use scirs2_core::ndarray::array;
    ///
    /// # fn example() -> Result<(), Box<dyn std::error::Error>> {
    /// let points = array![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
    /// let kdtree = KDTree::new(&points)?;
    ///
    /// // Count points within radius 0.7 of [0.5, 0.5]
    /// let count = kdtree.count_neighbors(&[0.5, 0.5], 0.7)?;
    /// assert_eq!(count, 4); // All points are within 0.7 units of [0.5, 0.5]
    /// # Ok(())
    /// # }
    /// ```
    pub fn count_neighbors(&self, point: &[T], radius: T) -> SpatialResult<usize> {
        if point.len() != self.ndim {
            return Err(SpatialError::DimensionError(format!(
                "Query point dimension ({}) does not match tree dimension ({})",
                point.len(),
                self.ndim
            )));
        }

        if radius < T::zero() {
            return Err(SpatialError::ValueError(
                "Radius must be non-negative".to_string(),
            ));
        }

        let mut count = 0;

        if let Some(root) = self.root {
            // If the radius is outside the bounds of the entire dataset, just return 0
            let bounds_dist = self.bounds.min_distance(point, &self.metric);
            if bounds_dist > radius {
                return Ok(0);
            }

            // Search recursively
            self.count_neighbors_recursive(root, point, radius, &mut count);
        }

        Ok(count)
    }

    /// Recursive helper for count_neighbors
    fn count_neighbors_recursive(
        &self,
        node_idx: usize,
        point: &[T],
        radius: T,
        count: &mut usize,
    ) {
        let node = &self.nodes[node_idx];
        let idx = node.idx;
        let axis = node.axis;

        // Calculate distance to current point
        let node_point = self.points.row(idx).to_vec();
        let dist = self.metric.distance(&node_point, point);

        // If point is within radius, increment count
        if dist <= radius {
            *count += 1;
        }

        // Determine which subtrees need to be searched
        let diff = point[axis] - node.value;

        // Always search the near subtree
        let (near, far) = if diff < T::zero() {
            (node.left, node.right)
        } else {
            (node.right, node.left)
        };

        if let Some(near_idx) = near {
            self.count_neighbors_recursive(near_idx, point, radius, count);
        }

        // Only search the far subtree if it could contain points within radius
        if diff.abs() <= radius {
            if let Some(far_idx) = far {
                self.count_neighbors_recursive(far_idx, point, radius, count);
            }
        }
    }

    /// Get the shape of the KD-Tree's point set
    ///
    /// # Returns
    ///
    /// * A tuple of (n_points, n_dimensions)
    pub fn shape(&self) -> (usize, usize) {
        (self.points.nrows(), self.ndim)
    }

    /// Get the number of points in the KD-Tree
    ///
    /// # Returns
    ///
    /// * The number of points
    pub fn npoints(&self) -> usize {
        self.points.nrows()
    }

    /// Get the dimensionality of the points in the KD-Tree
    ///
    /// # Returns
    ///
    /// * The dimensionality of the points
    pub fn ndim(&self) -> usize {
        self.ndim
    }

    /// Get the leaf size of the KD-Tree
    ///
    /// # Returns
    ///
    /// * The leaf size
    pub fn leafsize(&self) -> usize {
        self.leafsize
    }

    /// Get the bounds of the KD-Tree
    ///
    /// # Returns
    ///
    /// * The bounding rectangle of the entire dataset
    pub fn bounds(&self) -> &Rectangle<T> {
        &self.bounds
    }
}

#[cfg(test)]
mod tests {
    use super::{KDTree, Rectangle};
    use crate::distance::{
        ChebyshevDistance, Distance, EuclideanDistance, ManhattanDistance, MinkowskiDistance,
    };
    use approx::assert_relative_eq;
    use scirs2_core::ndarray::arr2;

    #[test]
    fn test_rectangle() {
        let mins = vec![0.0, 0.0];
        let maxes = vec![1.0, 1.0];
        let rect = Rectangle::new(mins, maxes);

        // Test contains
        assert!(rect.contains(&[0.5, 0.5]));
        assert!(rect.contains(&[0.0, 0.0]));
        assert!(rect.contains(&[1.0, 1.0]));
        assert!(!rect.contains(&[1.5, 0.5]));
        assert!(!rect.contains(&[0.5, 1.5]));

        // Test split
        let (left, right) = rect.split(0, 0.5);
        assert!(left.contains(&[0.25, 0.5]));
        assert!(!left.contains(&[0.75, 0.5]));
        assert!(!right.contains(&[0.25, 0.5]));
        assert!(right.contains(&[0.75, 0.5]));

        // Test min_distance
        let metric = EuclideanDistance::<f64>::new();
        assert_relative_eq!(rect.min_distance(&[0.5, 0.5], &metric), 0.0, epsilon = 1e-6);
        assert_relative_eq!(rect.min_distance(&[2.0, 0.5], &metric), 1.0, epsilon = 1e-6);
        assert_relative_eq!(
            rect.min_distance(&[2.0, 2.0], &metric),
            std::f64::consts::SQRT_2,
            epsilon = 1e-6
        );
    }

    #[test]
    fn test_kdtree_build() {
        let points = arr2(&[
            [2.0, 3.0],
            [5.0, 4.0],
            [9.0, 6.0],
            [4.0, 7.0],
            [8.0, 1.0],
            [7.0, 2.0],
        ]);

        let kdtree = KDTree::new(&points).expect("Operation failed");

        // Check that the tree has the correct number of nodes
        assert_eq!(kdtree.nodes.len(), points.nrows());

        // Check tree properties
        assert_eq!(kdtree.shape(), (6, 2));
        assert_eq!(kdtree.npoints(), 6);
        assert_eq!(kdtree.ndim(), 2);
        assert_eq!(kdtree.leafsize(), 16);

        // Check bounds
        assert_eq!(kdtree.bounds().mins(), &[2.0, 1.0]);
        assert_eq!(kdtree.bounds().maxes(), &[9.0, 7.0]);
    }

    #[test]
    fn test_kdtree_query() {
        let points = arr2(&[
            [2.0, 3.0],
            [5.0, 4.0],
            [9.0, 6.0],
            [4.0, 7.0],
            [8.0, 1.0],
            [7.0, 2.0],
        ]);

        let kdtree = KDTree::new(&points).expect("Operation failed");

        // Query for nearest neighbor to [3.0, 5.0]
        let (indices, distances) = kdtree.query(&[3.0, 5.0], 1).expect("Operation failed");
        assert_eq!(indices.len(), 1);
        assert_eq!(distances.len(), 1);

        // Calculate actual distances
        let query = [3.0, 5.0];
        let mut expected_dists = vec![];
        for i in 0..points.nrows() {
            let p = points.row(i).to_vec();
            let metric = EuclideanDistance::<f64>::new();
            expected_dists.push((i, metric.distance(&p, &query)));
        }
        expected_dists.sort_by(|a, b| a.1.partial_cmp(&b.1).unwrap_or(std::cmp::Ordering::Equal));

        // Verify we got one of the actual nearest neighbors (there might be ties)
        // Check that the distance matches the expected minimum distance
        assert_relative_eq!(distances[0], expected_dists[0].1, epsilon = 1e-6);

        // Verify the returned index is one of the points with minimum distance
        let min_dist = expected_dists[0].1;
        let valid_indices: Vec<usize> = expected_dists
            .iter()
            .filter(|(_, d)| (d - min_dist).abs() < 1e-6)
            .map(|(i, _)| *i)
            .collect();
        assert!(
            valid_indices.contains(&indices[0]),
            "Expected one of {:?} but got {}",
            valid_indices,
            indices[0]
        );
    }

    #[test]
    fn test_kdtree_query_k() {
        let points = arr2(&[
            [2.0, 3.0],
            [5.0, 4.0],
            [9.0, 6.0],
            [4.0, 7.0],
            [8.0, 1.0],
            [7.0, 2.0],
        ]);

        let kdtree = KDTree::new(&points).expect("Operation failed");

        // Query for 3 nearest neighbors to [3.0, 5.0]
        let (indices, distances) = kdtree.query(&[3.0, 5.0], 3).expect("Operation failed");
        assert_eq!(indices.len(), 3);
        assert_eq!(distances.len(), 3);

        // Calculate actual distances
        let query = [3.0, 5.0];
        let mut expected_dists = vec![];
        for i in 0..points.nrows() {
            let p = points.row(i).to_vec();
            let metric = EuclideanDistance::<f64>::new();
            expected_dists.push((i, metric.distance(&p, &query)));
        }
        expected_dists.sort_by(|a, b| a.1.partial_cmp(&b.1).unwrap_or(std::cmp::Ordering::Equal));

        // Verify we got the 3 actual nearest neighbors (for now, just check distances)
        let expected_indices: Vec<usize> = expected_dists.iter().take(3).map(|&(i, _)| i).collect();
        let expected_distances: Vec<f64> = expected_dists.iter().take(3).map(|&(_, d)| d).collect();

        // Check each returned index is in the expected set
        for i in &indices {
            assert!(expected_indices.contains(i));
        }

        // Check that distances are sorted
        assert!(distances[0] <= distances[1]);
        assert!(distances[1] <= distances[2]);

        // Check distance values match expected
        for i in 0..3 {
            assert_relative_eq!(distances[i], expected_distances[i], epsilon = 1e-6);
        }
    }

    #[test]
    fn test_kdtree_query_radius() {
        let points = arr2(&[
            [2.0, 3.0],
            [5.0, 4.0],
            [9.0, 6.0],
            [4.0, 7.0],
            [8.0, 1.0],
            [7.0, 2.0],
        ]);

        let kdtree = KDTree::new(&points).expect("Operation failed");

        // Query for points within radius 3.0 of [3.0, 5.0]
        let (indices, distances) = kdtree
            .query_radius(&[3.0, 5.0], 3.0)
            .expect("Operation failed");

        // Calculate expected results
        let query = [3.0, 5.0];
        let radius = 3.0;
        let mut expected_results = vec![];
        for i in 0..points.nrows() {
            let p = points.row(i).to_vec();
            let metric = EuclideanDistance::<f64>::new();
            let dist = metric.distance(&p, &query);
            if dist <= radius {
                expected_results.push((i, dist));
            }
        }
        expected_results.sort_by(|a, b| {
            match a.1.partial_cmp(&b.1).unwrap_or(std::cmp::Ordering::Equal) {
                std::cmp::Ordering::Equal => a.0.cmp(&b.0), // Use index as tiebreaker
                ord => ord,
            }
        });

        // Check that we got the expected number of points
        assert_eq!(indices.len(), expected_results.len());

        // Check that all returned points are within radius
        for i in 0..indices.len() {
            assert!(distances[i] <= radius + 1e-6);
        }

        // Check that the indices/distances pairs match expected results
        // Note: order might differ for equal distances
        let mut idx_dist_pairs: Vec<(usize, f64)> = indices
            .iter()
            .zip(distances.iter())
            .map(|(&i, &d)| (i, d))
            .collect();
        idx_dist_pairs.sort_by(|a, b| {
            match a.1.partial_cmp(&b.1).unwrap_or(std::cmp::Ordering::Equal) {
                std::cmp::Ordering::Equal => a.0.cmp(&b.0),
                ord => ord,
            }
        });

        for (actual, expected) in idx_dist_pairs.iter().zip(expected_results.iter()) {
            assert_eq!(actual.0, expected.0);
            assert_relative_eq!(actual.1, expected.1, epsilon = 1e-6);
        }
    }

    #[test]
    fn test_kdtree_count_neighbors() {
        let points = arr2(&[
            [2.0, 3.0],
            [5.0, 4.0],
            [9.0, 6.0],
            [4.0, 7.0],
            [8.0, 1.0],
            [7.0, 2.0],
        ]);

        let kdtree = KDTree::new(&points).expect("Operation failed");

        // Count points within radius 3.0 of [3.0, 5.0]
        let count = kdtree
            .count_neighbors(&[3.0, 5.0], 3.0)
            .expect("Operation failed");

        // Calculate actual count
        let query = [3.0, 5.0];
        let mut expected_count = 0;
        for i in 0..points.nrows() {
            let p = points.row(i).to_vec();
            let metric = EuclideanDistance::<f64>::new();
            let dist = metric.distance(&p, &query);
            if dist <= 3.0 {
                expected_count += 1;
            }
        }

        assert_eq!(count, expected_count);
    }

    #[test]
    fn test_kdtree_with_manhattan_distance() {
        let points = arr2(&[
            [2.0, 3.0],
            [5.0, 4.0],
            [9.0, 6.0],
            [4.0, 7.0],
            [8.0, 1.0],
            [7.0, 2.0],
        ]);

        let metric = ManhattanDistance::new();
        let kdtree = KDTree::with_metric(&points, metric).expect("Operation failed");

        // Query for nearest neighbor to [3.0, 5.0] using Manhattan distance
        let (indices, distances) = kdtree.query(&[3.0, 5.0], 1).expect("Operation failed");

        // Calculate actual distances using Manhattan distance
        let query = [3.0, 5.0];
        let mut expected_dists = vec![];
        for i in 0..points.nrows() {
            let p = points.row(i).to_vec();
            let m = ManhattanDistance::<f64>::new();
            expected_dists.push((i, m.distance(&p, &query)));
        }
        expected_dists.sort_by(|a, b| a.1.partial_cmp(&b.1).unwrap_or(std::cmp::Ordering::Equal));

        // Check that the distance matches the expected minimum distance
        assert_relative_eq!(distances[0], expected_dists[0].1, epsilon = 1e-6);

        // Verify the returned index is one of the points with minimum distance
        let min_dist = expected_dists[0].1;
        let valid_indices: Vec<usize> = expected_dists
            .iter()
            .filter(|(_, d)| (d - min_dist).abs() < 1e-6)
            .map(|(i, _)| *i)
            .collect();
        assert!(
            valid_indices.contains(&indices[0]),
            "Expected one of {:?} but got {}",
            valid_indices,
            indices[0]
        );
    }

    #[test]
    fn test_kdtree_with_chebyshev_distance() {
        let points = arr2(&[
            [2.0, 3.0],
            [5.0, 4.0],
            [9.0, 6.0],
            [4.0, 7.0],
            [8.0, 1.0],
            [7.0, 2.0],
        ]);

        let metric = ChebyshevDistance::new();
        let kdtree = KDTree::with_metric(&points, metric).expect("Operation failed");

        // Query for nearest neighbor to [3.0, 5.0] using Chebyshev distance
        let (indices, distances) = kdtree.query(&[3.0, 5.0], 1).expect("Operation failed");

        // Calculate actual distances using Chebyshev distance
        let query = [3.0, 5.0];
        let mut expected_dists = vec![];
        for i in 0..points.nrows() {
            let p = points.row(i).to_vec();
            let m = ChebyshevDistance::<f64>::new();
            expected_dists.push((i, m.distance(&p, &query)));
        }
        expected_dists.sort_by(|a, b| a.1.partial_cmp(&b.1).unwrap_or(std::cmp::Ordering::Equal));

        // Check that the distance matches the expected minimum distance
        assert_relative_eq!(distances[0], expected_dists[0].1, epsilon = 1e-6);

        // Verify the returned index is one of the points with minimum distance
        let min_dist = expected_dists[0].1;
        let valid_indices: Vec<usize> = expected_dists
            .iter()
            .filter(|(_, d)| (d - min_dist).abs() < 1e-6)
            .map(|(i, _)| *i)
            .collect();
        assert!(
            valid_indices.contains(&indices[0]),
            "Expected one of {:?} but got {}",
            valid_indices,
            indices[0]
        );
    }

    #[test]
    fn test_kdtree_with_minkowski_distance() {
        let points = arr2(&[
            [2.0, 3.0],
            [5.0, 4.0],
            [9.0, 6.0],
            [4.0, 7.0],
            [8.0, 1.0],
            [7.0, 2.0],
        ]);

        let metric = MinkowskiDistance::new(3.0);
        let kdtree = KDTree::with_metric(&points, metric).expect("Operation failed");

        // Query for nearest neighbor to [3.0, 5.0] using Minkowski distance (p=3)
        let (indices, distances) = kdtree.query(&[3.0, 5.0], 1).expect("Operation failed");

        // Calculate actual distances using Minkowski distance
        let query = [3.0, 5.0];
        let mut expected_dists = vec![];
        for i in 0..points.nrows() {
            let p = points.row(i).to_vec();
            let m = MinkowskiDistance::<f64>::new(3.0);
            expected_dists.push((i, m.distance(&p, &query)));
        }
        expected_dists.sort_by(|a, b| a.1.partial_cmp(&b.1).unwrap_or(std::cmp::Ordering::Equal));

        // Check that the distance matches the expected minimum distance
        assert_relative_eq!(distances[0], expected_dists[0].1, epsilon = 1e-6);

        // Verify the returned index is one of the points with minimum distance
        let min_dist = expected_dists[0].1;
        let valid_indices: Vec<usize> = expected_dists
            .iter()
            .filter(|(_, d)| (d - min_dist).abs() < 1e-6)
            .map(|(i, _)| *i)
            .collect();
        assert!(
            valid_indices.contains(&indices[0]),
            "Expected one of {:?} but got {}",
            valid_indices,
            indices[0]
        );
    }

    #[test]
    fn test_kdtree_with_custom_leaf_size() {
        let points = arr2(&[
            [2.0, 3.0],
            [5.0, 4.0],
            [9.0, 6.0],
            [4.0, 7.0],
            [8.0, 1.0],
            [7.0, 2.0],
        ]);

        // Use a very small leaf size to test that it works
        let leafsize = 1;
        let kdtree = KDTree::with_leaf_size(&points, leafsize).expect("Operation failed");

        assert_eq!(kdtree.leafsize(), 1);

        // Query for nearest neighbor to [3.0, 5.0]
        let (indices, distances) = kdtree.query(&[3.0, 5.0], 1).expect("Operation failed");

        // Calculate actual distances
        let query = [3.0, 5.0];
        let mut expected_dists = vec![];
        for i in 0..points.nrows() {
            let p = points.row(i).to_vec();
            let metric = EuclideanDistance::<f64>::new();
            expected_dists.push((i, metric.distance(&p, &query)));
        }
        expected_dists.sort_by(|a, b| a.1.partial_cmp(&b.1).unwrap_or(std::cmp::Ordering::Equal));

        // Verify we got one of the actual nearest neighbors (there might be ties)
        // Check that the distance matches the expected minimum distance
        assert_relative_eq!(distances[0], expected_dists[0].1, epsilon = 1e-6);

        // Verify the returned index is one of the points with minimum distance
        let min_dist = expected_dists[0].1;
        let valid_indices: Vec<usize> = expected_dists
            .iter()
            .filter(|(_, d)| (d - min_dist).abs() < 1e-6)
            .map(|(i, _)| *i)
            .collect();
        assert!(
            valid_indices.contains(&indices[0]),
            "Expected one of {:?} but got {}",
            valid_indices,
            indices[0]
        );
    }
}