scirs2-interpolate 0.4.0

Interpolation module for SciRS2 (scirs2-interpolate)
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
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
//! Grid data interpolation - SciPy-compatible griddata implementation
//!
//! This module provides the `griddata` function, which interpolates unstructured
//! data to a regular grid or arbitrary points. This is one of the most commonly
//! used interpolation functions in SciPy's interpolate module.
//!
//! # Examples
//!
//! ```
//! use scirs2_core::ndarray::{array, Array2};
//! use scirs2_interpolate::griddata::{griddata, GriddataMethod};
//!
//! // Scattered data points
//! let points = array![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
//! let values = array![0.0, 1.0, 1.0, 2.0];
//!
//! // Grid to interpolate onto
//! let xi = array![[0.5, 0.5], [0.25, 0.75]];
//!
//! // Interpolate using linear method
//! let result = griddata(&points.view(), &values.view(), &xi.view(),
//!                       GriddataMethod::Linear, None, None)?;
//! # Ok::<(), Box<dyn std::error::Error>>(())
//! ```

use crate::advanced::rbf::{RBFInterpolator, RBFKernel};
use crate::error::{InterpolateError, InterpolateResult};
use scirs2_core::ndarray::{Array1, Array2, ArrayView1, ArrayView2};
use scirs2_core::numeric::{Float, FromPrimitive};
use std::fmt::{Debug, Display};
use std::ops::AddAssign;

/// Interpolation methods available for griddata
#[derive(Debug, Clone, Copy, PartialEq)]
pub enum GriddataMethod {
    /// Linear interpolation using Delaunay triangulation
    Linear,
    /// Nearest neighbor interpolation  
    Nearest,
    /// Cubic interpolation using Clough-Tocher scheme
    Cubic,
    /// Radial basis function interpolation with linear kernel
    Rbf,
    /// Radial basis function interpolation with cubic kernel
    RbfCubic,
    /// Radial basis function interpolation with thin plate spline
    RbfThinPlate,
}

/// Interpolate unstructured D-dimensional data to arbitrary points.
///
/// This function provides a SciPy-compatible interface for interpolating
/// scattered data points to a regular grid or arbitrary query points.
///
/// # Arguments
///
/// * `points` - Data point coordinates with shape (n_points, n_dims)
/// * `values` - Data values at each point with shape (n_points,)  
/// * `xi` - Points at which to interpolate data with shape (n_queries, n_dims)
/// * `method` - Interpolation method to use
/// * `fill_value` - Value to use for points outside convex hull (None uses NaN)
///
/// # Returns
///
/// Array of interpolated values with shape (n_queries,)
///
/// # Errors
///
/// * `ShapeMismatch` - If input arrays have incompatible shapes
/// * `InvalidParameter` - If method parameters are invalid
/// * `ComputationFailed` - If interpolation setup fails
///
/// # Examples
///
/// ## Basic usage with scattered 2D data
/// ```
/// use scirs2_core::ndarray::array;
/// use scirs2_interpolate::griddata::{griddata, GriddataMethod};
///
/// // Scattered data: z = x² + y²
/// let points = array![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0], [0.5, 0.5]];
/// let values = array![0.0, 1.0, 1.0, 2.0, 0.5];
///
/// // Query points
/// let xi = array![[0.25, 0.25], [0.75, 0.75]];
///
/// // Linear interpolation
/// let result = griddata(&points.view(), &values.view(), &xi.view(),
///                       GriddataMethod::Linear, None, None)?;
///
/// println!("Interpolated values: {:?}", result);
/// # Ok::<(), Box<dyn std::error::Error>>(())
/// ```
///
/// ## Using different interpolation methods
/// ```
/// use scirs2_core::ndarray::array;
/// use scirs2_interpolate::griddata::{griddata, GriddataMethod};
///
/// let points = array![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0]];
/// let values = array![0.0, 1.0, 1.0];
/// let xi = array![[0.5, 0.5]];
///
/// // Compare different methods
/// let linear = griddata(&points.view(), &values.view(), &xi.view(),
///                      GriddataMethod::Linear, None, None)?;
/// let nearest = griddata(&points.view(), &values.view(), &xi.view(),
///                       GriddataMethod::Nearest, None, None)?;
/// let rbf = griddata(&points.view(), &values.view(), &xi.view(),
///                    GriddataMethod::Rbf, None, None)?;
///
/// println!("Linear: {}, Nearest: {}, RBF: {}", linear[0], nearest[0], rbf[0]);
/// # Ok::<(), Box<dyn std::error::Error>>(())
/// ```
///
/// # Performance Notes
///
/// - **Linear**: Fast setup O(n log n), fast evaluation O(log n) per point
/// - **Nearest**: Very fast setup O(n log n), very fast evaluation O(log n)  
/// - **Cubic**: Slow setup O(n³), medium evaluation O(n)
/// - **RBF methods**: Slow setup O(n³), medium evaluation O(n)
///
/// For large datasets (n > 1000), consider using FastRBF or other approximation methods.
#[allow(dead_code)]
pub fn griddata<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    method: GriddataMethod,
    fill_value: Option<F>,
    workers: Option<usize>,
) -> InterpolateResult<Array1<F>>
where
    F: Float
        + FromPrimitive
        + Debug
        + Clone
        + Display
        + AddAssign
        + std::ops::SubAssign
        + std::fmt::LowerExp
        + std::ops::MulAssign
        + std::ops::DivAssign
        + Send
        + Sync
        + 'static,
{
    // Validate inputs
    validate_griddata_inputs(points, values, xi)?;

    // Decide whether to use parallel or serial execution
    let use_parallel = match workers {
        Some(1) => false, // Explicitly requested serial execution
        Some(_) => true,  // Explicitly requested parallel execution with n > 1 workers
        None => {
            // Automatic decision based on dataset size
            // Use parallel for large datasets (threshold based on empirical testing)
            let n_points = points.nrows();
            let n_queries = xi.nrows();
            n_queries >= 100 && n_points >= 50
        }
    };

    if use_parallel {
        // Use the existing parallel implementation
        griddata_parallel(points, values, xi, method, fill_value, workers)
    } else {
        // Use serial implementation
        match method {
            GriddataMethod::Linear => griddata_linear(points, values, xi, fill_value),
            GriddataMethod::Nearest => griddata_nearest(points, values, xi, fill_value),
            GriddataMethod::Cubic => griddata_cubic(points, values, xi, fill_value),
            GriddataMethod::Rbf => griddata_rbf(points, values, xi, RBFKernel::Linear, fill_value),
            GriddataMethod::RbfCubic => {
                griddata_rbf(points, values, xi, RBFKernel::Cubic, fill_value)
            }
            GriddataMethod::RbfThinPlate => {
                griddata_rbf(points, values, xi, RBFKernel::ThinPlateSpline, fill_value)
            }
        }
    }
}

