scirs2-interpolate 0.4.1

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
//! Grid transformation and resampling utilities
//!
//! This module provides functions for transforming irregular data onto regular grids
//! and vice versa, as well as various grid manipulation utilities.

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

/// Grid transformation methods
#[derive(Debug, Clone, Copy, PartialEq)]
pub enum GridTransformMethod {
    /// Nearest neighbor assignment
    Nearest,
    /// Linear interpolation
    Linear,
    /// Cubic spline interpolation
    Cubic,
}

/// Create a grid of evenly spaced coordinates
///
/// # Arguments
///
/// * `bounds` - Min and max bounds for each dimension
/// * `shape` - Number of points in each dimension
///
/// # Returns
///
/// Vector of coordinate arrays, one for each dimension
///
/// # Examples
///
/// ```
/// use scirs2_interpolate::grid::create_regular_grid;
///
/// // Create a 10×20 grid from (0,0) to (1,2)
/// let grid_coords = create_regular_grid(
///     &[(0.0, 1.0), (0.0, 2.0)],
///     &[10, 20]
/// ).expect("Operation failed");
///
/// assert_eq!(grid_coords.len(), 2);
/// assert_eq!(grid_coords[0].len(), 10);
/// assert_eq!(grid_coords[1].len(), 20);
/// ```
#[allow(dead_code)]
pub fn create_regular_grid<F>(
    bounds: &[(F, F)],
    shape: &[usize],
) -> InterpolateResult<Vec<Array1<F>>>
where
    F: Float + FromPrimitive + Debug,
{
    if bounds.len() != shape.len() {
        return Err(InterpolateError::invalid_input(
            "bounds and shape must have the same length".to_string(),
        ));
    }

    let n_dims = bounds.len();
    let mut grid_coords = Vec::with_capacity(n_dims);

    for i in 0..n_dims {
        let (min, max) = bounds[i];
        let n_points = shape[i];

        if min >= max {
            return Err(InterpolateError::invalid_input(
                "min bound must be less than max bound".to_string(),
            ));
        }

        if n_points < 2 {
            return Err(InterpolateError::invalid_input(
                "grid shape must have at least 2 points in each dimension".to_string(),
            ));
        }

        let mut coords = Array1::zeros(n_points);

        if n_points == 1 {
            coords[0] = min;
        } else {
            let step = (max - min) / F::from_usize(n_points - 1).expect("Operation failed");
            for j in 0..n_points {
                coords[j] = min + F::from_usize(j).expect("Operation failed") * step;
            }
        }

        grid_coords.push(coords);
    }

    Ok(grid_coords)
}

/// Resample irregular (scattered) data onto a regular grid
///
/// # Arguments
///
/// * `points` - Coordinates of scattered data points (n_points × n_dimensions)
/// * `values` - Values at the scattered data points (n_points)
/// * `gridshape` - Shape of the output grid (number of points in each dimension)
/// * `grid_bounds` - Min and max bounds for each dimension of the grid
/// * `method` - Interpolation method to use
/// * `fill_value` - Value to use for grid points outside the convex hull of input points
///
/// # Returns
///
/// A tuple containing:
/// * The grid coordinates for each dimension (vector of arrays)
/// * The resampled values on the regular grid
#[allow(dead_code)]
pub fn resample_to_grid<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    gridshape: &[usize],
    grid_bounds: &[(F, F)],
    method: GridTransformMethod,
    fill_value: F,
) -> InterpolateResult<(Vec<Array1<F>>, ArrayD<F>)>
where
    F: Float
        + FromPrimitive
        + Debug
        + Clone
        + PartialOrd
        + Zero
        + 'static
        + std::fmt::Display
        + std::ops::AddAssign
        + std::ops::SubAssign
        + std::ops::MulAssign
        + std::ops::DivAssign
        + std::fmt::LowerExp
        + Send
        + Sync,
{
    if points.nrows() != values.len() {
        return Err(InterpolateError::invalid_input(
            "Number of points and values must match".to_string(),
        ));
    }

    if points.ncols() != grid_bounds.len() {
        return Err(InterpolateError::invalid_input(
            "Point dimensions must match grid _bounds dimensions".to_string(),
        ));
    }

    if grid_bounds.len() != gridshape.len() {
        return Err(InterpolateError::invalid_input(
            "Grid _bounds and shape dimensions must match".to_string(),
        ));
    }

    // Create the regular grid
    let grid_coords = create_regular_grid(grid_bounds, gridshape)?;

    // Create a multidimensional array with the specified shape
    let shape: Vec<usize> = gridshape.to_vec();
    let mut grid_values = ArrayD::from_elem(shape.clone(), fill_value);

    match method {
        GridTransformMethod::Nearest => {
            resample_nearest_neighbor(points, values, &grid_coords, &mut grid_values, fill_value)?;
        }
        GridTransformMethod::Linear => {
            resample_linear(points, values, &grid_coords, &mut grid_values, fill_value)?;
        }
        GridTransformMethod::Cubic => {
            resample_rbf(points, values, &grid_coords, &mut grid_values, fill_value)?;
        }
    }

    Ok((grid_coords, grid_values))
}

