Skip to main content

scirs2_interpolate/resampling/
mod.rs

1//! Grid resampling and extrapolation with configurable boundary modes.
2//!
3//! Provides:
4//! - 1-D resampling with linear, cubic spline, nearest-neighbour, and Lanczos methods
5//! - 2-D separable resampling on a regular grid
6//! - Scattered-to-grid conversion using inverse distance weighting (IDW)
7//! - Symbolic spline derivative (degree reduction)
8//! - Multiple extrapolation modes (Nearest, Linear, Polynomial, Reflection, Periodic, Zero, Constant)
9
10use crate::error::InterpolateError;
11
12// ─── ExtrapolationMode ───────────────────────────────────────────────────────
13
14/// How to handle queries outside the data domain.
15#[non_exhaustive]
16#[derive(Debug, Clone, PartialEq)]
17pub enum ExtrapolationMode {
18    /// Clamp to the nearest boundary value.
19    Nearest,
20    /// Linearly extrapolate using the slope at the boundary.
21    Linear,
22    /// Extrapolate with a polynomial of given degree fit to the last `degree+1` points.
23    Polynomial(usize),
24    /// Mirror / reflect the index about the boundary.
25    Reflection,
26    /// Wrap the index periodically.
27    Periodic,
28    /// Return zero outside the domain.
29    Zero,
30    /// Return a fixed constant value outside the domain.
31    Constant(f64),
32}
33
34// ─── ResamplingMethod ────────────────────────────────────────────────────────
35
36/// Interpolation method to use within the data domain.
37#[non_exhaustive]
38#[derive(Debug, Clone, PartialEq)]
39pub enum ResamplingMethod {
40    /// Bi/trilinear interpolation.
41    Linear,
42    /// Natural cubic spline.
43    CubicSpline,
44    /// Nearest-neighbour.
45    Nearest,
46    /// Lanczos windowed-sinc with the given number of lobes (a).
47    Lanczos(usize),
48}
49
50// ─── ResamplingConfig ────────────────────────────────────────────────────────
51
52/// Configuration for resampling operations.
53#[derive(Debug, Clone)]
54pub struct ResamplingConfig {
55    /// In-domain interpolation method.
56    pub method: ResamplingMethod,
57    /// Out-of-domain extrapolation strategy.
58    pub extrapolation: ExtrapolationMode,
59}
60
61impl Default for ResamplingConfig {
62    fn default() -> Self {
63        Self {
64            method: ResamplingMethod::Linear,
65            extrapolation: ExtrapolationMode::Nearest,
66        }
67    }
68}
69
70// ─── 1-D resampling ──────────────────────────────────────────────────────────
71
72/// Resample a 1-D signal from `(x_in, y_in)` to query points `x_out`.
73///
74/// `x_in` must be strictly increasing. `x_out` may be arbitrary.
75pub fn resample_1d(
76    x_in: &[f64],
77    y_in: &[f64],
78    x_out: &[f64],
79    config: &ResamplingConfig,
80) -> Result<Vec<f64>, InterpolateError> {
81    let n = x_in.len();
82    if n < 2 {
83        return Err(InterpolateError::InsufficientData(
84            "resample_1d requires at least 2 input points".to_string(),
85        ));
86    }
87    if n != y_in.len() {
88        return Err(InterpolateError::DimensionMismatch(format!(
89            "x_in length {} != y_in length {}",
90            n,
91            y_in.len()
92        )));
93    }
94
95    // Validate monotonicity
96    for i in 1..n {
97        if x_in[i] <= x_in[i - 1] {
98            return Err(InterpolateError::InvalidInput {
99                message: "x_in must be strictly increasing".to_string(),
100            });
101        }
102    }
103
104    // Precompute cubic spline coefficients if needed
105    let spline_coeffs: Option<Vec<[f64; 4]>> = match config.method {
106        ResamplingMethod::CubicSpline => Some(natural_cubic_spline_coeffs(x_in, y_in)?),
107        _ => None,
108    };
109
110    let x_min = x_in[0];
111    let x_max = x_in[n - 1];
112
113    let result: Result<Vec<f64>, InterpolateError> = x_out
114        .iter()
115        .map(|&xq| {
116            // Map possibly out-of-range xq to a resolved position
117            let xq_mapped = resolve_query(xq, x_min, x_max, &config.extrapolation);
118
119            match xq_mapped {
120                ResolvedQuery::InDomain(xr) => {
121                    interpolate_1d(x_in, y_in, xr, config, &spline_coeffs)
122                }
123                ResolvedQuery::Extrapolated(val) => Ok(val),
124                ResolvedQuery::ExtrapLinear(xr) => {
125                    // Linear extrapolation: use xr which may be outside domain
126                    interpolate_1d_linear_extrap(x_in, y_in, xr)
127                }
128                ResolvedQuery::ExtrapPolynomial(xr, deg) => {
129                    interpolate_1d_poly_extrap(x_in, y_in, xr, deg)
130                }
131            }
132        })
133        .collect();
134
135    result
136}
137
138// ─── Query resolution ────────────────────────────────────────────────────────
139
140enum ResolvedQuery {
141    InDomain(f64),
142    Extrapolated(f64),
143    ExtrapLinear(f64),
144    ExtrapPolynomial(f64, usize),
145}
146
147fn resolve_query(xq: f64, x_min: f64, x_max: f64, mode: &ExtrapolationMode) -> ResolvedQuery {
148    if xq >= x_min && xq <= x_max {
149        return ResolvedQuery::InDomain(xq);
150    }
151
152    match mode {
153        ExtrapolationMode::Nearest => ResolvedQuery::InDomain(xq.clamp(x_min, x_max)),
154        ExtrapolationMode::Linear => ResolvedQuery::ExtrapLinear(xq),
155        ExtrapolationMode::Polynomial(deg) => ResolvedQuery::ExtrapPolynomial(xq, *deg),
156        ExtrapolationMode::Reflection => {
157            let range = x_max - x_min;
158            if range < 1e-300 {
159                return ResolvedQuery::InDomain(x_min);
160            }
161            // Normalise to [0, 2*range) then reflect
162            let shifted = xq - x_min;
163            let period = 2.0 * range;
164            let t = shifted - (shifted / period).floor() * period;
165            let reflected = if t <= range { t } else { period - t };
166            ResolvedQuery::InDomain(x_min + reflected.clamp(0.0, range))
167        }
168        ExtrapolationMode::Periodic => {
169            let range = x_max - x_min;
170            if range < 1e-300 {
171                return ResolvedQuery::InDomain(x_min);
172            }
173            let shifted = xq - x_min;
174            let t = shifted - (shifted / range).floor() * range;
175            ResolvedQuery::InDomain(x_min + t.clamp(0.0, range))
176        }
177        ExtrapolationMode::Zero => ResolvedQuery::Extrapolated(0.0),
178        ExtrapolationMode::Constant(c) => ResolvedQuery::Extrapolated(*c),
179    }
180}
181
182// ─── 1-D interpolation methods ───────────────────────────────────────────────
183
184fn interpolate_1d(
185    x_in: &[f64],
186    y_in: &[f64],
187    xq: f64,
188    config: &ResamplingConfig,
189    spline_coeffs: &Option<Vec<[f64; 4]>>,
190) -> Result<f64, InterpolateError> {
191    let n = x_in.len();
192    let idx = binary_search_floor(x_in, xq);
193    let i = idx.min(n - 2);
194
195    match &config.method {
196        ResamplingMethod::Linear => {
197            let t = (xq - x_in[i]) / (x_in[i + 1] - x_in[i]);
198            Ok(y_in[i] * (1.0 - t) + y_in[i + 1] * t)
199        }
200        ResamplingMethod::Nearest => {
201            let i_near = if (xq - x_in[i]).abs() < (xq - x_in[(i + 1).min(n - 1)]).abs() {
202                i
203            } else {
204                (i + 1).min(n - 1)
205            };
206            Ok(y_in[i_near])
207        }
208        ResamplingMethod::CubicSpline => {
209            if let Some(coeffs) = spline_coeffs {
210                let dx = xq - x_in[i];
211                let [a, b, c, d] = coeffs[i];
212                Ok(a + b * dx + c * dx * dx + d * dx * dx * dx)
213            } else {
214                // Fallback to linear
215                let t = (xq - x_in[i]) / (x_in[i + 1] - x_in[i]);
216                Ok(y_in[i] * (1.0 - t) + y_in[i + 1] * t)
217            }
218        }
219        ResamplingMethod::Lanczos(a) => Ok(lanczos_interp(x_in, y_in, xq, *a)),
220    }
221}
222
223fn interpolate_1d_linear_extrap(
224    x_in: &[f64],
225    y_in: &[f64],
226    xq: f64,
227) -> Result<f64, InterpolateError> {
228    let n = x_in.len();
229    let x_min = x_in[0];
230    let x_max = x_in[n - 1];
231    if xq < x_min {
232        // Extrapolate left using first interval slope
233        let slope = (y_in[1] - y_in[0]) / (x_in[1] - x_in[0]);
234        Ok(y_in[0] + slope * (xq - x_min))
235    } else {
236        // Extrapolate right using last interval slope
237        let slope = (y_in[n - 1] - y_in[n - 2]) / (x_in[n - 1] - x_in[n - 2]);
238        Ok(y_in[n - 1] + slope * (xq - x_max))
239    }
240}
241
242fn interpolate_1d_poly_extrap(
243    x_in: &[f64],
244    y_in: &[f64],
245    xq: f64,
246    deg: usize,
247) -> Result<f64, InterpolateError> {
248    let n = x_in.len();
249    let x_min = x_in[0];
250    let pts = deg + 1;
251    // Pick boundary points
252    let (px, py): (Vec<f64>, Vec<f64>) = if xq < x_min {
253        // Use first `pts` points
254        let end = pts.min(n);
255        (x_in[..end].to_vec(), y_in[..end].to_vec())
256    } else {
257        // Use last `pts` points
258        let start = n.saturating_sub(pts);
259        (x_in[start..].to_vec(), y_in[start..].to_vec())
260    };
261
262    // Lagrange interpolation/extrapolation
263    Ok(lagrange_eval(&px, &py, xq))
264}
265
266// ─── Natural cubic spline coefficient computation ─────────────────────────────
267
268/// Compute natural cubic spline coefficients for `n-1` intervals.
269/// Returns coefficients `[a, b, c, d]` per interval such that
270/// `f(x) = a + b*(x-xi) + c*(x-xi)^2 + d*(x-xi)^3` for x in [xi, xi+1].
271fn natural_cubic_spline_coeffs(x: &[f64], y: &[f64]) -> Result<Vec<[f64; 4]>, InterpolateError> {
272    let n = x.len();
273    if n < 2 {
274        return Err(InterpolateError::InsufficientData(
275            "Need at least 2 points for spline".to_string(),
276        ));
277    }
278    let m = n - 1;
279    let mut h = vec![0.0f64; m];
280    for i in 0..m {
281        h[i] = x[i + 1] - x[i];
282        if h[i] <= 0.0 {
283            return Err(InterpolateError::InvalidInput {
284                message: "x must be strictly increasing".to_string(),
285            });
286        }
287    }
288
289    if n == 2 {
290        let b = (y[1] - y[0]) / h[0];
291        return Ok(vec![[y[0], b, 0.0, 0.0]]);
292    }
293
294    // Set up tridiagonal system for second derivatives σ
295    let mut alpha = vec![0.0f64; n];
296    for i in 1..m {
297        alpha[i] = 3.0 * ((y[i + 1] - y[i]) / h[i] - (y[i] - y[i - 1]) / h[i - 1]);
298    }
299
300    // Thomas algorithm (natural: σ_0 = σ_{n-1} = 0)
301    let mut l = vec![1.0f64; n];
302    let mut mu = vec![0.0f64; n];
303    let mut z = vec![0.0f64; n];
304
305    for i in 1..m {
306        l[i] = 2.0 * (x[i + 1] - x[i - 1]) - h[i - 1] * mu[i - 1];
307        if l[i].abs() < 1e-300 {
308            l[i] = 1e-300;
309        }
310        mu[i] = h[i] / l[i];
311        z[i] = (alpha[i] - h[i - 1] * z[i - 1]) / l[i];
312    }
313
314    let mut sigma = vec![0.0f64; n]; // second derivatives
315    for i in (1..m).rev() {
316        sigma[i] = z[i] - mu[i] * sigma[i + 1];
317    }
318
319    // Compute polynomial coefficients
320    let mut coeffs = Vec::with_capacity(m);
321    for i in 0..m {
322        let a = y[i];
323        let b = (y[i + 1] - y[i]) / h[i] - h[i] * (2.0 * sigma[i] + sigma[i + 1]) / 3.0;
324        let c = sigma[i];
325        let d = (sigma[i + 1] - sigma[i]) / (3.0 * h[i]);
326        coeffs.push([a, b, c, d]);
327    }
328
329    Ok(coeffs)
330}
331
332// ─── Lanczos windowed-sinc ────────────────────────────────────────────────────
333
334fn sinc(x: f64) -> f64 {
335    if x.abs() < 1e-12 {
336        1.0
337    } else {
338        let px = std::f64::consts::PI * x;
339        px.sin() / px
340    }
341}
342
343fn lanczos_kernel(x: f64, a: usize) -> f64 {
344    let af = a as f64;
345    if x.abs() >= af {
346        0.0
347    } else {
348        sinc(x) * sinc(x / af)
349    }
350}
351
352fn lanczos_interp(x_in: &[f64], y_in: &[f64], xq: f64, a: usize) -> f64 {
353    let n = x_in.len();
354    if n < 2 {
355        return y_in.first().copied().unwrap_or(0.0);
356    }
357    // Convert xq to fractional index in x_in (assume uniform spacing for simplicity,
358    // fallback to linear for non-uniform)
359    let i0 = binary_search_floor(x_in, xq);
360    let h = x_in[1] - x_in[0]; // approximate uniform step
361    if h.abs() < 1e-300 {
362        return y_in[i0.min(n - 1)];
363    }
364    let frac = (xq - x_in[i0.min(n - 1)]) / h;
365    let fi = i0 as f64 + frac;
366
367    let mut numer = 0.0f64;
368    let mut denom = 0.0f64;
369    let start = (fi as isize - a as isize).max(0) as usize;
370    let end = ((fi as isize + a as isize + 1) as usize).min(n);
371
372    for k in start..end {
373        let w = lanczos_kernel(fi - k as f64, a);
374        numer += w * y_in[k];
375        denom += w;
376    }
377
378    if denom.abs() < 1e-300 {
379        y_in[i0.min(n - 1)]
380    } else {
381        numer / denom
382    }
383}
384
385// ─── Helper: Lagrange interpolation ─────────────────────────────────────────
386
387fn lagrange_eval(px: &[f64], py: &[f64], xq: f64) -> f64 {
388    let n = px.len();
389    let mut result = 0.0f64;
390    for i in 0..n {
391        let mut li = 1.0f64;
392        for j in 0..n {
393            if i != j {
394                let denom = px[i] - px[j];
395                if denom.abs() < 1e-300 {
396                    continue;
397                }
398                li *= (xq - px[j]) / denom;
399            }
400        }
401        result += py[i] * li;
402    }
403    result
404}
405
406// ─── Helper: binary search (floor index) ─────────────────────────────────────
407
408/// Return index i such that x_in[i] <= xq < x_in[i+1].
409/// Clamps to [0, n-2].
410fn binary_search_floor(x_in: &[f64], xq: f64) -> usize {
411    let n = x_in.len();
412    if n == 0 {
413        return 0;
414    }
415    let mut lo = 0usize;
416    let mut hi = n - 1;
417    while lo + 1 < hi {
418        let mid = (lo + hi) / 2;
419        if x_in[mid] <= xq {
420            lo = mid;
421        } else {
422            hi = mid;
423        }
424    }
425    lo.min(n.saturating_sub(2))
426}
427
428// ─── 2-D resampling ──────────────────────────────────────────────────────────
429
430/// Resample a 2-D grid `grid[iy][ix]` from `(x_in, y_in)` axes to `(x_out, y_out)`.
431///
432/// Uses separable 1-D resampling: first along x, then along y.
433pub fn resample_2d(
434    grid: &[Vec<f64>],
435    x_in: &[f64],
436    y_in: &[f64],
437    x_out: &[f64],
438    y_out: &[f64],
439    config: &ResamplingConfig,
440) -> Result<Vec<Vec<f64>>, InterpolateError> {
441    let ny_in = y_in.len();
442    let nx_in = x_in.len();
443    if grid.len() != ny_in {
444        return Err(InterpolateError::DimensionMismatch(format!(
445            "grid has {} rows but y_in has {} elements",
446            grid.len(),
447            ny_in
448        )));
449    }
450    for (row_idx, row) in grid.iter().enumerate() {
451        if row.len() != nx_in {
452            return Err(InterpolateError::DimensionMismatch(format!(
453                "grid row {} has {} columns but x_in has {} elements",
454                row_idx,
455                row.len(),
456                nx_in
457            )));
458        }
459    }
460
461    // Step 1: resample along x for each input y row
462    // Produces intermediate grid of shape (ny_in, nx_out)
463    let mut intermediate: Vec<Vec<f64>> = Vec::with_capacity(ny_in);
464    for row in grid.iter() {
465        let resampled_row = resample_1d(x_in, row, x_out, config)?;
466        intermediate.push(resampled_row);
467    }
468
469    // Step 2: resample along y for each output x column
470    let nx_out = x_out.len();
471    let ny_out = y_out.len();
472    let mut output = vec![vec![0.0f64; nx_out]; ny_out];
473
474    for ix in 0..nx_out {
475        // Extract column from intermediate
476        let col: Vec<f64> = intermediate.iter().map(|row| row[ix]).collect();
477        let resampled_col = resample_1d(y_in, &col, y_out, config)?;
478        for iy in 0..ny_out {
479            output[iy][ix] = resampled_col[iy];
480        }
481    }
482
483    Ok(output)
484}
485
486// ─── Scattered to grid (IDW) ─────────────────────────────────────────────────
487
488/// Map scattered N-D points to a regular grid using inverse distance weighting.
489///
490/// `grid_ranges[d] = (min, max, n_points)` for dimension d.
491/// Returns a flattened array of shape `[n0 * n1 * ... * nd]` in row-major order.
492pub fn scattered_to_grid(
493    x: &[Vec<f64>],
494    y: &[f64],
495    grid_ranges: &[(f64, f64, usize)],
496    _config: &ResamplingConfig,
497) -> Result<Vec<f64>, InterpolateError> {
498    if x.is_empty() {
499        return Err(InterpolateError::InsufficientData(
500            "No scattered data points".to_string(),
501        ));
502    }
503    if x.len() != y.len() {
504        return Err(InterpolateError::DimensionMismatch(format!(
505            "x has {} rows but y has {} elements",
506            x.len(),
507            y.len()
508        )));
509    }
510    if grid_ranges.is_empty() {
511        return Err(InterpolateError::InvalidInput {
512            message: "grid_ranges must not be empty".to_string(),
513        });
514    }
515
516    let n_dims = grid_ranges.len();
517    let input_dims = x[0].len();
518    if input_dims != n_dims {
519        return Err(InterpolateError::DimensionMismatch(format!(
520            "x has {} dimensions but grid_ranges specifies {} dimensions",
521            input_dims, n_dims
522        )));
523    }
524
525    // Build grid axes
526    let axes: Vec<Vec<f64>> = grid_ranges
527        .iter()
528        .map(|&(lo, hi, n)| {
529            if n <= 1 {
530                vec![lo]
531            } else {
532                (0..n)
533                    .map(|i| lo + (hi - lo) * i as f64 / (n - 1) as f64)
534                    .collect()
535            }
536        })
537        .collect();
538
539    // Total grid size
540    let total: usize = axes.iter().map(|a| a.len()).product();
541    let mut result = vec![0.0f64; total];
542
543    // Enumerate multi-index
544    let shapes: Vec<usize> = axes.iter().map(|a| a.len()).collect();
545    let mut flat_idx = 0usize;
546
547    let mut multi = vec![0usize; n_dims];
548    loop {
549        // Build grid point coordinates
550        let gp: Vec<f64> = (0..n_dims).map(|d| axes[d][multi[d]]).collect();
551
552        // IDW with power p=2
553        let mut numer = 0.0f64;
554        let mut denom = 0.0f64;
555        for (xi, &yi) in x.iter().zip(y.iter()) {
556            let dist2: f64 = xi.iter().zip(gp.iter()).map(|(a, b)| (a - b).powi(2)).sum();
557            if dist2 < 1e-28 {
558                // Exact hit
559                numer = yi;
560                denom = 1.0;
561                break;
562            }
563            let w = 1.0 / dist2;
564            numer += w * yi;
565            denom += w;
566        }
567        result[flat_idx] = if denom > 1e-300 { numer / denom } else { 0.0 };
568
569        // Advance multi-index (row-major)
570        flat_idx += 1;
571        let mut carry = true;
572        for d in (0..n_dims).rev() {
573            if carry {
574                multi[d] += 1;
575                if multi[d] >= shapes[d] {
576                    multi[d] = 0;
577                } else {
578                    carry = false;
579                }
580            }
581        }
582        if carry {
583            break; // all indices wrapped
584        }
585    }
586
587    Ok(result)
588}
589
590// ─── SplineDerivative ────────────────────────────────────────────────────────
591
592/// A piecewise polynomial (spline) on a set of knot intervals.
593///
594/// Each segment `[knots[i], knots[i+1]]` is represented by a polynomial
595/// of degree `degree` with coefficients `coefficients[i]` stored in
596/// *ascending order* (coefficient of x^0 first).
597#[derive(Debug, Clone)]
598pub struct SplineDerivative {
599    /// Polynomial coefficients per segment, `coefficients[i][k]` = coeff of `(x - knots[i])^k`.
600    pub coefficients: Vec<Vec<f64>>,
601    /// Knot values (segment boundaries), length = n_segments + 1.
602    pub knots: Vec<f64>,
603    /// Polynomial degree of each segment.
604    pub degree: usize,
605}
606
607impl SplineDerivative {
608    /// Create a new spline from coefficients, knots, and degree.
609    pub fn new(
610        coefficients: Vec<Vec<f64>>,
611        knots: Vec<f64>,
612        degree: usize,
613    ) -> Result<Self, InterpolateError> {
614        if knots.len() < 2 {
615            return Err(InterpolateError::InsufficientData(
616                "SplineDerivative needs at least 2 knots".to_string(),
617            ));
618        }
619        let n_seg = knots.len() - 1;
620        if coefficients.len() != n_seg {
621            return Err(InterpolateError::DimensionMismatch(format!(
622                "Expected {} coefficient vectors for {} segments, got {}",
623                n_seg,
624                n_seg,
625                coefficients.len()
626            )));
627        }
628        Ok(Self {
629            coefficients,
630            knots,
631            degree,
632        })
633    }
634
635    /// Differentiate this spline, returning a new spline of `degree - 1`.
636    pub fn differentiate(spline: &SplineDerivative) -> Result<Self, InterpolateError> {
637        if spline.degree == 0 {
638            return Err(InterpolateError::InvalidOperation(
639                "Cannot differentiate a degree-0 spline".to_string(),
640            ));
641        }
642        let new_degree = spline.degree - 1;
643        let new_coeffs: Vec<Vec<f64>> = spline
644            .coefficients
645            .iter()
646            .map(|seg_coeffs| {
647                // Differentiate polynomial: d/dx (c_k (x-x_i)^k) = k * c_k * (x-x_i)^(k-1)
648                // Result has one fewer coefficient
649                let n = seg_coeffs.len().min(spline.degree + 1);
650                (1..n)
651                    .map(|k| k as f64 * seg_coeffs[k])
652                    .collect::<Vec<f64>>()
653            })
654            .collect();
655
656        Self::new(new_coeffs, spline.knots.clone(), new_degree)
657    }
658
659    /// Evaluate the spline at point `x`.
660    pub fn evaluate(&self, x: f64) -> Result<f64, InterpolateError> {
661        let n = self.knots.len();
662        if n < 2 {
663            return Err(InterpolateError::InsufficientData(
664                "No segments to evaluate".to_string(),
665            ));
666        }
667
668        // Find segment
669        let seg = if x <= self.knots[0] {
670            0
671        } else if x >= self.knots[n - 1] {
672            n - 2
673        } else {
674            binary_search_floor(&self.knots, x)
675        };
676
677        let dx = x - self.knots[seg];
678        let coeffs = &self.coefficients[seg];
679        // Horner's method
680        let mut val = 0.0f64;
681        for &c in coeffs.iter().rev() {
682            val = val * dx + c;
683        }
684        Ok(val)
685    }
686}
687
688// ─── Grid Resampling Convenience Functions (WS227) ──────────────────────────
689
690/// Resample scattered 1-D data onto a uniform grid of `n_grid_points`.
691///
692/// Returns `(grid_x, grid_y)` where `grid_x` is evenly spaced over
693/// `[min(scattered_x), max(scattered_x)]`.
694pub fn resample_to_regular(
695    scattered_x: &[f64],
696    scattered_y: &[f64],
697    n_grid_points: usize,
698    config: &ResamplingConfig,
699) -> Result<(Vec<f64>, Vec<f64>), InterpolateError> {
700    if scattered_x.len() < 2 {
701        return Err(InterpolateError::InsufficientData(
702            "resample_to_regular requires at least 2 input points".to_string(),
703        ));
704    }
705    if scattered_x.len() != scattered_y.len() {
706        return Err(InterpolateError::DimensionMismatch(format!(
707            "scattered_x len {} != scattered_y len {}",
708            scattered_x.len(),
709            scattered_y.len()
710        )));
711    }
712    if n_grid_points < 2 {
713        return Err(InterpolateError::InvalidInput {
714            message: "n_grid_points must be >= 2".to_string(),
715        });
716    }
717
718    // Sort the input data by x.
719    let mut pairs: Vec<(f64, f64)> = scattered_x
720        .iter()
721        .copied()
722        .zip(scattered_y.iter().copied())
723        .collect();
724    pairs.sort_by(|a, b| a.0.partial_cmp(&b.0).unwrap_or(std::cmp::Ordering::Equal));
725
726    // Deduplicate by x (keep last y for each unique x).
727    let mut sorted_x: Vec<f64> = Vec::with_capacity(pairs.len());
728    let mut sorted_y: Vec<f64> = Vec::with_capacity(pairs.len());
729    for &(px, py) in &pairs {
730        if let Some(&last_x) = sorted_x.last() {
731            if (px - last_x).abs() < 1e-15_f64 {
732                // Replace y for duplicate x.
733                if let Some(ly) = sorted_y.last_mut() {
734                    *ly = py;
735                }
736                continue;
737            }
738        }
739        sorted_x.push(px);
740        sorted_y.push(py);
741    }
742
743    if sorted_x.len() < 2 {
744        return Err(InterpolateError::InsufficientData(
745            "After deduplication, fewer than 2 unique x values remain".to_string(),
746        ));
747    }
748
749    let x_min = sorted_x[0];
750    let x_max = sorted_x[sorted_x.len() - 1];
751    let step = (x_max - x_min) / (n_grid_points - 1) as f64;
752    let grid_x: Vec<f64> = (0..n_grid_points)
753        .map(|i| x_min + i as f64 * step)
754        .collect();
755
756    let grid_y = resample_1d(&sorted_x, &sorted_y, &grid_x, config)?;
757    Ok((grid_x, grid_y))
758}
759
760/// Resample 1-D data onto arbitrary target x-coordinates.
761///
762/// `data_x` must be sortable; it will be sorted internally.
763pub fn resample_to_irregular(
764    data_x: &[f64],
765    data_y: &[f64],
766    target_x: &[f64],
767    config: &ResamplingConfig,
768) -> Result<Vec<f64>, InterpolateError> {
769    if data_x.len() < 2 {
770        return Err(InterpolateError::InsufficientData(
771            "resample_to_irregular requires at least 2 input points".to_string(),
772        ));
773    }
774    if data_x.len() != data_y.len() {
775        return Err(InterpolateError::DimensionMismatch(format!(
776            "data_x len {} != data_y len {}",
777            data_x.len(),
778            data_y.len()
779        )));
780    }
781
782    // Sort input by x.
783    let mut pairs: Vec<(f64, f64)> = data_x.iter().copied().zip(data_y.iter().copied()).collect();
784    pairs.sort_by(|a, b| a.0.partial_cmp(&b.0).unwrap_or(std::cmp::Ordering::Equal));
785
786    let sorted_x: Vec<f64> = pairs.iter().map(|p| p.0).collect();
787    let sorted_y: Vec<f64> = pairs.iter().map(|p| p.1).collect();
788
789    resample_1d(&sorted_x, &sorted_y, target_x, config)
790}
791
792/// Resample 2-D scattered data onto a regular nx × ny grid.
793///
794/// Uses inverse-distance weighting to map scattered points to a grid covering
795/// `[min_x, max_x] × [min_y, max_y]`.
796///
797/// Returns a `grid_ny × grid_nx` nested `Vec<Vec<f64>>` in row-major order.
798pub fn resample_scattered_2d(
799    scattered_xy: &[(f64, f64)],
800    values: &[f64],
801    grid_nx: usize,
802    grid_ny: usize,
803) -> Result<Vec<Vec<f64>>, InterpolateError> {
804    if scattered_xy.is_empty() {
805        return Err(InterpolateError::InsufficientData(
806            "No scattered data points for 2D resampling".to_string(),
807        ));
808    }
809    if scattered_xy.len() != values.len() {
810        return Err(InterpolateError::DimensionMismatch(format!(
811            "scattered_xy len {} != values len {}",
812            scattered_xy.len(),
813            values.len()
814        )));
815    }
816    if grid_nx < 2 || grid_ny < 2 {
817        return Err(InterpolateError::InvalidInput {
818            message: "grid_nx and grid_ny must each be >= 2".to_string(),
819        });
820    }
821
822    let x_vals: Vec<f64> = scattered_xy.iter().map(|p| p.0).collect();
823    let y_vals: Vec<f64> = scattered_xy.iter().map(|p| p.1).collect();
824
825    let x_min = x_vals.iter().copied().fold(f64::INFINITY, f64::min);
826    let x_max = x_vals.iter().copied().fold(f64::NEG_INFINITY, f64::max);
827    let y_min = y_vals.iter().copied().fold(f64::INFINITY, f64::min);
828    let y_max = y_vals.iter().copied().fold(f64::NEG_INFINITY, f64::max);
829
830    // Prevent degenerate grids.
831    let dx = if (x_max - x_min).abs() < 1e-15 {
832        1.0
833    } else {
834        (x_max - x_min) / (grid_nx - 1) as f64
835    };
836    let dy = if (y_max - y_min).abs() < 1e-15 {
837        1.0
838    } else {
839        (y_max - y_min) / (grid_ny - 1) as f64
840    };
841
842    let mut grid = vec![vec![0.0f64; grid_nx]; grid_ny];
843
844    for iy in 0..grid_ny {
845        let gy = y_min + iy as f64 * dy;
846        for ix in 0..grid_nx {
847            let gx = x_min + ix as f64 * dx;
848
849            // IDW with power 2
850            let mut numer = 0.0_f64;
851            let mut denom = 0.0_f64;
852            let mut exact_hit = false;
853            for (idx, &(sx, sy)) in scattered_xy.iter().enumerate() {
854                let dist2 = (sx - gx).powi(2) + (sy - gy).powi(2);
855                if dist2 < 1e-28 {
856                    grid[iy][ix] = values[idx];
857                    exact_hit = true;
858                    break;
859                }
860                let w = 1.0 / dist2;
861                numer += w * values[idx];
862                denom += w;
863            }
864            if !exact_hit {
865                grid[iy][ix] = if denom > 1e-300 { numer / denom } else { 0.0 };
866            }
867        }
868    }
869
870    Ok(grid)
871}
872
873/// Downsample a 1-D signal by keeping every `factor`-th point.
874///
875/// The first point is always kept. Returns `(x_out, y_out)`.
876pub fn downsample(
877    x: &[f64],
878    y: &[f64],
879    factor: usize,
880) -> Result<(Vec<f64>, Vec<f64>), InterpolateError> {
881    if x.len() != y.len() {
882        return Err(InterpolateError::DimensionMismatch(format!(
883            "x len {} != y len {}",
884            x.len(),
885            y.len()
886        )));
887    }
888    if factor == 0 {
889        return Err(InterpolateError::InvalidInput {
890            message: "downsample factor must be >= 1".to_string(),
891        });
892    }
893    if x.is_empty() {
894        return Ok((Vec::new(), Vec::new()));
895    }
896
897    let x_out: Vec<f64> = x.iter().copied().step_by(factor).collect();
898    let y_out: Vec<f64> = y.iter().copied().step_by(factor).collect();
899    Ok((x_out, y_out))
900}
901
902/// Upsample a 1-D signal by inserting `factor - 1` interpolated points
903/// between each original pair.
904///
905/// Uses the configured resampling method (default: linear).
906pub fn upsample(
907    x: &[f64],
908    y: &[f64],
909    factor: usize,
910    config: &ResamplingConfig,
911) -> Result<(Vec<f64>, Vec<f64>), InterpolateError> {
912    if x.len() != y.len() {
913        return Err(InterpolateError::DimensionMismatch(format!(
914            "x len {} != y len {}",
915            x.len(),
916            y.len()
917        )));
918    }
919    if factor == 0 {
920        return Err(InterpolateError::InvalidInput {
921            message: "upsample factor must be >= 1".to_string(),
922        });
923    }
924    if x.len() < 2 {
925        return Ok((x.to_vec(), y.to_vec()));
926    }
927
928    // Sort input by x.
929    let mut pairs: Vec<(f64, f64)> = x.iter().copied().zip(y.iter().copied()).collect();
930    pairs.sort_by(|a, b| a.0.partial_cmp(&b.0).unwrap_or(std::cmp::Ordering::Equal));
931
932    let sorted_x: Vec<f64> = pairs.iter().map(|p| p.0).collect();
933    let sorted_y: Vec<f64> = pairs.iter().map(|p| p.1).collect();
934
935    let n = sorted_x.len();
936    // Total output points: (n - 1) * factor + 1
937    let n_out = (n - 1) * factor + 1;
938    let mut x_out = Vec::with_capacity(n_out);
939
940    for i in 0..(n - 1) {
941        let x0 = sorted_x[i];
942        let x1 = sorted_x[i + 1];
943        for j in 0..factor {
944            let t = j as f64 / factor as f64;
945            x_out.push(x0 + t * (x1 - x0));
946        }
947    }
948    // Include the last point.
949    x_out.push(sorted_x[n - 1]);
950
951    let y_out = resample_1d(&sorted_x, &sorted_y, &x_out, config)?;
952    Ok((x_out, y_out))
953}
954
955// ─── Tests ───────────────────────────────────────────────────────────────────
956
957#[cfg(test)]
958mod tests {
959    use super::*;
960
961    #[test]
962    fn test_resample_1d_linear_identity() {
963        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
964        let y: Vec<f64> = x.clone();
965        let config = ResamplingConfig {
966            method: ResamplingMethod::Linear,
967            extrapolation: ExtrapolationMode::Nearest,
968        };
969        let x_out: Vec<f64> = (0..10).map(|i| i as f64 * 0.5 + 0.5).collect();
970        let result = resample_1d(&x, &y, &x_out, &config).expect("resample");
971        for (got, expected) in result.iter().zip(x_out.iter()) {
972            // y = x, so linear interpolation should be exact
973            let x_clamped = expected.clamp(x[0], x[x.len() - 1]);
974            assert!(
975                (got - x_clamped).abs() < 1e-10,
976                "Linear identity failed: got={got}, expected={x_clamped}"
977            );
978        }
979    }
980
981    #[test]
982    fn test_extrapolation_nearest_boundary() {
983        let x: Vec<f64> = vec![0.0, 1.0, 2.0, 3.0];
984        let y: Vec<f64> = vec![10.0, 20.0, 30.0, 40.0];
985        let config = ResamplingConfig {
986            method: ResamplingMethod::Linear,
987            extrapolation: ExtrapolationMode::Nearest,
988        };
989        let x_out = vec![-1.0, 5.0];
990        let result = resample_1d(&x, &y, &x_out, &config).expect("resample");
991        assert!(
992            (result[0] - 10.0).abs() < 1e-10,
993            "Left boundary clamped: {}",
994            result[0]
995        );
996        assert!(
997            (result[1] - 40.0).abs() < 1e-10,
998            "Right boundary clamped: {}",
999            result[1]
1000        );
1001    }
1002
1003    #[test]
1004    fn test_extrapolation_zero() {
1005        let x: Vec<f64> = vec![0.0, 1.0, 2.0];
1006        let y: Vec<f64> = vec![1.0, 2.0, 3.0];
1007        let config = ResamplingConfig {
1008            method: ResamplingMethod::Linear,
1009            extrapolation: ExtrapolationMode::Zero,
1010        };
1011        let x_out = vec![-1.0, 5.0];
1012        let result = resample_1d(&x, &y, &x_out, &config).expect("resample");
1013        assert!((result[0] - 0.0).abs() < 1e-10);
1014        assert!((result[1] - 0.0).abs() < 1e-10);
1015    }
1016
1017    #[test]
1018    fn test_extrapolation_constant() {
1019        let x: Vec<f64> = vec![0.0, 1.0, 2.0];
1020        let y: Vec<f64> = vec![1.0, 2.0, 3.0];
1021        let config = ResamplingConfig {
1022            method: ResamplingMethod::Linear,
1023            extrapolation: ExtrapolationMode::Constant(99.0),
1024        };
1025        let x_out = vec![-5.0, 10.0];
1026        let result = resample_1d(&x, &y, &x_out, &config).expect("resample");
1027        assert!((result[0] - 99.0).abs() < 1e-10);
1028        assert!((result[1] - 99.0).abs() < 1e-10);
1029    }
1030
1031    #[test]
1032    fn test_extrapolation_periodic() {
1033        let x: Vec<f64> = vec![0.0, 1.0, 2.0, 3.0];
1034        let y: Vec<f64> = vec![0.0, 1.0, 2.0, 3.0]; // y = x in domain
1035        let config = ResamplingConfig {
1036            method: ResamplingMethod::Linear,
1037            extrapolation: ExtrapolationMode::Periodic,
1038        };
1039        // x=3.5 should map to x=0.5 (period = 3)
1040        let result = resample_1d(&x, &y, &[3.5], &config).expect("periodic");
1041        assert!(
1042            (result[0] - 0.5).abs() < 0.2,
1043            "Periodic wrap: got {} expected ~0.5",
1044            result[0]
1045        );
1046    }
1047
1048    #[test]
1049    fn test_extrapolation_linear() {
1050        let x: Vec<f64> = vec![0.0, 1.0, 2.0, 3.0];
1051        let y: Vec<f64> = vec![0.0, 2.0, 4.0, 6.0]; // y = 2x
1052        let config = ResamplingConfig {
1053            method: ResamplingMethod::Linear,
1054            extrapolation: ExtrapolationMode::Linear,
1055        };
1056        // x=4.0 should extrapolate to 8.0
1057        let result = resample_1d(&x, &y, &[4.0], &config).expect("linear extrap");
1058        assert!(
1059            (result[0] - 8.0).abs() < 1e-8,
1060            "Linear extrapolation: got {} expected 8.0",
1061            result[0]
1062        );
1063    }
1064
1065    #[test]
1066    fn test_cubic_spline_resample() {
1067        let x: Vec<f64> = (0..6).map(|i| i as f64).collect();
1068        let y: Vec<f64> = x.iter().map(|&xi| xi * xi).collect(); // y = x²
1069        let config = ResamplingConfig {
1070            method: ResamplingMethod::CubicSpline,
1071            extrapolation: ExtrapolationMode::Nearest,
1072        };
1073        let x_out = vec![0.5, 1.5, 2.5, 3.5];
1074        let result = resample_1d(&x, &y, &x_out, &config).expect("cubic");
1075        for (xq, &yq) in x_out.iter().zip(result.iter()) {
1076            let exact = xq * xq;
1077            // Natural cubic spline imposes zero second derivatives at boundaries,
1078            // which introduces a boundary error for polynomial data.
1079            // Tolerance is relaxed near boundaries (x<1 or x>4).
1080            let tol = if *xq < 1.0 || *xq > 4.0 { 0.2 } else { 0.05 };
1081            assert!(
1082                (yq - exact).abs() < tol,
1083                "Cubic spline on y=x² at x={xq}: got {yq}, expected {exact}"
1084            );
1085        }
1086    }
1087
1088    #[test]
1089    fn test_nearest_method() {
1090        let x: Vec<f64> = vec![0.0, 1.0, 2.0];
1091        let y: Vec<f64> = vec![10.0, 20.0, 30.0];
1092        let config = ResamplingConfig {
1093            method: ResamplingMethod::Nearest,
1094            extrapolation: ExtrapolationMode::Nearest,
1095        };
1096        let result = resample_1d(&x, &y, &[0.3, 0.7], &config).expect("nearest");
1097        assert!(
1098            (result[0] - 10.0).abs() < 1e-10,
1099            "Nearest left: {}",
1100            result[0]
1101        );
1102        assert!(
1103            (result[1] - 20.0).abs() < 1e-10,
1104            "Nearest right: {}",
1105            result[1]
1106        );
1107    }
1108
1109    #[test]
1110    fn test_scattered_to_grid_2d() {
1111        // Simple scattered points: y = x1 + x2
1112        let x: Vec<Vec<f64>> = vec![
1113            vec![0.0, 0.0],
1114            vec![1.0, 0.0],
1115            vec![0.0, 1.0],
1116            vec![1.0, 1.0],
1117            vec![0.5, 0.5],
1118        ];
1119        let y: Vec<f64> = x.iter().map(|xi| xi[0] + xi[1]).collect();
1120        let grid_ranges = vec![(0.0, 1.0, 3), (0.0, 1.0, 3)];
1121        let config = ResamplingConfig::default();
1122        let result = scattered_to_grid(&x, &y, &grid_ranges, &config).expect("scattered_to_grid");
1123        assert_eq!(result.len(), 9); // 3×3
1124                                     // All values should be >= 0 and <= 2
1125        for &v in &result {
1126            assert!(v >= -0.1 && v <= 2.1, "Value out of expected range: {v}");
1127        }
1128    }
1129
1130    #[test]
1131    fn test_spline_derivative_differentiation() {
1132        // Quadratic: y = 3x² + 2x + 1 on [0, 2]
1133        // Coefficients: [1, 2, 3] for (x-0)^0, (x-0)^1, (x-0)^2
1134        let spline = SplineDerivative::new(vec![vec![1.0, 2.0, 3.0]], vec![0.0, 2.0], 2)
1135            .expect("create spline");
1136
1137        // Derivative: 2 + 6x => coeffs [2, 6]
1138        let deriv = SplineDerivative::differentiate(&spline).expect("differentiate");
1139        assert_eq!(deriv.degree, 1);
1140
1141        // Evaluate at x=1: should be 2 + 6*1 = 8
1142        let val = deriv.evaluate(1.0).expect("evaluate");
1143        assert!(
1144            (val - 8.0).abs() < 1e-10,
1145            "Derivative at x=1: got {val}, expected 8.0"
1146        );
1147    }
1148
1149    #[test]
1150    fn test_spline_evaluate() {
1151        // Linear: y = 2x + 1 on [0, 1] and y = 3x - 0 on [1, 2]
1152        let spline =
1153            SplineDerivative::new(vec![vec![1.0, 2.0], vec![3.0, 0.0]], vec![0.0, 1.0, 2.0], 1)
1154                .expect("create spline");
1155        let v0 = spline.evaluate(0.5).expect("eval");
1156        assert!((v0 - 2.0).abs() < 1e-10, "Got {v0}"); // 1 + 2*0.5 = 2
1157    }
1158
1159    #[test]
1160    fn test_resample_2d() {
1161        // 3×3 grid with values row + col
1162        let grid: Vec<Vec<f64>> = (0..3)
1163            .map(|i| (0..3).map(|j| (i + j) as f64).collect())
1164            .collect();
1165        let x_in: Vec<f64> = vec![0.0, 1.0, 2.0];
1166        let y_in: Vec<f64> = vec![0.0, 1.0, 2.0];
1167        let x_out: Vec<f64> = vec![0.5, 1.0, 1.5];
1168        let y_out: Vec<f64> = vec![0.5, 1.0, 1.5];
1169        let config = ResamplingConfig::default();
1170        let result =
1171            resample_2d(&grid, &x_in, &y_in, &x_out, &y_out, &config).expect("resample_2d");
1172        assert_eq!(result.len(), 3);
1173        assert_eq!(result[0].len(), 3);
1174        // At (0.5, 0.5): value should be ~1.0
1175        assert!(
1176            (result[0][0] - 1.0).abs() < 0.1,
1177            "2D resample at (0.5,0.5): got {}",
1178            result[0][0]
1179        );
1180    }
1181
1182    #[test]
1183    fn test_reflection_extrapolation() {
1184        let x: Vec<f64> = vec![0.0, 1.0, 2.0, 3.0];
1185        let y: Vec<f64> = vec![0.0, 1.0, 4.0, 9.0];
1186        let config = ResamplingConfig {
1187            method: ResamplingMethod::Linear,
1188            extrapolation: ExtrapolationMode::Reflection,
1189        };
1190        // x=-0.5 should reflect to 0.5 within [0,3]
1191        let result = resample_1d(&x, &y, &[-0.5], &config).expect("reflection");
1192        // Should be in range since reflected
1193        assert!(
1194            result[0].is_finite(),
1195            "Reflection should produce finite value"
1196        );
1197    }
1198
1199    // ── WS227 Grid Resampling tests ──────────────────────────────────────
1200
1201    #[test]
1202    fn test_resample_to_regular_roundtrip() {
1203        // y = 2x on [0, 4] — resample to 5-point grid, then evaluate at original.
1204        let x: Vec<f64> = vec![0.0, 1.0, 2.0, 3.0, 4.0];
1205        let y: Vec<f64> = x.iter().map(|&xi| 2.0 * xi).collect();
1206        let config = ResamplingConfig::default();
1207
1208        let (grid_x, grid_y) =
1209            resample_to_regular(&x, &y, 9, &config).expect("resample_to_regular");
1210        assert_eq!(grid_x.len(), 9);
1211        assert_eq!(grid_y.len(), 9);
1212
1213        // At each grid point the value should be close to 2*x.
1214        for (gx, gy) in grid_x.iter().zip(grid_y.iter()) {
1215            let expected = 2.0 * gx;
1216            assert!(
1217                (gy - expected).abs() < 0.1,
1218                "resample_to_regular roundtrip: at x={gx} got {gy}, expected {expected}"
1219            );
1220        }
1221    }
1222
1223    #[test]
1224    fn test_resample_to_irregular() {
1225        let x: Vec<f64> = vec![0.0, 1.0, 2.0, 3.0, 4.0];
1226        let y: Vec<f64> = x.iter().map(|&xi| xi * xi).collect(); // y = x^2
1227        let config = ResamplingConfig::default();
1228
1229        let targets = vec![0.5, 1.5, 2.5, 3.5];
1230        let result =
1231            resample_to_irregular(&x, &y, &targets, &config).expect("resample_to_irregular");
1232        assert_eq!(result.len(), 4);
1233
1234        // Linear interpolation on x^2: approximate values.
1235        for (i, &tgt) in targets.iter().enumerate() {
1236            let exact = tgt * tgt;
1237            // Linear interp on x^2 has some error, but should be reasonable.
1238            assert!(
1239                (result[i] - exact).abs() < 1.0,
1240                "resample_to_irregular at x={tgt}: got {}, expected ~{exact}",
1241                result[i]
1242            );
1243        }
1244    }
1245
1246    #[test]
1247    fn test_resample_scattered_2d_grid_covers_domain() {
1248        let scattered: Vec<(f64, f64)> =
1249            vec![(0.0, 0.0), (1.0, 0.0), (0.0, 1.0), (1.0, 1.0), (0.5, 0.5)];
1250        let values: Vec<f64> = scattered.iter().map(|&(x, y)| x + y).collect();
1251
1252        let grid = resample_scattered_2d(&scattered, &values, 3, 3).expect("scattered 2d");
1253        assert_eq!(grid.len(), 3);
1254        assert_eq!(grid[0].len(), 3);
1255
1256        // All values should be in [0, 2] (since x+y ranges from 0 to 2).
1257        for row in &grid {
1258            for &v in row {
1259                assert!(v >= -0.1 && v <= 2.1, "2D grid value out of range: {v}");
1260            }
1261        }
1262    }
1263
1264    #[test]
1265    fn test_downsample() {
1266        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
1267        let y: Vec<f64> = x.iter().map(|&xi| xi * xi).collect();
1268
1269        let (dx, dy) = downsample(&x, &y, 3).expect("downsample");
1270        // With factor 3: indices 0, 3, 6, 9
1271        assert_eq!(dx.len(), 4);
1272        assert!((dx[0] - 0.0).abs() < 1e-12);
1273        assert!((dx[1] - 3.0).abs() < 1e-12);
1274        assert!((dx[2] - 6.0).abs() < 1e-12);
1275        assert!((dx[3] - 9.0).abs() < 1e-12);
1276        // Values should match y = x^2
1277        assert!((dy[1] - 9.0).abs() < 1e-12);
1278    }
1279
1280    #[test]
1281    fn test_upsample_preserves_function() {
1282        let x: Vec<f64> = vec![0.0, 1.0, 2.0, 3.0];
1283        let y: Vec<f64> = vec![0.0, 2.0, 4.0, 6.0]; // y = 2x
1284        let config = ResamplingConfig::default();
1285
1286        let (ux, uy) = upsample(&x, &y, 3, &config).expect("upsample");
1287        // (n-1)*factor + 1 = 3*3 + 1 = 10
1288        assert_eq!(ux.len(), 10);
1289        assert_eq!(uy.len(), 10);
1290
1291        // For y=2x, all upsampled values should be very close to 2*x.
1292        for (xi, yi) in ux.iter().zip(uy.iter()) {
1293            let expected = 2.0 * xi;
1294            assert!(
1295                (yi - expected).abs() < 0.1,
1296                "upsample: at x={xi} got {yi}, expected {expected}"
1297            );
1298        }
1299    }
1300
1301    #[test]
1302    fn test_downsample_factor_1_identity() {
1303        let x: Vec<f64> = vec![0.0, 1.0, 2.0];
1304        let y: Vec<f64> = vec![1.0, 2.0, 3.0];
1305        let (dx, dy) = downsample(&x, &y, 1).expect("factor 1");
1306        assert_eq!(dx.len(), 3);
1307        assert_eq!(dy.len(), 3);
1308    }
1309
1310    #[test]
1311    fn test_downsample_factor_zero_error() {
1312        let result = downsample(&[1.0], &[1.0], 0);
1313        assert!(result.is_err());
1314    }
1315
1316    #[test]
1317    fn test_upsample_factor_zero_error() {
1318        let config = ResamplingConfig::default();
1319        let result = upsample(&[1.0, 2.0], &[1.0, 2.0], 0, &config);
1320        assert!(result.is_err());
1321    }
1322
1323    #[test]
1324    fn test_resample_to_regular_too_few_points() {
1325        let config = ResamplingConfig::default();
1326        let result = resample_to_regular(&[1.0], &[1.0], 5, &config);
1327        assert!(result.is_err());
1328    }
1329}