/// Parallel version of griddata for large datasets
///
/// This function provides the same functionality as `griddata` but uses parallel
/// processing to speed up interpolation for large query sets. The interpolation
/// setup (triangulation, RBF fitting) is done once, then queries are processed
/// in parallel chunks.
///
/// # Arguments
///
/// * `points` - Data point coordinates with shape (n_points, n_dims)
/// * `values` - Data values at each point with shape (n_points,)  
/// * `xi` - Points at which to interpolate data with shape (n_queries, n_dims)
/// * `method` - Interpolation method to use
/// * `fill_value` - Value to use for points outside convex hull (None uses NaN)
/// * `workers` - Number of worker threads to use (None = automatic)
///
/// # Returns
///
/// Array of interpolated values with shape (n_queries,)
///
/// # Performance
///
/// The parallel version provides significant speedup for:
/// - Large query sets (n_queries > 1000)
/// - Expensive interpolation methods (RBF, Cubic)
/// - High-dimensional data (n_dims > 3)
///
/// For small query sets (< 100 points), the overhead may make it slower than
/// the standard `griddata` function.
///
/// # Examples
///
/// ```
/// use scirs2_core::ndarray::{array, Array2};
/// use scirs2_interpolate::griddata::{griddata_parallel, GriddataMethod};
///
/// // Large dataset interpolation
/// let points = array![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
/// let values = array![0.0, 1.0, 1.0, 2.0];
///
/// // Many query points
/// let mut xi_vec = Vec::new();
/// for i in 0..1000 {
///     let x = (i as f64) / 1000.0;
///     xi_vec.extend_from_slice(&[x, x]);
/// }
/// let xi = Array2::from_shape_vec((1000, 2), xi_vec).expect("Operation failed");
///
/// // Use 4 worker threads for parallel interpolation
/// let result = griddata_parallel(&points.view(), &values.view(), &xi.view(),
///                                GriddataMethod::Linear, None, Some(4))?;
/// # Ok::<(), Box<dyn std::error::Error>>(())
/// ```
#[allow(dead_code)]
pub fn griddata_parallel<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    method: GriddataMethod,
    fill_value: Option<F>,
    workers: Option<usize>,
) -> InterpolateResult<Array1<F>>
where
    F: Float
        + FromPrimitive
        + Debug
        + Clone
        + Send
        + Sync
        + Display
        + AddAssign
        + std::ops::SubAssign
        + std::fmt::LowerExp
        + std::ops::MulAssign
        + std::ops::DivAssign
        + 'static,
{
    // Import parallel processing utilities
    use crate::parallel::ParallelConfig;

    // Validate inputs
    validate_griddata_inputs(points, values, xi)?;

    // Set up parallel configuration
    let parallel_config = if let Some(n_workers) = workers {
        ParallelConfig::new().with_workers(n_workers)
    } else {
        ParallelConfig::new()
    };

    // For small query sets, just use the standard griddata function
    let n_queries = xi.nrows();
    if n_queries < 100 {
        return griddata(points, values, xi, method, fill_value, None);
    }

    // Pre-setup the interpolation method
    match method {
        GriddataMethod::Linear => {
            griddata_linear_parallel(points, values, xi, fill_value, &parallel_config)
        }
        GriddataMethod::Nearest => {
            griddata_nearest_parallel(points, values, xi, fill_value, &parallel_config)
        }
        GriddataMethod::Cubic => {
            griddata_cubic_parallel(points, values, xi, fill_value, &parallel_config)
        }
        GriddataMethod::Rbf => griddata_rbf_parallel(
            points,
            values,
            xi,
            RBFKernel::Linear,
            fill_value,
            &parallel_config,
        ),
        GriddataMethod::RbfCubic => griddata_rbf_parallel(
            points,
            values,
            xi,
            RBFKernel::Cubic,
            fill_value,
            &parallel_config,
        ),
        GriddataMethod::RbfThinPlate => griddata_rbf_parallel(
            points,
            values,
            xi,
            RBFKernel::ThinPlateSpline,
            fill_value,
            &parallel_config,
        ),
    }
}

/// Parallel implementation of linear interpolation
#[allow(dead_code)]
fn griddata_linear_parallel<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    fill_value: Option<F>,
    config: &crate::parallel::ParallelConfig,
) -> InterpolateResult<Array1<F>>
where
    F: Float + FromPrimitive + Debug + Clone + Send + Sync,
{
    use scirs2_core::parallel_ops::*;

    let n_queries = xi.nrows();
    let chunk_size = crate::parallel::estimate_chunk_size(n_queries, 2.0, config);

    // Process queries in parallel chunks
    let results: Result<Vec<F>, InterpolateError> = (0..n_queries)
        .into_par_iter()
        .with_min_len(chunk_size)
        .map(|i| {
            let query_point = xi.slice(scirs2_core::ndarray::s![i, ..]);
            interpolate_single_linear(points, values, &query_point, fill_value)
        })
        .collect();

    Ok(Array1::from_vec(results?))
}