/// Resample using nearest neighbor interpolation
#[allow(dead_code)]
fn resample_nearest_neighbor<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    grid_coords: &[Array1<F>],
    grid_values: &mut ArrayD<F>,
    fill_value: F,
) -> InterpolateResult<()>
where
    F: Float + FromPrimitive + Debug + Clone + PartialOrd + Zero,
{
    let n_dims = grid_coords.len();

    // For each grid point, find the nearest data point
    let gridshape: Vec<usize> = grid_coords.iter().map(|coord| coord.len()).collect();

    // Generate all grid point coordinates
    let mut indices = vec![0; n_dims];

    loop {
        // Convert indices to actual coordinates
        let mut grid_point = vec![F::zero(); n_dims];
        for (dim, &idx) in indices.iter().enumerate() {
            grid_point[dim] = grid_coords[dim][idx];
        }

        // Find nearest data point
        let mut min_dist_sq = F::infinity();
        let mut nearest_value = fill_value;

        for i in 0..points.nrows() {
            let mut dist_sq = F::zero();
            for j in 0..n_dims {
                let diff = points[[i, j]] - grid_point[j];
                dist_sq = dist_sq + diff * diff;
            }

            if dist_sq < min_dist_sq {
                min_dist_sq = dist_sq;
                nearest_value = values[i];
            }
        }

        // Set the grid _value
        grid_values[&indices[..]] = nearest_value;

        // Increment indices
        if !increment_indices(&mut indices, &gridshape) {
            break;
        }
    }

    Ok(())
}

/// Resample using RBF interpolation for smooth results
#[allow(dead_code)]
fn resample_rbf<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    grid_coords: &[Array1<F>],
    grid_values: &mut ArrayD<F>,
    fill_value: F,
) -> InterpolateResult<()>
where
    F: Float
        + FromPrimitive
        + Debug
        + Clone
        + PartialOrd
        + Zero
        + 'static
        + std::fmt::Display
        + std::ops::AddAssign
        + std::ops::SubAssign
        + std::ops::MulAssign
        + std::ops::DivAssign
        + std::fmt::LowerExp
        + Send
        + Sync,
{
    // Create RBF interpolator
    let rbf = RBFInterpolator::new(
        points,
        values,
        RBFKernel::Gaussian,
        F::from_f64(1.0).unwrap_or_else(|| F::one()),
    )?;

    let n_dims = grid_coords.len();
    let gridshape: Vec<usize> = grid_coords.iter().map(|coord| coord.len()).collect();
    let mut indices = vec![0; n_dims];

    loop {
        // Convert indices to actual coordinates
        let mut grid_point = Array1::zeros(n_dims);
        for (dim, &idx) in indices.iter().enumerate() {
            grid_point[dim] = grid_coords[dim][idx];
        }

        // Evaluate RBF at this grid point
        let interp_value = match rbf.interpolate(&grid_point.view().insert_axis(Axis(0))) {
            Ok(val) => val[0],
            Err(_) => fill_value,
        };

        grid_values[&indices[..]] = interp_value;

        // Increment indices
        if !increment_indices(&mut indices, &gridshape) {
            break;
        }
    }

    Ok(())
}