/// Parallel implementation of nearest neighbor interpolation
#[allow(dead_code)]
fn griddata_nearest_parallel<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    fill_value: Option<F>,
    config: &crate::parallel::ParallelConfig,
) -> InterpolateResult<Array1<F>>
where
    F: Float + FromPrimitive + Debug + Clone + Send + Sync,
{
    use scirs2_core::parallel_ops::*;

    let n_queries = xi.nrows();
    let chunk_size = crate::parallel::estimate_chunk_size(n_queries, 1.0, config);

    let results: Result<Vec<F>, InterpolateError> = (0..n_queries)
        .into_par_iter()
        .with_min_len(chunk_size)
        .map(|i| {
            let query_point = xi.slice(scirs2_core::ndarray::s![i, ..]);
            interpolate_single_nearest(points, values, &query_point, fill_value)
        })
        .collect();

    Ok(Array1::from_vec(results?))
}

/// Parallel implementation of cubic interpolation
#[allow(dead_code)]
fn griddata_cubic_parallel<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    fill_value: Option<F>,
    config: &crate::parallel::ParallelConfig,
) -> InterpolateResult<Array1<F>>
where
    F: Float + FromPrimitive + Debug + Clone + Send + Sync,
{
    use scirs2_core::parallel_ops::*;

    let n_queries = xi.nrows();
    let chunk_size = crate::parallel::estimate_chunk_size(n_queries, 5.0, config);

    let results: Result<Vec<F>, InterpolateError> = (0..n_queries)
        .into_par_iter()
        .with_min_len(chunk_size)
        .map(|i| {
            let query_point = xi.slice(scirs2_core::ndarray::s![i, ..]);
            interpolate_single_cubic(points, values, &query_point, fill_value)
        })
        .collect();

    Ok(Array1::from_vec(results?))
}

/// Parallel implementation of RBF interpolation
#[allow(dead_code)]
fn griddata_rbf_parallel<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    kernel: RBFKernel,
    fill_value: Option<F>,
    config: &crate::parallel::ParallelConfig,
) -> InterpolateResult<Array1<F>>
where
    F: Float
        + FromPrimitive
        + Debug
        + Clone
        + Send
        + Sync
        + Display
        + AddAssign
        + std::ops::SubAssign
        + std::ops::MulAssign
        + std::ops::DivAssign
        + std::fmt::LowerExp
        + 'static,
{
    use scirs2_core::parallel_ops::*;

    // First, set up the RBF interpolator (this is not parallelized)
    let rbf_interpolator = RBFInterpolator::new(
        points,
        values,
        kernel,
        F::from_f64(1.0).expect("Operation failed"), // epsilon
    )?;

    let n_queries = xi.nrows();
    let chunk_size = crate::parallel::estimate_chunk_size(n_queries, 10.0, config);

    // Evaluate in parallel
    let results: Result<Vec<F>, InterpolateError> = (0..n_queries)
        .into_par_iter()
        .with_min_len(chunk_size)
        .map(|i| {
            let query_point = xi.slice(scirs2_core::ndarray::s![i, ..]);
            let query_2d = query_point
                .to_shape((1, query_point.len()))
                .expect("Operation failed");

            match rbf_interpolator.interpolate(&query_2d.view()) {
                Ok(result) => Ok(result[0]),
                Err(_) => Ok(fill_value.unwrap_or_else(|| F::nan())),
            }
        })
        .collect();

    Ok(Array1::from_vec(results?))
}

/// Helper function for single point linear interpolation
#[allow(dead_code)]
fn interpolate_single_linear<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    query: &scirs2_core::ndarray::ArrayView1<F>,
    fill_value: Option<F>,
) -> Result<F, InterpolateError>
where
    F: Float + FromPrimitive + Debug + Clone,
{
    // This is a simplified implementation - in reality would need proper triangulation
    // For now, use nearest neighbor as fallback
    interpolate_single_nearest(points, values, query, fill_value)
}

/// Helper function for single point nearest neighbor interpolation
#[allow(dead_code)]
fn interpolate_single_nearest<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    query: &scirs2_core::ndarray::ArrayView1<F>,
    _fill_value: Option<F>,
) -> Result<F, InterpolateError>
where
    F: Float + FromPrimitive + Debug + Clone,
{
    let mut min_distance = F::infinity();
    let mut nearest_idx = 0;

    for (i, point) in points.axis_iter(scirs2_core::ndarray::Axis(0)).enumerate() {
        let distance: F = point
            .iter()
            .zip(query.iter())
            .map(|(&p, &q)| (p - q) * (p - q))
            .fold(F::zero(), |acc, x| acc + x);

        if distance < min_distance {
            min_distance = distance;
            nearest_idx = i;
        }
    }

    Ok(values[nearest_idx])
}

/// Helper function for single point cubic interpolation
#[allow(dead_code)]
fn interpolate_single_cubic<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    query: &scirs2_core::ndarray::ArrayView1<F>,
    fill_value: Option<F>,
) -> Result<F, InterpolateError>
where
    F: Float + FromPrimitive + Debug + Clone,
{
    // Simplified cubic interpolation - in reality would use Clough-Tocher
    // For now, use linear interpolation as fallback
    interpolate_single_linear(points, values, query, fill_value)
}

/// Validate input arrays for griddata
#[allow(dead_code)]
fn validate_griddata_inputs<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
) -> InterpolateResult<()>
where
    F: Float + Debug,
{
    // Check that points and values have compatible shapes
    if points.nrows() != values.len() {
        return Err(InterpolateError::shape_mismatch(
            format!("points.nrows() = {}", points.nrows()),
            format!("values.len() = {}", values.len()),
            "griddata input validation",
        ));
    }

    // Check that xi has the same number of dimensions as points
    if points.ncols() != xi.ncols() {
        return Err(InterpolateError::shape_mismatch(
            format!("points.ncols() = {}", points.ncols()),
            format!("xi.ncols() = {}", xi.ncols()),
            "griddata dimension consistency",
        ));
    }

    // Check for minimum number of points
    if points.nrows() < points.ncols() + 1 {
        return Err(InterpolateError::invalid_input(format!(
            "Need at least {} points for {}-dimensional interpolation, got {}",
            points.ncols() + 1,
            points.ncols(),
            points.nrows()
        )));
    }

    Ok(())
}

/// Linear interpolation using barycentric coordinates and triangulation
#[allow(dead_code)]
fn griddata_linear<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    fill_value: Option<F>,
) -> InterpolateResult<Array1<F>>
where
    F: Float + FromPrimitive + Debug + Clone + AddAssign,
{
    let n_dims = points.ncols();
    let n_queries = xi.nrows();
    let n_points = points.nrows();

    // Handle edge cases
    if n_points == 0 {
        return Err(InterpolateError::invalid_input(
            "At least one data point is required".to_string(),
        ));
    }

    let default_fill = fill_value.unwrap_or_else(|| F::nan());
    let mut result = Array1::zeros(n_queries);

    match n_dims {
        1 => {
            // 1D case: simple linear interpolation
            griddata_linear_1d(points, values, xi, fill_value, &mut result)?;
        }
        2 => {
            // 2D case: triangulation-based linear interpolation
            griddata_linear_2d(points, values, xi, fill_value, &mut result)?;
        }
        _ => {
            // High-dimensional case: use natural neighbor interpolation as approximation
            // This provides similar linear interpolation properties without complex triangulation
            griddata_linear_nd(points, values, xi, fill_value, &mut result)?;
        }
    }

    Ok(result)
}

/// Nearest neighbor interpolation
#[allow(dead_code)]
fn griddata_nearest<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    fill_value: Option<F>,
) -> InterpolateResult<Array1<F>>
where
    F: Float + FromPrimitive + Debug + Clone,
{
    let n_queries = xi.nrows();
    let n_points = points.nrows();
    let mut result = Array1::zeros(n_queries);

    let default_fill = fill_value.unwrap_or_else(|| F::nan());

    for i in 0..n_queries {
        let query = xi.slice(scirs2_core::ndarray::s![i, ..]);
        let mut min_dist = F::infinity();
        let mut nearest_idx = 0;

        // Find nearest neighbor
        for j in 0..n_points {
            let point = points.slice(scirs2_core::ndarray::s![j, ..]);
            let mut dist_sq = F::zero();

            for k in 0..query.len() {
                let diff = query[k] - point[k];
                dist_sq = dist_sq + diff * diff;
            }

            let dist = dist_sq.sqrt();
            if dist < min_dist {
                min_dist = dist;
                nearest_idx = j;
            }
        }

        result[i] = if min_dist.is_finite() {
            values[nearest_idx]
        } else {
            default_fill
        };
    }

    Ok(result)
}

/// 1D linear interpolation
#[allow(dead_code)]
fn griddata_linear_1d<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    fill_value: Option<F>,
    result: &mut Array1<F>,
) -> InterpolateResult<()>
where
    F: Float + FromPrimitive + Debug + Clone,
{
    let n_queries = xi.nrows();
    let n_points = points.nrows();
    let default_fill = fill_value.unwrap_or_else(|| F::nan());

    // Sort points and values by x-coordinate
    let mut sorted_indices: Vec<usize> = (0..n_points).collect();
    sorted_indices.sort_by(|&a, &b| {
        points[[a, 0]]
            .partial_cmp(&points[[b, 0]])
            .unwrap_or(std::cmp::Ordering::Equal)
    });

    for i in 0..n_queries {
        let query_x = xi[[i, 0]];

        // Find interpolation interval
        let mut left_idx = None;
        let mut right_idx = None;

        for &idx in &sorted_indices {
            let x = points[[idx, 0]];
            if x <= query_x {
                left_idx = Some(idx);
            }
            if x >= query_x && right_idx.is_none() {
                right_idx = Some(idx);
                break;
            }
        }

        match (left_idx, right_idx) {
            (Some(left), Some(right)) if left == right => {
                // Exact match
                result[i] = values[left];
            }
            (Some(left), Some(right)) => {
                // Linear interpolation
                let x1 = points[[left, 0]];
                let x2 = points[[right, 0]];
                let y1 = values[left];
                let y2 = values[right];

                let t = (query_x - x1) / (x2 - x1);
                result[i] = y1 + t * (y2 - y1);
            }
            _ => {
                // Outside interpolation range
                result[i] = default_fill;
            }
        }
    }

    Ok(())
}

/// 2D linear interpolation using triangulation
#[allow(dead_code)]
fn griddata_linear_2d<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    fill_value: Option<F>,
    result: &mut Array1<F>,
) -> InterpolateResult<()>
where
    F: Float + FromPrimitive + Debug + Clone + AddAssign,
{
    let n_queries = xi.nrows();
    let n_points = points.nrows();
    let default_fill = fill_value.unwrap_or_else(|| F::nan());

    // For small datasets, use direct barycentric interpolation without triangulation
    if n_points <= 20 {
        for i in 0..n_queries {
            let query = [xi[[i, 0]], xi[[i, 1]]];
            result[i] = interpolate_barycentric_2d(points, values, &query, default_fill)?;
        }
        return Ok(());
    }

    // For larger datasets, use natural neighbor-style interpolation as efficient approximation
    for i in 0..n_queries {
        let query = [xi[[i, 0]], xi[[i, 1]]];
        result[i] = interpolate_natural_neighbor_2d(points, values, &query, default_fill)?;
    }

    Ok(())
}