/// Linear interpolation for grid resampling (simplified implementation)
#[allow(dead_code)]
fn resample_linear<F>(
    points: &ArrayView2<F>,
    values: &ArrayView1<F>,
    grid_coords: &[Array1<F>],
    grid_values: &mut ArrayD<F>,
    fill_value: F,
) -> InterpolateResult<()>
where
    F: Float + FromPrimitive + Debug + Clone + PartialOrd + Zero,
{
    // For simplicity, fall back to nearest neighbor for multidimensional case
    // A full implementation would use multilinear interpolation
    resample_nearest_neighbor(points, values, grid_coords, grid_values, fill_value)
}

/// Helper function to increment multi-dimensional indices
#[allow(dead_code)]
fn increment_indices(indices: &mut [usize], shape: &[usize]) -> bool {
    for i in (0..indices.len()).rev() {
        indices[i] += 1;
        if indices[i] < shape[i] {
            return true;
        }
        indices[i] = 0;
    }
    false
}

/// Resample a regular grid to another regular grid with different resolution or bounds
///
/// # Arguments
///
/// * `src_coords` - Source grid coordinates (vector of arrays, one per dimension)
/// * `src_values` - Source grid values
/// * `dst_coords` - Destination grid coordinates (vector of arrays, one per dimension)
/// * `method` - Interpolation method to use
/// * `fill_value` - Value to use for grid points outside the source grid
///
/// # Returns
///
/// Resampled values on the destination grid
#[allow(dead_code)]
pub fn resample_grid_to_grid<F, D>(
    src_coords: &[Array1<F>],
    src_values: &scirs2_core::ndarray::ArrayView<F, D>,
    dst_coords: &[Array1<F>],
    method: GridTransformMethod,
    fill_value: F,
) -> InterpolateResult<ArrayD<F>>
where
    F: Float + FromPrimitive + Debug + Clone + PartialOrd + Zero + 'static,
    D: scirs2_core::ndarray::Dimension,
{
    if src_coords.len() != dst_coords.len() {
        return Err(InterpolateError::invalid_input(
            "Source and destination must have same number of dimensions".to_string(),
        ));
    }

    let _n_dims = src_coords.len(); // Reserved for future use

    // Verify source coordinates match source _values shape
    for (i, coord) in src_coords.iter().enumerate() {
        if coord.len() != src_values.shape()[i] {
            return Err(InterpolateError::invalid_input(format!(
                "Source coordinate dimension {} length doesn't match _values shape",
                i
            )));
        }
    }

    // Create destination grid shape
    let dstshape: Vec<usize> = dst_coords.iter().map(|coord| coord.len()).collect();
    let mut dst_values = ArrayD::from_elem(dstshape.clone(), fill_value);

    match method {
        GridTransformMethod::Nearest => {
            grid_to_grid_nearest(
                src_coords,
                src_values,
                dst_coords,
                &mut dst_values,
                fill_value,
            )?;
        }
        GridTransformMethod::Linear => {
            grid_to_grid_linear(
                src_coords,
                src_values,
                dst_coords,
                &mut dst_values,
                fill_value,
            )?;
        }
        GridTransformMethod::Cubic => {
            // For cubic, we'll use linear interpolation as it's more stable for grids
            grid_to_grid_linear(
                src_coords,
                src_values,
                dst_coords,
                &mut dst_values,
                fill_value,
            )?;
        }
    }

    Ok(dst_values)
}

/// Convert multi-dimensional indices to linear index
#[allow(dead_code)]
fn ravel_multi_index(indices: &[usize], shape: &[usize]) -> usize {
    let mut linear_idx = 0;
    let mut stride = 1;

    for i in (0..indices.len()).rev() {
        linear_idx += indices[i] * stride;
        stride *= shape[i];
    }

    linear_idx
}