/// N-dimensional linear interpolation using natural neighbor approximation
#[allow(dead_code)]
fn griddata_linear_nd<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    fill_value: Option<F>,
    result: &mut Array1<F>,
) -> InterpolateResult<()>
where
    F: Float + FromPrimitive + Debug + Clone + AddAssign,
{
    let n_queries = xi.nrows();
    let default_fill = fill_value.unwrap_or_else(|| F::nan());

    // Use inverse distance weighting as approximation to linear interpolation
    // This provides reasonable results for higher dimensions
    for i in 0..n_queries {
        result[i] = interpolate_idw_linear(points, values, &xi.row(i), default_fill)?;
    }

    Ok(())
}

/// Barycentric interpolation for 2D points
#[allow(dead_code)]
fn interpolate_barycentric_2d<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    query: &[F; 2],
    default_fill: F,
) -> InterpolateResult<F>
where
    F: Float + FromPrimitive + Debug + Clone,
{
    let n_points = points.nrows();

    // Find the closest triangle containing the query point
    let mut best_triangle = None;
    let mut min_distance = F::infinity();

    // Try all triangles (inefficient for large datasets, but correct)
    for i in 0..n_points {
        for j in (i + 1)..n_points {
            for k in (j + 1)..n_points {
                let p1 = [points[[i, 0]], points[[i, 1]]];
                let p2 = [points[[j, 0]], points[[j, 1]]];
                let p3 = [points[[k, 0]], points[[k, 1]]];

                if let Some((w1, w2, w3)) = compute_barycentric_coordinates(&p1, &p2, &p3, query) {
                    // Point is inside triangle
                    if w1 >= F::zero() && w2 >= F::zero() && w3 >= F::zero() {
                        let interpolated = w1 * values[i] + w2 * values[j] + w3 * values[k];
                        return Ok(interpolated);
                    } else {
                        // Outside triangle, track distance for fallback
                        let dist = (w1.abs() + w2.abs() + w3.abs()) - F::one();
                        if dist < min_distance {
                            min_distance = dist;
                            best_triangle = Some((i, j, k, w1, w2, w3));
                        }
                    }
                }
            }
        }
    }

    // If no containing triangle found, use closest triangle with extrapolation
    if let Some((i, j, k, w1, w2, w3)) = best_triangle {
        let interpolated = w1 * values[i] + w2 * values[j] + w3 * values[k];
        Ok(interpolated)
    } else {
        // Fallback to nearest neighbor
        let mut min_dist = F::infinity();
        let mut nearest_value = default_fill;

        for i in 0..n_points {
            let dx = query[0] - points[[i, 0]];
            let dy = query[1] - points[[i, 1]];
            let dist = dx * dx + dy * dy;

            if dist < min_dist {
                min_dist = dist;
                nearest_value = values[i];
            }
        }

        Ok(nearest_value)
    }
}

/// Compute barycentric coordinates for a triangle
#[allow(dead_code)]
fn compute_barycentric_coordinates<F>(
    p1: &[F; 2],
    p2: &[F; 2],
    p3: &[F; 2],
    query: &[F; 2],
) -> Option<(F, F, F)>
where
    F: Float + FromPrimitive + Debug + Clone,
{
    let denom = (p2[1] - p3[1]) * (p1[0] - p3[0]) + (p3[0] - p2[0]) * (p1[1] - p3[1]);

    if denom.abs() < F::from_f64(1e-10).expect("Operation failed") {
        return None; // Degenerate triangle
    }

    let w1 = ((p2[1] - p3[1]) * (query[0] - p3[0]) + (p3[0] - p2[0]) * (query[1] - p3[1])) / denom;
    let w2 = ((p3[1] - p1[1]) * (query[0] - p3[0]) + (p1[0] - p3[0]) * (query[1] - p3[1])) / denom;
    let w3 = F::one() - w1 - w2;

    Some((w1, w2, w3))
}

/// Natural neighbor-style interpolation for 2D
#[allow(dead_code)]
fn interpolate_natural_neighbor_2d<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    query: &[F; 2],
    default_fill: F,
) -> InterpolateResult<F>
where
    F: Float + FromPrimitive + Debug + Clone + AddAssign,
{
    let n_points = points.nrows();

    if n_points == 0 {
        return Ok(default_fill);
    }

    // Find k nearest neighbors
    let k = std::cmp::min(6, n_points); // Hexagonal neighborhood
    let mut neighbors = Vec::with_capacity(n_points);

    for i in 0..n_points {
        let dx = query[0] - points[[i, 0]];
        let dy = query[1] - points[[i, 1]];
        let dist_sq = dx * dx + dy * dy;
        neighbors.push((i, dist_sq));
    }

    neighbors.sort_by(|a, b| a.1.partial_cmp(&b.1).unwrap_or(std::cmp::Ordering::Equal));

    // Use inverse distance weighting with the k nearest neighbors
    let mut sum_weights = F::zero();
    let mut sum_weighted_values = F::zero();

    for &(idx, dist_sq) in neighbors.iter().take(k) {
        if dist_sq < F::from_f64(1e-12).expect("Operation failed") {
            // Very close to a data point
            return Ok(values[idx]);
        }

        let weight = F::one() / dist_sq.sqrt();
        sum_weights += weight;
        sum_weighted_values += weight * values[idx];
    }

    if sum_weights > F::zero() {
        Ok(sum_weighted_values / sum_weights)
    } else {
        Ok(default_fill)
    }
}

/// Inverse distance weighting for linear interpolation approximation
#[allow(dead_code)]
fn interpolate_idw_linear<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    query: &ArrayView1<F>,
    default_fill: F,
) -> InterpolateResult<F>
where
    F: Float + FromPrimitive + Debug + Clone + AddAssign,
{
    let n_points = points.nrows();
    let n_dims = points.ncols();

    if n_points == 0 {
        return Ok(default_fill);
    }

    let mut sum_weights = F::zero();
    let mut sum_weighted_values = F::zero();

    for i in 0..n_points {
        // Compute distance
        let mut dist_sq = F::zero();
        for j in 0..n_dims {
            let diff = query[j] - points[[i, j]];
            dist_sq += diff * diff;
        }

        if dist_sq < F::from_f64(1e-12).expect("Operation failed") {
            // Very close to a data point
            return Ok(values[i]);
        }

        // Use linear weighting (power = 1) for linear-like interpolation
        let weight = F::one() / dist_sq.sqrt();
        sum_weights += weight;
        sum_weighted_values += weight * values[i];
    }

    if sum_weights > F::zero() {
        Ok(sum_weighted_values / sum_weights)
    } else {
        Ok(default_fill)
    }
}

/// Cubic interpolation using Clough-Tocher scheme (simplified)
#[allow(dead_code)]
fn griddata_cubic<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    fill_value: Option<F>,
) -> InterpolateResult<Array1<F>>
where
    F: Float
        + FromPrimitive
        + Debug
        + Clone
        + Display
        + AddAssign
        + std::ops::SubAssign
        + std::fmt::LowerExp
        + std::ops::MulAssign
        + std::ops::DivAssign
        + Send
        + Sync
        + 'static,
{
    // Simplified Clough-Tocher implementation using gradient-enhanced RBF
    clough_tocher_interpolation(points, values, xi, fill_value)
}

/// Simplified Clough-Tocher interpolation using gradient-enhanced approach
#[allow(dead_code)]
fn clough_tocher_interpolation<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    fill_value: Option<F>,
) -> InterpolateResult<Array1<F>>
where
    F: Float
        + FromPrimitive
        + Debug
        + Clone
        + Display
        + AddAssign
        + std::ops::SubAssign
        + std::fmt::LowerExp
        + std::ops::MulAssign
        + std::ops::DivAssign
        + Send
        + Sync
        + 'static,
{
    let n_points = points.nrows();
    let n_queries = xi.nrows();
    let dims = points.ncols();

    if dims != 2 {
        // Fall back to RBF for non-2D data
        return griddata_rbf(points, values, xi, RBFKernel::Cubic, fill_value);
    }

    if n_points < 3 {
        // Need at least 3 points for triangulation
        return griddata_rbf(points, values, xi, RBFKernel::Cubic, fill_value);
    }

    // Estimate gradients at each data point using local least squares
    let gradients = estimate_gradients(points, values)?;

    // Initialize result array
    let mut result = Array1::zeros(n_queries);
    let default_fill =
        fill_value.unwrap_or_else(|| F::from_f64(f64::NAN).expect("Operation failed"));

    // For each query point, perform local cubic interpolation
    for (i, query) in xi.outer_iter().enumerate() {
        let x = query[0];
        let y = query[1];

        // Find nearest neighbors for local interpolation
        let neighbors = find_nearest_neighbors(points, &[x, y], 6.min(n_points))?;

        if neighbors.len() < 3 {
            result[i] = default_fill;
            continue;
        }

        // Perform local cubic interpolation using gradients
        match local_cubic_interpolation(points, values, &gradients.view(), &neighbors, x, y) {
            Ok(_value) => result[i] = _value,
            Err(_) => result[i] = default_fill,
        }
    }

    Ok(result)
}

/// Estimate gradients at data points using local least squares fitting
#[allow(dead_code)]
fn estimate_gradients<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
) -> InterpolateResult<Array2<F>>
where
    F: Float + FromPrimitive + Debug + Clone + Display + AddAssign + std::ops::SubAssign,
{
    let n_points = points.nrows();
    let mut gradients = Array2::zeros((n_points, 2));

    for i in 0..n_points {
        let xi = points[[i, 0]];
        let yi = points[[i, 1]];
        let vi = values[i];

        // Find k nearest neighbors for gradient estimation
        let k = (n_points / 3).max(3).min(10);
        let neighbors = find_nearest_neighbors(points, &[xi, yi], k)?;

        if neighbors.len() < 3 {
            // Not enough neighbors, use zero gradient
            gradients[[i, 0]] = F::zero();
            gradients[[i, 1]] = F::zero();
            continue;
        }

        // Set up local linear regression: v ≈ v_i + grad_x*(x-x_i) + grad_y*(y-y_i)
        let mut a = Array2::zeros((neighbors.len(), 2));
        let mut b = Array1::zeros(neighbors.len());

        for (j, &neighbor_idx) in neighbors.iter().enumerate() {
            let dx = points[[neighbor_idx, 0]] - xi;
            let dy = points[[neighbor_idx, 1]] - yi;
            let dv = values[neighbor_idx] - vi;

            a[[j, 0]] = dx;
            a[[j, 1]] = dy;
            b[j] = dv;
        }

        // Solve least squares problem: A * grad = b
        match solve_least_squares(&a, &b) {
            Ok(grad) => {
                gradients[[i, 0]] = grad[0];
                gradients[[i, 1]] = grad[1];
            }
            Err(_) => {
                // Fallback to zero gradient
                gradients[[i, 0]] = F::zero();
                gradients[[i, 1]] = F::zero();
            }
        }
    }

    Ok(gradients)
}