/// Grid-to-grid resampling using nearest neighbor
#[allow(dead_code)]
fn grid_to_grid_nearest<F, D>(
    src_coords: &[Array1<F>],
    src_values: &scirs2_core::ndarray::ArrayView<F, D>,
    dst_coords: &[Array1<F>],
    dst_values: &mut ArrayD<F>,
    fill_value: F,
) -> InterpolateResult<()>
where
    F: Float + FromPrimitive + Debug + Clone + PartialOrd + Zero,
    D: scirs2_core::ndarray::Dimension,
{
    let n_dims = src_coords.len();
    let dstshape: Vec<usize> = dst_coords.iter().map(|coord| coord.len()).collect();
    let mut indices = vec![0; n_dims];

    loop {
        // Get destination grid point coordinates
        let mut dst_point = vec![F::zero(); n_dims];
        for (dim, &idx) in indices.iter().enumerate() {
            dst_point[dim] = dst_coords[dim][idx];
        }

        // Find nearest source grid indices
        let mut src_indices = vec![0; n_dims];
        let valid = true;

        for dim in 0..n_dims {
            let coord = dst_point[dim];
            let src_coord = &src_coords[dim];

            // Find nearest index in source coordinates
            let mut best_idx = 0;
            let mut min_dist = (src_coord[0] - coord).abs();

            for (i, &src_val) in src_coord.iter().enumerate() {
                let dist = (src_val - coord).abs();
                if dist < min_dist {
                    min_dist = dist;
                    best_idx = i;
                }
            }

            src_indices[dim] = best_idx;
        }

        // Get the source _value and assign to destination
        if valid {
            // Use linear indexing for all cases to avoid generic dimension issues
            let linear_idx = ravel_multi_index(&src_indices, src_values.shape());
            let src_value = src_values.as_slice().expect("Operation failed")[linear_idx];
            let dst_linear_idx = ravel_multi_index(&indices, &dstshape);
            dst_values.as_slice_mut().expect("Operation failed")[dst_linear_idx] = src_value;
        } else {
            let dst_linear_idx = ravel_multi_index(&indices, &dstshape);
            dst_values.as_slice_mut().expect("Operation failed")[dst_linear_idx] = fill_value;
        }

        // Increment indices
        if !increment_indices(&mut indices, &dstshape) {
            break;
        }
    }

    Ok(())
}

/// Grid-to-grid resampling using linear interpolation
#[allow(dead_code)]
fn grid_to_grid_linear<F, D>(
    src_coords: &[Array1<F>],
    src_values: &scirs2_core::ndarray::ArrayView<F, D>,
    dst_coords: &[Array1<F>],
    dst_values: &mut ArrayD<F>,
    fill_value: F,
) -> InterpolateResult<()>
where
    F: Float + FromPrimitive + Debug + Clone + PartialOrd + Zero,
    D: scirs2_core::ndarray::Dimension,
{
    let n_dims = src_coords.len();
    let dstshape: Vec<usize> = dst_coords.iter().map(|coord| coord.len()).collect();
    let mut indices = vec![0; n_dims];

    loop {
        // Get destination grid point coordinates
        let mut dst_point = vec![F::zero(); n_dims];
        for (dim, &idx) in indices.iter().enumerate() {
            dst_point[dim] = dst_coords[dim][idx];
        }

        // Perform multilinear interpolation
        let interpolated_value =
            multilinear_interpolate(src_coords, src_values, &dst_point, fill_value)?;

        dst_values[&indices[..]] = interpolated_value;

        // Increment indices
        if !increment_indices(&mut indices, &dstshape) {
            break;
        }
    }

    Ok(())
}

/// Perform multilinear interpolation at a single point
#[allow(dead_code)]
fn multilinear_interpolate<F, D>(
    coords: &[Array1<F>],
    values: &scirs2_core::ndarray::ArrayView<F, D>,
    point: &[F],
    fill_value: F,
) -> InterpolateResult<F>
where
    F: Float + FromPrimitive + Debug + Clone + PartialOrd + Zero,
    D: scirs2_core::ndarray::Dimension,
{
    let n_dims = coords.len();

    // Find bounding grid cells for each dimension
    let mut lower_indices = vec![0; n_dims];
    let mut upper_indices = vec![0; n_dims];
    let mut weights = vec![F::zero(); n_dims];

    for dim in 0..n_dims {
        let coord_array = &coords[dim];
        let target = point[dim];

        // Find the interval containing the target
        let mut found = false;
        for i in 0..coord_array.len() - 1 {
            if target >= coord_array[i] && target <= coord_array[i + 1] {
                lower_indices[dim] = i;
                upper_indices[dim] = i + 1;

                // Calculate interpolation weight
                let dx = coord_array[i + 1] - coord_array[i];
                if dx.abs() > F::zero() {
                    weights[dim] = (target - coord_array[i]) / dx;
                } else {
                    weights[dim] = F::zero();
                }
                found = true;
                break;
            }
        }

        if !found {
            // Point is outside grid bounds
            return Ok(fill_value);
        }
    }

    // Perform multilinear interpolation
    // For N dimensions, we need 2^N corner values
    let n_corners = 1 << n_dims; // 2^n_dims
    let mut result = F::zero();

    for corner in 0..n_corners {
        let mut corner_indices = vec![0; n_dims];
        let mut corner_weight = F::one();

        for dim in 0..n_dims {
            if (corner >> dim) & 1 == 0 {
                corner_indices[dim] = lower_indices[dim];
                corner_weight = corner_weight * (F::one() - weights[dim]);
            } else {
                corner_indices[dim] = upper_indices[dim];
                corner_weight = corner_weight * weights[dim];
            }
        }

        // Get _value at this corner using linear indexing
        let linear_idx = ravel_multi_index(&corner_indices, values.shape());
        let corner_value = values.as_slice().expect("Operation failed")[linear_idx];
        result = result + corner_weight * corner_value;
    }

    Ok(result)
}

/// Maps values from a regular grid to arbitrary points using interpolation
///
/// # Arguments
///
/// * `grid_coords` - Grid coordinates (vector of arrays, one per dimension)
/// * `grid_values` - Values on the regular grid
/// * `query_points` - Points at which to evaluate (n_points × n_dimensions)
/// * `method` - Interpolation method to use
/// * `fill_value` - Value to use for query points outside the grid
///
/// # Returns
///
/// Values at the query points
#[allow(dead_code)]
pub fn map_grid_to_points<F, D>(
    grid_coords: &[Array1<F>],
    grid_values: &scirs2_core::ndarray::ArrayView<F, D>,
    query_points: &ArrayView2<F>,
    method: GridTransformMethod,
    fill_value: F,
) -> InterpolateResult<Array1<F>>
where
    F: Float + FromPrimitive + Debug + Clone + PartialOrd + Zero + 'static,
    D: scirs2_core::ndarray::Dimension,
{
    let n_points = query_points.nrows();
    let n_dims = query_points.ncols();

    if grid_coords.len() != n_dims {
        return Err(InterpolateError::invalid_input(
            "Grid coordinates and query point dimensions must match".to_string(),
        ));
    }

    // Verify grid coordinates match grid _values shape
    for (i, coord) in grid_coords.iter().enumerate() {
        if coord.len() != grid_values.shape()[i] {
            return Err(InterpolateError::invalid_input(format!(
                "Grid coordinate dimension {} length doesn't match _values shape",
                i
            )));
        }
    }

    let mut result = Array1::zeros(n_points);

    for i in 0..n_points {
        let query_point: Vec<F> = query_points.row(i).to_vec();

        let interpolated_value = match method {
            GridTransformMethod::Nearest => {
                grid_nearest_neighbor(grid_coords, grid_values, &query_point, fill_value)?
            }
            GridTransformMethod::Linear => {
                multilinear_interpolate(grid_coords, grid_values, &query_point, fill_value)?
            }
            GridTransformMethod::Cubic => {
                // For cubic, we fall back to linear for stability
                multilinear_interpolate(grid_coords, grid_values, &query_point, fill_value)?
            }
        };

        result[i] = interpolated_value;
    }

    Ok(result)
}

/// Find the nearest grid point value
#[allow(dead_code)]
fn grid_nearest_neighbor<F, D>(
    grid_coords: &[Array1<F>],
    grid_values: &scirs2_core::ndarray::ArrayView<F, D>,
    query_point: &[F],
    _fill_value: F,
) -> InterpolateResult<F>
where
    F: Float + FromPrimitive + Debug + Clone + PartialOrd + Zero,
    D: scirs2_core::ndarray::Dimension,
{
    let n_dims = grid_coords.len();
    let mut nearest_indices = vec![0; n_dims];

    for dim in 0..n_dims {
        let coord_array = &grid_coords[dim];
        let target = query_point[dim];

        // Find nearest index
        let mut best_idx = 0;
        let mut min_dist = (coord_array[0] - target).abs();

        for (i, &coord_val) in coord_array.iter().enumerate() {
            let dist = (coord_val - target).abs();
            if dist < min_dist {
                min_dist = dist;
                best_idx = i;
            }
        }

        nearest_indices[dim] = best_idx;
    }

    let linear_idx = ravel_multi_index(&nearest_indices, grid_values.shape());
    Ok(grid_values.as_slice().expect("Operation failed")[linear_idx])
}