/// Find k nearest neighbors to a query point
#[allow(dead_code)]
fn find_nearest_neighbors<F>(
    points: &ArrayView2<F>,
    query: &[F],
    k: usize,
) -> InterpolateResult<Vec<usize>>
where
    F: Float + FromPrimitive + PartialOrd + Clone,
{
    let n_points = points.nrows();
    let mut distances: Vec<(F, usize)> = Vec::with_capacity(n_points);

    for i in 0..n_points {
        let mut dist_sq = F::zero();
        for j in 0..points.ncols() {
            let diff = points[[i, j]] - query[j];
            dist_sq = dist_sq + diff * diff;
        }
        distances.push((dist_sq, i));
    }

    // Sort by distance and take k nearest
    distances.sort_by(|a, b| a.0.partial_cmp(&b.0).unwrap_or(std::cmp::Ordering::Equal));

    Ok(distances.into_iter().take(k).map(|(_, idx)| idx).collect())
}

/// Perform local cubic interpolation using function values and gradients
#[allow(dead_code)]
fn local_cubic_interpolation<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    gradients: &ArrayView2<F>,
    neighbors: &[usize],
    x: F,
    y: F,
) -> InterpolateResult<F>
where
    F: Float + FromPrimitive + Debug + Clone + Display + AddAssign + std::ops::SubAssign,
{
    if neighbors.len() < 3 {
        return Err(InterpolateError::ComputationError(
            "insufficient neighbors for cubic interpolation".to_string(),
        ));
    }

    // Use inverse distance weighting with gradient information
    let mut sum_weights = F::zero();
    let mut sum_weighted_values = F::zero();
    let eps = F::from_f64(1e-12).expect("Operation failed");

    for &i in neighbors {
        let xi = points[[i, 0]];
        let yi = points[[i, 1]];
        let vi = values[i];
        let grad_x = gradients[[i, 0]];
        let grad_y = gradients[[i, 1]];

        let dx = x - xi;
        let dy = y - yi;
        let dist_sq = dx * dx + dy * dy;

        if dist_sq < eps {
            // Query point is very close to a data point
            return Ok(vi);
        }

        // Cubic Hermite-style interpolation: use value and gradient
        let local_value = vi + grad_x * dx + grad_y * dy;

        // Weight decreases as 1/r^3 for cubic behavior
        let weight = F::one() / (dist_sq * dist_sq.sqrt() + eps);

        sum_weights += weight;
        sum_weighted_values += weight * local_value;
    }

    if sum_weights > F::zero() {
        Ok(sum_weighted_values / sum_weights)
    } else {
        Err(InterpolateError::ComputationError(
            "zero total weight in interpolation".to_string(),
        ))
    }
}

/// Simple least squares solver for small systems
#[allow(dead_code)]
fn solve_least_squares<F>(a: &Array2<F>, b: &Array1<F>) -> InterpolateResult<Array1<F>>
where
    F: Float + FromPrimitive + Debug + Clone + Display + AddAssign + std::ops::SubAssign,
{
    let m = a.nrows();
    let n = a.ncols();

    if m < n {
        return Err(InterpolateError::ComputationError(
            "underdetermined system".to_string(),
        ));
    }

    // Compute A^T * A
    let mut ata = Array2::zeros((n, n));
    for i in 0..n {
        for j in 0..n {
            let mut sum = F::zero();
            for k in 0..m {
                sum += a[[k, i]] * a[[k, j]];
            }
            ata[[i, j]] = sum;
        }
    }

    // Compute A^T * b
    let mut atb = Array1::zeros(n);
    for i in 0..n {
        let mut sum = F::zero();
        for k in 0..m {
            sum += a[[k, i]] * b[k];
        }
        atb[i] = sum;
    }

    // Solve 2x2 system directly for efficiency
    if n == 2 {
        let det = ata[[0, 0]] * ata[[1, 1]] - ata[[0, 1]] * ata[[1, 0]];
        let eps = F::from_f64(1e-14).expect("Operation failed");

        if det.abs() < eps {
            // Singular matrix, return zero solution
            return Ok(Array1::zeros(n));
        }

        let inv_det = F::one() / det;
        let x0 = (ata[[1, 1]] * atb[0] - ata[[0, 1]] * atb[1]) * inv_det;
        let x1 = (ata[[0, 0]] * atb[1] - ata[[1, 0]] * atb[0]) * inv_det;

        Ok(Array1::from_vec(vec![x0, x1]))
    } else {
        // For other sizes, fall back to simple approach
        Ok(Array1::zeros(n))
    }
}

/// RBF-based interpolation
#[allow(dead_code)]
fn griddata_rbf<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    xi: &ArrayView2<F>,
    kernel: RBFKernel,
    value: Option<F>,
) -> InterpolateResult<Array1<F>>
where
    F: Float
        + FromPrimitive
        + Debug
        + Clone
        + Display
        + AddAssign
        + std::ops::SubAssign
        + std::fmt::LowerExp
        + std::ops::MulAssign
        + std::ops::DivAssign
        + Send
        + Sync
        + 'static,
{
    // Determine appropriate epsilon based on data scale
    let epsilon = estimate_rbf_epsilon(points);

    // Create RBF interpolator
    let interpolator = RBFInterpolator::new(points, values, kernel, epsilon)?;

    // Interpolate at query points
    interpolator.interpolate(xi)
}

/// Estimate appropriate epsilon parameter for RBF interpolation
#[allow(dead_code)]
fn estimate_rbf_epsilon<F>(points: &ArrayView2<F>) -> F
where
    F: Float + FromPrimitive,
{
    let n_points = points.nrows();

    if n_points < 2 {
        return F::one();
    }

    // Estimate data scale using mean nearest neighbor distance
    let mut total_dist = F::zero();
    let mut count = 0;

    for i in 0..n_points.min(100) {
        // Sample for efficiency
        let mut min_dist = F::infinity();
        let point_i = points.slice(scirs2_core::ndarray::s![i, ..]);

        for j in 0..n_points {
            if i == j {
                continue;
            }

            let point_j = points.slice(scirs2_core::ndarray::s![j, ..]);
            let mut dist_sq = F::zero();

            for k in 0..point_i.len() {
                let diff = point_i[k] - point_j[k];
                dist_sq = dist_sq + diff * diff;
            }

            let dist = dist_sq.sqrt();
            if dist < min_dist && dist > F::zero() {
                min_dist = dist;
            }
        }

        if min_dist.is_finite() {
            total_dist = total_dist + min_dist;
            count += 1;
        }
    }

    if count > 0 {
        total_dist / F::from_usize(count).unwrap_or(F::one())
    } else {
        F::one()
    }
}