/// Efficient grid coordinate range checking
#[allow(dead_code)]
fn point_in_grid_bounds<F>(_gridcoords: &[Array1<F>], point: &[F]) -> bool
where
    F: Float + PartialOrd,
{
    for (dim, coord_array) in _gridcoords.iter().enumerate() {
        let target = point[dim];
        let min_coord = coord_array[0];
        let max_coord = coord_array[coord_array.len() - 1];

        if target < min_coord || target > max_coord {
            return false;
        }
    }
    true
}

/// Create a tensor product grid from coordinate arrays
///
/// This function creates all combinations of coordinates from the input arrays,
/// useful for creating meshgrids for evaluation.
#[allow(dead_code)]
pub fn create_meshgrid<F>(coords: &[Array1<F>]) -> InterpolateResult<Array2<F>>
where
    F: Float + FromPrimitive + Debug + Clone + Zero,
{
    let n_dims = coords.len();
    if n_dims == 0 {
        return Err(InterpolateError::invalid_input(
            "At least one coordinate array required".to_string(),
        ));
    }

    // Calculate total number of grid points
    let mut total_points = 1;
    for coord in coords {
        total_points *= coord.len();
    }

    let mut result = Array2::zeros((total_points, n_dims));
    let shapes: Vec<usize> = coords.iter().map(|c| c.len()).collect();
    let mut indices = vec![0; n_dims];

    for row in 0..total_points {
        // Set coordinates for this grid point
        for (dim, &idx) in indices.iter().enumerate() {
            result[[row, dim]] = coords[dim][idx];
        }

        // Increment multi-dimensional indices
        increment_indices(&mut indices, &shapes);
    }

    Ok(result)
}

/// Calculate grid spacing for each dimension
#[allow(dead_code)]
pub fn calculate_grid_spacing<F>(coords: &[Array1<F>]) -> InterpolateResult<Vec<F>>
where
    F: Float + FromPrimitive + Debug + Clone,
{
    let mut spacings = Vec::with_capacity(coords.len());

    for coord in coords {
        if coord.len() < 2 {
            return Err(InterpolateError::invalid_input(
                "Grid coordinates must have at least 2 points".to_string(),
            ));
        }

        // Calculate average spacing (assumes roughly uniform grid)
        let total_range = coord[coord.len() - 1] - coord[0];
        let n_intervals = F::from_usize(coord.len() - 1).expect("Operation failed");
        let avg_spacing = total_range / n_intervals;

        spacings.push(avg_spacing);
    }

    Ok(spacings)
}

#[cfg(test)]
mod tests {
    use super::*;
    use approx::assert_abs_diff_eq;
    // テスト関数で使用するモジュール

    #[test]
    fn test_create_regular_grid() {
        // 1D grid
        let grid_1d = create_regular_grid(&[(0.0, 1.0)], &[5]).expect("Operation failed");

        assert_eq!(grid_1d.len(), 1);
        assert_eq!(grid_1d[0].len(), 5);
        assert_abs_diff_eq!(grid_1d[0][0], 0.0, epsilon = 1e-10);
        assert_abs_diff_eq!(grid_1d[0][4], 1.0, epsilon = 1e-10);

        // 2D grid
        let grid_2d =
            create_regular_grid(&[(0.0, 1.0), (-1.0, 1.0)], &[3, 5]).expect("Operation failed");

        assert_eq!(grid_2d.len(), 2);
        assert_eq!(grid_2d[0].len(), 3);
        assert_eq!(grid_2d[1].len(), 5);
        assert_abs_diff_eq!(grid_2d[0][0], 0.0, epsilon = 1e-10);
        assert_abs_diff_eq!(grid_2d[0][2], 1.0, epsilon = 1e-10);
        assert_abs_diff_eq!(grid_2d[1][0], -1.0, epsilon = 1e-10);
        assert_abs_diff_eq!(grid_2d[1][4], 1.0, epsilon = 1e-10);
    }
}