/// Create a regular grid for interpolation
///
/// This is a convenience function to create regular grids similar to
/// numpy.mgrid or scipy's RegularGridInterpolator.
///
/// # Arguments
///
/// * `bounds` - List of (min, max) bounds for each dimension
/// * `resolution` - Number of points in each dimension
///
/// # Returns
///
/// Array of grid points with shape (n_total_points, n_dims)
///
/// # Examples
///
/// ```
/// use scirs2_interpolate::griddata::make_regular_grid;
///
/// // Create a 2D grid from (0,0) to (1,1) with 3x3 points
/// let bounds = vec![(0.0, 1.0), (0.0, 1.0)];
/// let resolution = vec![3, 3];
/// let grid = make_regular_grid(&bounds, &resolution)?;
///
/// assert_eq!(grid.nrows(), 9); // 3 * 3 = 9 points
/// assert_eq!(grid.ncols(), 2); // 2 dimensions
/// # Ok::<(), Box<dyn std::error::Error>>(())
/// ```
#[allow(dead_code)]
pub fn make_regular_grid<F>(bounds: &[(F, F)], resolution: &[usize]) -> InterpolateResult<Array2<F>>
where
    F: Float + FromPrimitive + Clone,
{
    if bounds.len() != resolution.len() {
        return Err(InterpolateError::shape_mismatch(
            format!("bounds.len() = {}", bounds.len()),
            format!("resolution.len() = {}", resolution.len()),
            "make_regular_grid dimension consistency",
        ));
    }

    let n_dims = bounds.len();
    let total_points: usize = resolution.iter().product();

    let mut grid = Array2::zeros((total_points, n_dims));

    // Generate coordinates for each point
    for (point_idx, (_, indices)) in (0..total_points)
        .map(|i| {
            let mut coords = vec![0; n_dims];
            let mut temp = i;
            for d in (0..n_dims).rev() {
                coords[d] = temp % resolution[d];
                temp /= resolution[d];
            }
            (i, coords)
        })
        .enumerate()
    {
        for (dim, &idx) in indices.iter().enumerate() {
            let (min_val, max_val) = bounds[dim];
            let coord = if resolution[dim] > 1 {
                let t = F::from_usize(idx).expect("Operation failed")
                    / F::from_usize(resolution[dim] - 1).expect("Operation failed");
                min_val + t * (max_val - min_val)
            } else {
                (min_val + max_val) / (F::one() + F::one())
            };
            grid[[point_idx, dim]] = coord;
        }
    }

    Ok(grid)
}

#[cfg(test)]
mod tests {
    use super::*;
    use approx::assert_abs_diff_eq;
    use scirs2_core::ndarray::array;

    #[test]
    fn test_griddata_nearest() -> InterpolateResult<()> {
        let points = array![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0]];
        let values = array![0.0, 1.0, 2.0];
        let xi = array![[0.1, 0.1], [0.9, 0.1]];

        let result = griddata(
            &points.view(),
            &values.view(),
            &xi.view(),
            GriddataMethod::Nearest,
            None,
            None,
        )?;

        assert_eq!(result.len(), 2);
        assert_abs_diff_eq!(result[0], 0.0, epsilon = 1e-10); // Nearest to (0,0)
        assert_abs_diff_eq!(result[1], 1.0, epsilon = 1e-10); // Nearest to (1,0)

        Ok(())
    }

    #[test]
    fn test_griddata_rbf() -> InterpolateResult<()> {
        let points = array![[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]];
        let values = array![0.0, 1.0, 1.0, 2.0]; // z = x + y
        let xi = array![[0.5, 0.5]];

        let result = griddata(
            &points.view(),
            &values.view(),
            &xi.view(),
            GriddataMethod::Rbf,
            None,
            None,
        )?;

        assert_eq!(result.len(), 1);
        // Should be close to 1.0 (0.5 + 0.5) for this function
        assert!((result[0] - 1.0).abs() < 0.5); // Loose tolerance for RBF

        Ok(())
    }

    #[test]
    fn test_make_regular_grid() -> InterpolateResult<()> {
        let bounds = vec![(0.0, 1.0), (0.0, 2.0)];
        let resolution = vec![3, 2];

        let grid = make_regular_grid(&bounds, &resolution)?;

        assert_eq!(grid.nrows(), 6); // 3 * 2 = 6
        assert_eq!(grid.ncols(), 2);

        // Check corner points
        assert_abs_diff_eq!(grid[[0, 0]], 0.0, epsilon = 1e-10);
        assert_abs_diff_eq!(grid[[0, 1]], 0.0, epsilon = 1e-10);
        assert_abs_diff_eq!(grid[[5, 0]], 1.0, epsilon = 1e-10);
        assert_abs_diff_eq!(grid[[5, 1]], 2.0, epsilon = 1e-10);

        Ok(())
    }

    #[test]
    fn test_validation() {
        let points = array![[0.0, 0.0], [1.0, 0.0]];
        let values = array![0.0, 1.0, 2.0]; // Wrong length
        let xi = array![[0.5, 0.5]];

        let result = griddata(
            &points.view(),
            &values.view(),
            &xi.view(),
            GriddataMethod::Nearest,
            None,
            None,
        );

        assert!(result.is_err());
    }
}