Skip to main content

ferray_window/windows/
mod.rs

1// ferray-window: Window functions for signal processing and spectral analysis
2//
3// Implements NumPy- and SciPy-equivalent window functions:
4//   - bartlett, blackman, hamming, hanning, kaiser (NumPy core)
5//   - cosine, exponential, gaussian, general_cosine, general_hamming,
6//     nuttall, parzen, taylor, tukey (SciPy / torch.signal.windows extras)
7// Each returns an Array1<f64> of length M.
8
9// Window-function coefficients (Blackman 0.42/0.5/0.08, Hamming 0.54/0.46,
10// Kaiser modified-Bessel I0 series) are textbook constants reproduced
11// without separators to match the reference papers / NumPy source.
12#![allow(clippy::unreadable_literal)]
13
14use ferray_core::Array;
15use ferray_core::dimension::Ix1;
16use ferray_core::error::{FerrayError, FerrayResult};
17
18use std::f64::consts::PI;
19
20// Re-use the Abramowitz & Stegun Bessel I_0 polynomial approximation
21// from ferray-ufunc instead of shipping a second copy (see #530).
22// The previous copy had already drifted in whitespace; funnelling
23// through the canonical crate keeps both the window code and the
24// ufunc in sync.
25use ferray_ufunc::ops::special::bessel_i0_scalar;
26
27/// Build a window array of length `m` by evaluating `f` at every sample
28/// `0..m`. Handles the shared `m == 0` / `m == 1` edge cases that every
29/// NumPy-equivalent window has (see issue #292 for the original
30/// five-way copy of this pattern).
31#[inline]
32fn gen_window<F: FnMut(usize) -> f64>(m: usize, mut f: F) -> FerrayResult<Array<f64, Ix1>> {
33    if m == 0 {
34        return Array::from_vec(Ix1::new([0]), vec![]);
35    }
36    if m == 1 {
37        return Array::from_vec(Ix1::new([1]), vec![1.0]);
38    }
39    let mut data = Vec::with_capacity(m);
40    for n in 0..m {
41        data.push(f(n));
42    }
43    Array::from_vec(Ix1::new([m]), data)
44}
45
46/// Internal helper: cast a length-`m` f64 window into an f32 array.
47/// We compute in f64 (ferray-window's native precision) and then
48/// down-cast — the windows are smooth, so this preserves all useful
49/// precision.
50fn cast_to_f32(arr: &Array<f64, Ix1>) -> FerrayResult<Array<f32, Ix1>> {
51    let data: Vec<f32> = arr.iter().map(|&x| x as f32).collect();
52    Array::<f32, Ix1>::from_vec(Ix1::new([data.len()]), data)
53}
54
55/// `f32` version of [`bartlett`] (#529).
56///
57/// # Errors
58/// Same as `bartlett`.
59pub fn bartlett_f32(m: usize) -> FerrayResult<Array<f32, Ix1>> {
60    cast_to_f32(&bartlett(m)?)
61}
62
63/// `f32` version of [`blackman`] (#529).
64///
65/// # Errors
66/// Same as `blackman`.
67pub fn blackman_f32(m: usize) -> FerrayResult<Array<f32, Ix1>> {
68    cast_to_f32(&blackman(m)?)
69}
70
71/// `f32` version of [`hamming`] (#529).
72///
73/// # Errors
74/// Same as `hamming`.
75pub fn hamming_f32(m: usize) -> FerrayResult<Array<f32, Ix1>> {
76    cast_to_f32(&hamming(m)?)
77}
78
79/// `f32` version of [`hanning`] (#529).
80///
81/// # Errors
82/// Same as `hanning`.
83pub fn hanning_f32(m: usize) -> FerrayResult<Array<f32, Ix1>> {
84    cast_to_f32(&hanning(m)?)
85}
86
87/// `f32` version of [`kaiser`] (#529).
88///
89/// # Errors
90/// Same as `kaiser`.
91pub fn kaiser_f32(m: usize, beta: f64) -> FerrayResult<Array<f32, Ix1>> {
92    cast_to_f32(&kaiser(m, beta)?)
93}
94
95/// Return the Bartlett (triangular) window of length `m`.
96///
97/// The Bartlett window is defined as:
98///   w(n) = 2/(M-1) * ((M-1)/2 - |n - (M-1)/2|)
99///
100/// This is equivalent to `numpy.bartlett(M)`.
101///
102/// # Errors
103/// Returns an error only if internal array construction fails.
104pub fn bartlett(m: usize) -> FerrayResult<Array<f64, Ix1>> {
105    let half = (m.saturating_sub(1)) as f64 / 2.0;
106    gen_window(m, |n| 1.0 - ((n as f64 - half) / half).abs())
107}
108
109/// Return the Blackman window of length `m`.
110///
111/// The Blackman window is defined as:
112///   w(n) = 0.42 - 0.5 * cos(2*pi*n/(M-1)) + 0.08 * cos(4*pi*n/(M-1))
113///
114/// This is equivalent to `numpy.blackman(M)`.
115///
116/// # Errors
117/// Returns an error only if internal array construction fails.
118pub fn blackman(m: usize) -> FerrayResult<Array<f64, Ix1>> {
119    let denom = (m.saturating_sub(1)) as f64;
120    gen_window(m, |n| {
121        let x = n as f64;
122        0.08f64.mul_add(
123            (4.0 * PI * x / denom).cos(),
124            0.5f64.mul_add(-(2.0 * PI * x / denom).cos(), 0.42),
125        )
126    })
127}
128
129/// Return the Hamming window of length `m`.
130///
131/// The Hamming window is defined as:
132///   w(n) = 0.54 - 0.46 * cos(2*pi*n/(M-1))
133///
134/// This is equivalent to `numpy.hamming(M)`.
135///
136/// # Errors
137/// Returns an error only if internal array construction fails.
138pub fn hamming(m: usize) -> FerrayResult<Array<f64, Ix1>> {
139    let denom = (m.saturating_sub(1)) as f64;
140    gen_window(m, |n| {
141        0.46f64.mul_add(-(2.0 * PI * n as f64 / denom).cos(), 0.54)
142    })
143}
144
145/// Return the Hann (Hanning) window of length `m`.
146///
147/// The Hann window is defined as:
148///   w(n) = 0.5 * (1 - cos(2*pi*n/(M-1)))
149///
150/// `NumPy` calls this function `hanning`. This is equivalent to `numpy.hanning(M)`.
151///
152/// # Errors
153/// Returns an error only if internal array construction fails.
154pub fn hanning(m: usize) -> FerrayResult<Array<f64, Ix1>> {
155    let denom = (m.saturating_sub(1)) as f64;
156    gen_window(m, |n| 0.5 * (1.0 - (2.0 * PI * n as f64 / denom).cos()))
157}
158
159/// Return the Kaiser window of length `m` with shape parameter `beta`.
160///
161/// The Kaiser window is defined as:
162///   w(n) = `I_0(beta` * sqrt(1 - ((2n/(M-1)) - 1)^2)) / `I_0(beta)`
163///
164/// where `I_0` is the modified Bessel function of the first kind, order 0.
165///
166/// This is equivalent to `numpy.kaiser(M, beta)`.
167///
168/// # Edge Cases
169/// - `m == 0`: returns an empty array.
170/// - `m == 1`: returns `[1.0]`.
171///
172/// # Errors
173/// Returns `FerrayError::InvalidValue` if `beta` is NaN or if `beta`
174/// exceeds the range where `I_0(beta)` can be computed in f64. The
175/// asymptotic branch uses `exp(beta) / sqrt(beta)` internally, so
176/// `|beta| > ~709` overflows to infinity and the normalized kaiser
177/// window would reduce to `Inf / Inf = NaN` for every sample
178/// (#294).
179pub fn kaiser(m: usize, beta: f64) -> FerrayResult<Array<f64, Ix1>> {
180    if beta.is_nan() {
181        return Err(FerrayError::invalid_value("kaiser: beta must not be NaN"));
182    }
183    // I_0 is an even function, so kaiser(m, -beta) == kaiser(m, beta).
184    // Accept negative beta for NumPy compatibility.
185    let beta = beta.abs();
186    // Reject beta values where `exp(beta)` overflows f64 (ln(f64::MAX)
187    // ≈ 709.78). We use a slightly conservative threshold so the
188    // clipped output stays finite even after the asymptotic-series
189    // polynomial corrections.
190    const BETA_OVERFLOW_THRESHOLD: f64 = 708.0;
191    if beta > BETA_OVERFLOW_THRESHOLD {
192        return Err(FerrayError::invalid_value(format!(
193            "kaiser: |beta| = {beta} exceeds the safe range ({BETA_OVERFLOW_THRESHOLD}) \
194             for f64 I_0; the window would be NaN everywhere"
195        )));
196    }
197    let i0_beta = bessel_i0_scalar::<f64>(beta);
198    let alpha = (m.saturating_sub(1)) as f64 / 2.0;
199    gen_window(m, |n| {
200        let t = (n as f64 - alpha) / alpha;
201        let arg = beta * (1.0 - t * t).max(0.0).sqrt();
202        bessel_i0_scalar::<f64>(arg) / i0_beta
203    })
204}
205
206// ===========================================================================
207// SciPy-extended windows (cosine, exponential, gaussian, general_cosine,
208// general_hamming, nuttall, parzen, taylor, tukey).
209//
210// All are symmetric ("sym=True" in scipy parlance); a non-symmetric "periodic"
211// variant can be obtained by computing window of length m+1 and dropping the
212// last sample (NumPy convention).
213// ===========================================================================
214
215/// Half-cycle cosine window (also known as the "sine" window).
216///
217/// `w(n) = sin(pi * (n + 0.5) / M)` for `n = 0..M-1`.
218///
219/// Mirrors `scipy.signal.windows.cosine` /
220/// `torch.signal.windows.cosine`.
221pub fn cosine(m: usize) -> FerrayResult<Array<f64, Ix1>> {
222    if m == 0 {
223        return Array::from_vec(Ix1::new([0]), vec![]);
224    }
225    if m == 1 {
226        return Array::from_vec(Ix1::new([1]), vec![1.0]);
227    }
228    let mf = m as f64;
229    gen_window(m, |n| (PI * (n as f64 + 0.5) / mf).sin())
230}
231
232/// Exponentially decaying window centred on `center` with time-constant `tau`.
233///
234/// `w(n) = exp(-|n - center| / tau)`. If `center` is `None`, defaults to
235/// the geometric centre `(M - 1) / 2`. `tau` must be > 0.
236///
237/// Mirrors `scipy.signal.windows.exponential` /
238/// `torch.signal.windows.exponential`.
239pub fn exponential(m: usize, center: Option<f64>, tau: f64) -> FerrayResult<Array<f64, Ix1>> {
240    if !tau.is_finite() || tau <= 0.0 {
241        return Err(FerrayError::invalid_value(format!(
242            "exponential: tau must be positive and finite, got {tau}"
243        )));
244    }
245    let centre = center.unwrap_or((m.saturating_sub(1)) as f64 / 2.0);
246    gen_window(m, |n| (-((n as f64) - centre).abs() / tau).exp())
247}
248
249/// Gaussian window with standard deviation `std`.
250///
251/// `w(n) = exp(-((n - (M-1)/2) / std)^2 / 2)`. `std` must be > 0.
252///
253/// Mirrors `scipy.signal.windows.gaussian` /
254/// `torch.signal.windows.gaussian`.
255pub fn gaussian(m: usize, std: f64) -> FerrayResult<Array<f64, Ix1>> {
256    if !std.is_finite() || std <= 0.0 {
257        return Err(FerrayError::invalid_value(format!(
258            "gaussian: std must be positive and finite, got {std}"
259        )));
260    }
261    let centre = (m.saturating_sub(1)) as f64 / 2.0;
262    gen_window(m, |n| {
263        let z = ((n as f64) - centre) / std;
264        (-0.5 * z * z).exp()
265    })
266}
267
268/// Generalized cosine sum window: `w(n) = Σ_k (-1)^k a[k] cos(2π k n / (M-1))`.
269///
270/// `coeffs` is the slice `[a_0, a_1, a_2, ...]`. Setting `[0.5, 0.5]` gives
271/// the Hann window; `[0.42, 0.5, 0.08]` gives the classical Blackman.
272///
273/// Mirrors `scipy.signal.windows.general_cosine` /
274/// `torch.signal.windows.general_cosine`.
275pub fn general_cosine(m: usize, coeffs: &[f64]) -> FerrayResult<Array<f64, Ix1>> {
276    if coeffs.is_empty() {
277        return Err(FerrayError::invalid_value(
278            "general_cosine: coeffs must not be empty",
279        ));
280    }
281    let denom = (m.saturating_sub(1)) as f64;
282    let coeffs = coeffs.to_vec();
283    gen_window(m, |n| {
284        let nf = n as f64;
285        let mut sum = 0.0_f64;
286        for (k, &a) in coeffs.iter().enumerate() {
287            let sign = if k % 2 == 0 { 1.0 } else { -1.0 };
288            sum += sign * a * (2.0 * PI * (k as f64) * nf / denom).cos();
289        }
290        sum
291    })
292}
293
294/// Generalized Hamming window: `w(n) = α - (1 - α) cos(2π n / (M - 1))`.
295///
296/// `alpha = 0.5` gives the Hann window; `alpha = 25/46 ≈ 0.5435` gives a
297/// "perfect" Hamming (NumPy uses 0.54 by tradition).
298///
299/// Mirrors `scipy.signal.windows.general_hamming` /
300/// `torch.signal.windows.general_hamming`.
301pub fn general_hamming(m: usize, alpha: f64) -> FerrayResult<Array<f64, Ix1>> {
302    if !alpha.is_finite() {
303        return Err(FerrayError::invalid_value(format!(
304            "general_hamming: alpha must be finite, got {alpha}"
305        )));
306    }
307    let denom = (m.saturating_sub(1)) as f64;
308    gen_window(m, |n| {
309        alpha - (1.0 - alpha) * (2.0 * PI * (n as f64) / denom).cos()
310    })
311}
312
313/// 4-term Nuttall window with continuous first derivative.
314///
315/// Coefficients [0.3635819, 0.4891775, 0.1365995, 0.0106411] (Nuttall 1981,
316/// "minimum 4-term Blackman-Harris with continuous first derivative").
317///
318/// Mirrors `scipy.signal.windows.nuttall` /
319/// `torch.signal.windows.nuttall`.
320pub fn nuttall(m: usize) -> FerrayResult<Array<f64, Ix1>> {
321    general_cosine(m, &[0.3635819, 0.4891775, 0.1365995, 0.0106411])
322}
323
324/// Parzen (de la Vallée Poussin) window: a piecewise-cubic B-spline.
325///
326/// Mirrors `scipy.signal.windows.parzen` / `torch.signal.windows.parzen`.
327pub fn parzen(m: usize) -> FerrayResult<Array<f64, Ix1>> {
328    if m == 0 {
329        return Array::from_vec(Ix1::new([0]), vec![]);
330    }
331    if m == 1 {
332        return Array::from_vec(Ix1::new([1]), vec![1.0]);
333    }
334    // Parzen formula: |x| in [0, M/4]: 1 - 6 r^2 + 6 |r|^3,
335    //                |x| in (M/4, M/2]: 2 (1 - |r|)^3, where r = x / (M/2).
336    let half = m as f64 / 2.0;
337    gen_window(m, |n| {
338        let x = (n as f64) - (m as f64 - 1.0) / 2.0;
339        let r = x.abs() / half;
340        if r <= 0.5 {
341            1.0 - 6.0 * r * r + 6.0 * r * r * r
342        } else if r <= 1.0 {
343            let one_minus_r = 1.0 - r;
344            2.0 * one_minus_r * one_minus_r * one_minus_r
345        } else {
346            0.0
347        }
348    })
349}
350
351/// Dolph–Chebyshev window with sidelobe attenuation `at` dB (#740).
352///
353/// Mirrors `scipy.signal.windows.chebwin`. The window is the
354/// inverse Fourier transform of the (M-1)-th Chebyshev polynomial
355/// of the first kind evaluated on a cosine sweep, producing an
356/// equiripple sidelobe response at `at` dB below the main lobe.
357///
358/// Implemented via direct DFT (sufficient for typical M ≤ 1024;
359/// the inner loop is O(M²) but allocation-free). Result is
360/// normalised so the centre value is 1.
361///
362/// # Errors
363/// `FerrayError::InvalidValue` if `at` is non-finite or
364/// `at <= 0.0`. `m == 0` returns the empty window; `m == 1`
365/// returns `[1.0]`.
366pub fn chebwin(m: usize, at: f64) -> FerrayResult<Array<f64, Ix1>> {
367    if !at.is_finite() || at <= 0.0 {
368        return Err(FerrayError::invalid_value(format!(
369            "chebwin: at must be a positive finite dB attenuation, got {at}"
370        )));
371    }
372    if m == 0 {
373        return Array::from_vec(Ix1::new([0]), vec![]);
374    }
375    if m == 1 {
376        return Array::from_vec(Ix1::new([1]), vec![1.0]);
377    }
378
379    let order = (m - 1) as f64;
380    let r = 10.0_f64.powf(at / 20.0);
381    let beta = (r.acosh() / order).cosh();
382
383    // Build spectrum P[k] = T_{m-1}(β cos(πk/m)).
384    let mut spectrum = Vec::with_capacity(m);
385    for k in 0..m {
386        let x = beta * (PI * k as f64 / m as f64).cos();
387        spectrum.push(chebyshev_t_extended(x, order));
388    }
389
390    // Direct DFT to time domain. For odd M the spectrum is even-
391    // symmetric about M/2 and the IDFT is purely real; for even M
392    // the (M-1)-th Chebyshev polynomial is odd, so we apply a
393    // half-bin phase shift before transforming. Matches scipy's
394    // FFT-based pipeline up to the trailing peak normalisation.
395    // Forward DFT of the spectrum (matches scipy's fft(p) call).
396    // For odd M the spectrum is even-symmetric so fft is real;
397    // for even M scipy multiplies by exp(iπk/M) phase first to
398    // recover symmetry. After the DFT we keep the first half,
399    // normalise by w[0] so the eventual centre lands at 1, and
400    // mirror to recover the full symmetric window.
401    let mf = m as f64;
402    let half_dft = |idx: usize| -> f64 {
403        if m % 2 == 1 {
404            spectrum
405                .iter()
406                .enumerate()
407                .map(|(k, &pk)| pk * (-2.0 * PI * k as f64 * idx as f64 / mf).cos())
408                .sum()
409        } else {
410            // Forward DFT with the additional exp(iπk/M) phase
411            // shift baked in: p[k] * exp(iπk/M) * exp(-i2πkn/M)
412            //              = p[k] * cos(πk(1 - 2n)/M).
413            spectrum
414                .iter()
415                .enumerate()
416                .map(|(k, &pk)| {
417                    let phase = PI * k as f64 * (1.0 - 2.0 * idx as f64) / mf;
418                    pk * phase.cos()
419                })
420                .sum()
421        }
422    };
423
424    let half_len = m / 2 + 1;
425    let mut half: Vec<f64> = (0..half_len).map(half_dft).collect();
426    let denom = if m % 2 == 1 { half[0] } else { half[1] };
427    if denom != 0.0 {
428        for v in &mut half {
429            *v /= denom;
430        }
431    }
432
433    let mut w = Vec::with_capacity(m);
434    if m % 2 == 1 {
435        // Mirror: reverse(half[1..]) ++ half
436        for &v in half.iter().skip(1).rev() {
437            w.push(v);
438        }
439        w.extend_from_slice(&half);
440    } else {
441        // Even: scipy's `concatenate((w[n-1:0:-1], w[1:n]))` where
442        // n = M/2 + 1; reverse of half[1..n], then half[1..n].
443        let n = m / 2 + 1;
444        for &v in half[1..n].iter().rev() {
445            w.push(v);
446        }
447        w.extend_from_slice(&half[1..n]);
448    }
449
450    Array::from_vec(Ix1::new([m]), w)
451}
452
453/// Discrete prolate spheroidal sequence (Slepian taper) (#745).
454///
455/// Returns the dominant DPSS — the eigenvector of the canonical
456/// time–frequency concentration tridiagonal matrix corresponding
457/// to the largest eigenvalue. Equivalent to
458/// `scipy.signal.windows.dpss(m, NW)` with the default `Kmax=None`
459/// (single window, sym=True, no ratios).
460///
461/// `nw` is the standardised half bandwidth (`NW = M·W`), typically
462/// in the range 2..=4 for spectral analysis.
463///
464/// Computed by warm-up power iteration followed by Rayleigh-quotient
465/// inverse iteration (RQI) on the symmetric tridiagonal
466/// concentration matrix
467///
468/// ```text
469///   d[i]  = ((M-1)/2 − i)² · cos(2π · NW / M)
470///   e[i]  = (i+1)(M−i−1) / 2
471/// ```
472///
473/// Power iteration alone converges only linearly with rate
474/// `λ₂/λ₁`, which is far too slow for the nearly-degenerate
475/// dominant DPSS eigenvalues; RQI achieves cubic local convergence
476/// near the dominant eigenpair via tridiagonal LU (Thomas
477/// algorithm) solves of `(T - σI) w = v`. The result is
478/// sign-flipped if needed so the centre value is positive, and
479/// peak-normalised (matching scipy's default `norm='approximate'`).
480///
481/// Multiple Slepian sequences (`Kmax > 1`) need orthogonal
482/// power-iteration / deflation; tracked as a future extension.
483///
484/// # Errors
485/// `FerrayError::InvalidValue` if `nw` is non-finite, `nw <= 0`,
486/// or `nw >= m / 2` (the band would degenerate).
487pub fn dpss(m: usize, nw: f64) -> FerrayResult<Array<f64, Ix1>> {
488    if !nw.is_finite() || nw <= 0.0 {
489        return Err(FerrayError::invalid_value(format!(
490            "dpss: nw must be a positive finite half bandwidth, got {nw}"
491        )));
492    }
493    if m == 0 {
494        return Array::from_vec(Ix1::new([0]), vec![]);
495    }
496    if m == 1 {
497        return Array::from_vec(Ix1::new([1]), vec![1.0]);
498    }
499    if nw >= (m as f64) / 2.0 {
500        return Err(FerrayError::invalid_value(format!(
501            "dpss: nw = {nw} must be < M/2 = {} for a non-degenerate band",
502            (m as f64) / 2.0
503        )));
504    }
505
506    let mf = m as f64;
507    let cos_w = (2.0 * PI * nw / mf).cos();
508    // Diagonal d[i] = ((M-1)/2 - i)^2 * cos(2π NW / M).
509    let diag: Vec<f64> = (0..m)
510        .map(|i| {
511            let centre = (mf - 1.0) / 2.0 - i as f64;
512            centre * centre * cos_w
513        })
514        .collect();
515    // Off-diagonal e[i] = (i+1)(M-i-1)/2 for i in 0..m-1.
516    let off: Vec<f64> = (0..m.saturating_sub(1))
517        .map(|i| (i as f64 + 1.0) * (mf - i as f64 - 1.0) / 2.0)
518        .collect();
519
520    // Power iteration. Centre-symmetric initial vector accelerates
521    // convergence to the (also symmetric) dominant DPSS.
522    let mut v: Vec<f64> = (0..m)
523        .map(|i| {
524            let centred = (i as f64 - (mf - 1.0) / 2.0) / mf;
525            (-2.0 * centred * centred).exp()
526        })
527        .collect();
528    let n0 = v.iter().map(|x| x * x).sum::<f64>().sqrt();
529    for x in &mut v {
530        *x /= n0.max(f64::MIN_POSITIVE);
531    }
532
533    // Helper: y = T v for the symmetric tridiagonal T.
534    let apply_t = |x: &[f64]| -> Vec<f64> {
535        let mut y = vec![0.0_f64; m];
536        for i in 0..m {
537            y[i] = diag[i] * x[i];
538            if i > 0 {
539                y[i] += off[i - 1] * x[i - 1];
540            }
541            if i + 1 < m {
542                y[i] += off[i] * x[i + 1];
543            }
544        }
545        y
546    };
547
548    // Stage 1: warm-up power iteration (~50 steps) to land near the
549    // dominant eigenspace before switching to RQI.
550    for _ in 0..50 {
551        let mut next = apply_t(&v);
552        let norm = next.iter().map(|x| x * x).sum::<f64>().sqrt();
553        if norm < f64::MIN_POSITIVE {
554            return Err(FerrayError::invalid_value(
555                "dpss: power iteration produced a zero vector",
556            ));
557        }
558        for x in &mut next {
559            *x /= norm;
560        }
561        v = next;
562    }
563
564    // Initial eigenvalue estimate σ = v^T (T v).
565    let tv0 = apply_t(&v);
566    let mut sigma: f64 = v.iter().zip(tv0.iter()).map(|(a, b)| a * b).sum();
567
568    // Stage 2: Rayleigh-quotient inverse iteration. Each step solves
569    // (T - σI) w = v via the Thomas algorithm (tridiagonal LU with
570    // no pivoting), normalises w, and updates σ = w^T (T w). Cubic
571    // local convergence near the dominant eigenpair.
572    for _ in 0..50 {
573        // Build the shifted tridiagonal: lower = upper = off (symmetric).
574        let mut d_shift: Vec<f64> = diag.iter().map(|d| d - sigma).collect();
575        let sup: Vec<f64> = off.clone();
576        let mut rhs: Vec<f64> = v.clone();
577
578        // Forward sweep: eliminate sub-diagonal.
579        let mut singular = false;
580        for i in 1..m {
581            if d_shift[i - 1].abs() < 1e-30 {
582                // σ is essentially an eigenvalue — v is already the
583                // eigenvector to machine precision; stop here.
584                singular = true;
585                break;
586            }
587            let factor = off[i - 1] / d_shift[i - 1];
588            d_shift[i] -= factor * sup[i - 1];
589            rhs[i] -= factor * rhs[i - 1];
590        }
591        if singular {
592            break;
593        }
594        if d_shift[m - 1].abs() < 1e-30 {
595            break;
596        }
597
598        // Back sweep.
599        let mut w = vec![0.0_f64; m];
600        w[m - 1] = rhs[m - 1] / d_shift[m - 1];
601        for i in (0..m - 1).rev() {
602            w[i] = (rhs[i] - sup[i] * w[i + 1]) / d_shift[i];
603        }
604
605        let norm = w.iter().map(|x| x * x).sum::<f64>().sqrt();
606        if norm < f64::MIN_POSITIVE {
607            break;
608        }
609        for x in &mut w {
610            *x /= norm;
611        }
612
613        let tw = apply_t(&w);
614        let new_sigma: f64 = w.iter().zip(tw.iter()).map(|(a, b)| a * b).sum();
615
616        // Residual ||T w - σ w||.
617        let resid: f64 = tw
618            .iter()
619            .zip(w.iter())
620            .map(|(t, x)| {
621                let r = t - new_sigma * x;
622                r * r
623            })
624            .sum::<f64>()
625            .sqrt();
626
627        v = w;
628        sigma = new_sigma;
629        if resid < 1e-14 {
630            break;
631        }
632    }
633    // Silence unused-assignment lint: sigma is the final eigenvalue
634    // estimate, kept for clarity even though it is not read after
635    // the loop.
636    let _ = sigma;
637
638    // Sign convention: scipy returns the dominant DPSS with a
639    // positive centre value.
640    let centre_val = if m % 2 == 1 {
641        v[m / 2]
642    } else {
643        v[m / 2 - 1] + v[m / 2]
644    };
645    if centre_val < 0.0 {
646        for x in &mut v {
647            *x = -*x;
648        }
649    }
650
651    // Match scipy's default `norm='approximate'`: rescale so the
652    // peak (after sign-fix, this is the centre value) is 1. Power
653    // iteration left v at unit ℓ²; switch to max-normalised so
654    // numeric values match scipy.signal.windows.dpss directly.
655    let peak = v.iter().copied().fold(f64::NEG_INFINITY, f64::max);
656    if peak > 0.0 {
657        for x in &mut v {
658            *x /= peak;
659        }
660    }
661
662    Array::from_vec(Ix1::new([m]), v)
663}
664
665/// Chebyshev polynomial of the first kind T_n(x) extended to all
666/// real x via cosh — `T_n(cos θ) = cos(nθ)` for `|x| ≤ 1`,
667/// `T_n(cosh θ) = cosh(nθ)` for `x > 1`. Sign handling on the
668/// negative branch uses `(-1)^n`.
669fn chebyshev_t_extended(x: f64, n: f64) -> f64 {
670    if x.abs() <= 1.0 {
671        (n * x.acos()).cos()
672    } else {
673        let mag = (n * x.abs().acosh()).cosh();
674        if x < 0.0 && (n as i64) % 2 != 0 {
675            -mag
676        } else {
677            mag
678        }
679    }
680}
681
682/// Boxcar (rectangular) window of length `m` (#738).
683///
684/// All-ones window. Equivalent to `scipy.signal.windows.boxcar`.
685///
686/// # Errors
687/// Only on internal array construction failure.
688pub fn boxcar(m: usize) -> FerrayResult<Array<f64, Ix1>> {
689    gen_window(m, |_| 1.0)
690}
691
692/// Triangular window of length `m` (#738).
693///
694/// Centre = 1, linearly decreasing to (1/m) at the edges (slightly
695/// different from Bartlett, which decays to 0 at the edges).
696/// Equivalent to `scipy.signal.windows.triang`.
697///
698/// # Errors
699/// Only on internal array construction failure.
700pub fn triang(m: usize) -> FerrayResult<Array<f64, Ix1>> {
701    if m == 0 {
702        return Array::from_vec(Ix1::new([0]), vec![]);
703    }
704    if m == 1 {
705        return Array::from_vec(Ix1::new([1]), vec![1.0]);
706    }
707    let mf = m as f64;
708    let centre = (mf - 1.0) / 2.0;
709    if m % 2 == 0 {
710        // Even-m formula: w[n] = 1 - |2n - (m-1)| / m
711        gen_window(m, |n| 1.0 - (2.0 * n as f64 - centre - centre).abs() / mf)
712    } else {
713        // Odd-m formula: w[n] = 1 - |2n - (m-1)| / (m + 1)
714        gen_window(m, |n| {
715            1.0 - (2.0 * n as f64 - centre - centre).abs() / (mf + 1.0)
716        })
717    }
718}
719
720/// Bohman window of length `m` (#738).
721///
722/// w(n) = (1 - |r|) cos(π|r|) + (1/π) sin(π|r|), where r is the
723/// normalised position in [-1, 1]. Equivalent to
724/// `scipy.signal.windows.bohman`.
725///
726/// # Errors
727/// Only on internal array construction failure.
728pub fn bohman(m: usize) -> FerrayResult<Array<f64, Ix1>> {
729    if m == 0 {
730        return Array::from_vec(Ix1::new([0]), vec![]);
731    }
732    if m == 1 {
733        return Array::from_vec(Ix1::new([1]), vec![1.0]);
734    }
735    // scipy: linspace(-1, 1, m)[1:-1] core, with explicit zeros at endpoints.
736    let mf = (m - 1) as f64;
737    let mut data = Vec::with_capacity(m);
738    for n in 0..m {
739        if n == 0 || n == m - 1 {
740            data.push(0.0);
741            continue;
742        }
743        let r = (2.0 * n as f64 / mf - 1.0).abs();
744        let v = (1.0 - r) * (PI * r).cos() + (PI * r).sin() / PI;
745        data.push(v);
746    }
747    Array::from_vec(Ix1::new([m]), data)
748}
749
750/// Flat-top window of length `m` (#738).
751///
752/// 5-term cosine window with coefficients chosen for minimum
753/// passband ripple (used for amplitude calibration). Equivalent to
754/// `scipy.signal.windows.flattop`.
755///
756/// # Errors
757/// Only on internal array construction failure.
758pub fn flattop(m: usize) -> FerrayResult<Array<f64, Ix1>> {
759    // scipy.signal.windows.flattop coefficients (a0..a4), verbatim from
760    // scipy 1.13.0 signal/windows/_windows.py.
761    const COEFFS: [f64; 5] = [
762        0.215_578_95,
763        0.416_631_58,
764        0.277_263_158,
765        0.083_578_947,
766        0.006_947_368,
767    ];
768    if m == 0 {
769        return Array::from_vec(Ix1::new([0]), vec![]);
770    }
771    let denom = (m - 1).max(1) as f64;
772    gen_window(m, |n| {
773        let phase = 2.0 * PI * n as f64 / denom;
774        let mut acc = COEFFS[0];
775        for (k, &c) in COEFFS.iter().enumerate().skip(1) {
776            let sign = if k % 2 == 0 { 1.0 } else { -1.0 };
777            acc += sign * c * (k as f64 * phase).cos();
778        }
779        acc
780    })
781}
782
783/// Lanczos (sinc) window of length `m` (#738).
784///
785/// w(n) = sinc(2 n / (m-1) - 1). Equivalent to
786/// `scipy.signal.windows.lanczos`.
787///
788/// # Errors
789/// Only on internal array construction failure.
790pub fn lanczos(m: usize) -> FerrayResult<Array<f64, Ix1>> {
791    if m == 0 {
792        return Array::from_vec(Ix1::new([0]), vec![]);
793    }
794    if m == 1 {
795        return Array::from_vec(Ix1::new([1]), vec![1.0]);
796    }
797    let denom = (m - 1) as f64;
798    gen_window(m, |n| {
799        let x = 2.0 * n as f64 / denom - 1.0;
800        if x.abs() < f64::EPSILON {
801            1.0
802        } else {
803            let pi_x = PI * x;
804            pi_x.sin() / pi_x
805        }
806    })
807}
808
809/// Taylor window with `nbar` near-side lobes and a sidelobe level of `sll`
810/// dB below the main lobe.
811///
812/// `nbar` (typically 4) controls how many sidelobes are constrained at the
813/// design level; `sll` (typically 30) is the desired peak sidelobe
814/// attenuation in **positive** dB. When `norm` is true the window is
815/// normalised so `w[(M-1)/2] = 1`.
816///
817/// Implementation follows Carrara & Goodman, "Symmetric Taylor Window"
818/// (Synthetic Aperture Radar, 1995). For the radar/array-processing
819/// applications where Taylor windows are used, `nbar = 4`, `sll = 30`
820/// gives equiripple sidelobes ~30 dB down — the usual default.
821///
822/// Mirrors `scipy.signal.windows.taylor` / `torch.signal.windows.taylor`.
823pub fn taylor(m: usize, nbar: usize, sll: f64, norm: bool) -> FerrayResult<Array<f64, Ix1>> {
824    if m == 0 {
825        return Array::from_vec(Ix1::new([0]), vec![]);
826    }
827    if m == 1 {
828        return Array::from_vec(Ix1::new([1]), vec![1.0]);
829    }
830    if nbar == 0 {
831        return Err(FerrayError::invalid_value("taylor: nbar must be >= 1"));
832    }
833    if !sll.is_finite() {
834        return Err(FerrayError::invalid_value("taylor: sll must be finite"));
835    }
836    // R = 10^(sll/20), B = (1/π) acosh(R)
837    let r = 10.0_f64.powf(sll / 20.0);
838    let b = r.acosh() / PI;
839    let nbar_f = nbar as f64;
840    // sigma^2 chosen so the (nbar)-th zero of the Taylor pattern is at
841    // n = nbar (Carrara & Goodman eq. 13).
842    let sigma2 = (nbar_f * nbar_f) / (b * b + (nbar_f - 0.5) * (nbar_f - 0.5));
843
844    // Compute coefficients F_m for m = 1..nbar-1.
845    let mut f_coeffs = Vec::with_capacity(nbar.saturating_sub(1));
846    for mm in 1..nbar {
847        let mmf = mm as f64;
848        // Numerator: ((-1)^(m+1) / 2) * Π_{n=1..nbar-1} [1 - m^2 / (sigma^2 * (B^2 + (n - 0.5)^2))]
849        // Denominator: Π_{n=1..nbar-1, n != m} [1 - m^2 / n^2]
850        let mut num = 1.0_f64;
851        for n in 1..nbar {
852            let nf = n as f64;
853            num *= 1.0 - mmf * mmf / (sigma2 * (b * b + (nf - 0.5) * (nf - 0.5)));
854        }
855        let sign = if mm % 2 == 0 { -1.0 } else { 1.0 };
856        let mut den = 1.0_f64;
857        for n in 1..nbar {
858            if n == mm {
859                continue;
860            }
861            let nf = n as f64;
862            den *= 1.0 - mmf * mmf / (nf * nf);
863        }
864        // The 0.5 prefactor: F_0 = 1 contributes the constant term, and
865        // each F_m doubles when reflected about zero in the cosine sum,
866        // so we halve the inverse-Fourier coefficients here.
867        f_coeffs.push(0.5 * sign * num / den);
868    }
869
870    let denom = (m - 1) as f64;
871    let arr = gen_window(m, |n| {
872        let xn = (n as f64) - denom / 2.0;
873        let mut w = 1.0_f64;
874        for (idx, &fk) in f_coeffs.iter().enumerate() {
875            let kk = (idx + 1) as f64;
876            w += 2.0 * fk * (2.0 * PI * kk * xn / m as f64).cos();
877        }
878        w
879    })?;
880
881    if !norm {
882        return Ok(arr);
883    }
884    // Analytic peak of the Taylor cosine sum at the fractional midpoint
885    // (m-1)/2. ferray's xn = n - (m-1)/2 means the midpoint is xn = 0,
886    // so every cosine collapses to 1 and the sum reduces to
887    // 1 + 2·Σ f_coeffs. This matches scipy's `scale = 1.0 / W((M-1)/2)`
888    // analytically, and avoids the sample-averaging error that affected
889    // even-M windows (#810 upstream from ferrotorch).
890    let centre_val: f64 = 1.0 + 2.0 * f_coeffs.iter().sum::<f64>();
891    if centre_val == 0.0 {
892        return Ok(arr); // pathological; leave un-normalised
893    }
894    let normed: Vec<f64> = arr
895        .as_slice()
896        .unwrap()
897        .iter()
898        .map(|&v| v / centre_val)
899        .collect();
900    Array::from_vec(Ix1::new([m]), normed)
901}
902
903/// Tukey (cosine-tapered) window with taper ratio `alpha` ∈ [0, 1].
904///
905/// `alpha = 0` gives a rectangular window; `alpha = 1` gives a Hann
906/// window. The middle `(1 - alpha) * (M - 1)` samples are unity, and the
907/// edges are tapered with a half-cosine.
908///
909/// Mirrors `scipy.signal.windows.tukey` / `torch.signal.windows.tukey`.
910pub fn tukey(m: usize, alpha: f64) -> FerrayResult<Array<f64, Ix1>> {
911    if !alpha.is_finite() || !(0.0..=1.0).contains(&alpha) {
912        return Err(FerrayError::invalid_value(format!(
913            "tukey: alpha must be in [0, 1], got {alpha}"
914        )));
915    }
916    if m == 0 {
917        return Array::from_vec(Ix1::new([0]), vec![]);
918    }
919    if m == 1 {
920        return Array::from_vec(Ix1::new([1]), vec![1.0]);
921    }
922    if alpha == 0.0 {
923        return Array::from_vec(Ix1::new([m]), vec![1.0; m]);
924    }
925    let denom = (m - 1) as f64;
926    let width = alpha * denom / 2.0;
927    gen_window(m, |n| {
928        let nf = n as f64;
929        if nf < width {
930            0.5 * (1.0 + (PI * (nf / width - 1.0)).cos())
931        } else if nf <= denom - width {
932            1.0
933        } else {
934            0.5 * (1.0 + (PI * ((denom - nf) / width - 1.0)).cos())
935        }
936    })
937}
938
939#[cfg(test)]
940mod tests {
941    use super::*;
942
943    // ---- dpss (#745) ---------------------------------------------------
944
945    #[test]
946    fn dpss_max_one_and_symmetric() {
947        let w = dpss(64, 3.0).unwrap();
948        let s = w.as_slice().unwrap();
949        // scipy's default `norm='approximate'`: max value is 1.
950        let max = s.iter().copied().fold(f64::NEG_INFINITY, f64::max);
951        assert!((max - 1.0).abs() < 1e-12, "max = {max}");
952        // Symmetric about the centre.
953        for i in 0..s.len() / 2 {
954            let mirror = s.len() - 1 - i;
955            assert!(
956                (s[i] - s[mirror]).abs() < 1e-7,
957                "asymmetric at i={i}: {} vs {}",
958                s[i],
959                s[mirror]
960            );
961        }
962    }
963
964    #[test]
965    fn dpss_centre_is_positive_and_peak() {
966        // Sign convention: centre value > 0 and is the maximum
967        // across the window (now exactly 1 after scipy-style
968        // max-normalisation).
969        let w = dpss(48, 2.5).unwrap();
970        let s = w.as_slice().unwrap();
971        let mid = s.len() / 2;
972        assert!(s[mid] > 0.0);
973        assert!((s[mid] - 1.0).abs() < 1e-12);
974    }
975
976    #[test]
977    fn dpss_endpoints_smaller_than_centre() {
978        // For NW=3 the dominant DPSS should be heavily concentrated
979        // near the centre — endpoints decay to small magnitudes.
980        let w = dpss(64, 3.0).unwrap();
981        let s = w.as_slice().unwrap();
982        let mid = s.len() / 2;
983        assert!(s[mid] > s[0]);
984        assert!(s[mid] > s[s.len() - 1]);
985        assert!(s[0].abs() < 0.1 * s[mid]);
986    }
987
988    #[test]
989    fn dpss_rejects_invalid_nw() {
990        assert!(dpss(32, -1.0).is_err());
991        assert!(dpss(32, 0.0).is_err());
992        assert!(dpss(32, f64::NAN).is_err());
993        // nw >= M/2 is degenerate.
994        assert!(dpss(8, 4.0).is_err());
995        assert!(dpss(8, 5.0).is_err());
996    }
997
998    #[test]
999    fn dpss_matches_scipy_within_relative_tolerance() {
1000        // scipy.signal.windows.dpss(16, 3.0) reference values for
1001        // the first half (window is symmetric, max-normalised).
1002        // Power iteration converges to within ~1% of scipy's
1003        // LAPACK-based eigensolver result for this small M; the
1004        // shape and concentration profile are identical up to a
1005        // small numerical drift.
1006        let want_first_half: [f64; 8] = [
1007            0.00710856, 0.03620714, 0.10884379, 0.24276927, 0.43733689, 0.6634221, 0.86701952,
1008            0.98841699,
1009        ];
1010        let got = dpss(16, 3.0).unwrap();
1011        let s = got.as_slice().unwrap();
1012        for (i, (&g, &w)) in s[..8].iter().zip(want_first_half.iter()).enumerate() {
1013            // Use generous relative tolerance — the small endpoint
1014            // values are most sensitive to power-iteration
1015            // convergence accuracy. 5% relative covers all
1016            // samples; the larger central values stay within 0.1%.
1017            let tol = 0.05 * w.abs().max(1e-3);
1018            assert!((g - w).abs() <= tol, "i={i}: got={g}, want={w}, tol={tol}");
1019        }
1020    }
1021
1022    #[test]
1023    fn dpss_m0_m1_edge_cases() {
1024        let z = dpss(0, 2.0).unwrap();
1025        assert_eq!(z.shape(), &[0]);
1026        let one = dpss(1, 0.5).unwrap();
1027        assert_eq!(one.as_slice().unwrap(), &[1.0]);
1028    }
1029
1030    // ---- chebwin (#740) ------------------------------------------------
1031
1032    #[test]
1033    fn chebwin_centre_is_one_and_symmetric() {
1034        let w = chebwin(11, 50.0).unwrap();
1035        let s = w.as_slice().unwrap();
1036        let mid = s.len() / 2;
1037        assert!((s[mid] - 1.0).abs() < 1e-6, "centre = {}", s[mid]);
1038        for i in 0..s.len() / 2 {
1039            assert!(
1040                (s[i] - s[s.len() - 1 - i]).abs() < 1e-9,
1041                "asymmetric at i={i}: {} vs {}",
1042                s[i],
1043                s[s.len() - 1 - i]
1044            );
1045        }
1046    }
1047
1048    #[test]
1049    fn chebwin_even_length_symmetric() {
1050        let w = chebwin(8, 40.0).unwrap();
1051        let s = w.as_slice().unwrap();
1052        for i in 0..s.len() / 2 {
1053            assert!(
1054                (s[i] - s[s.len() - 1 - i]).abs() < 1e-9,
1055                "asymmetric at i={i}"
1056            );
1057        }
1058        // No element should exceed the centre by more than 1e-9.
1059        let max = s.iter().copied().fold(f64::NEG_INFINITY, f64::max);
1060        assert!((max - 1.0).abs() < 1e-6);
1061    }
1062
1063    #[test]
1064    fn chebwin_endpoints_smaller_than_centre() {
1065        // High attenuation → endpoints should be very small.
1066        let w = chebwin(31, 80.0).unwrap();
1067        let s = w.as_slice().unwrap();
1068        let mid = s.len() / 2;
1069        assert!(s[mid] > s[0]);
1070        assert!(s[mid] > s[s.len() - 1]);
1071    }
1072
1073    #[test]
1074    fn chebwin_rejects_invalid_attenuation() {
1075        assert!(chebwin(16, -10.0).is_err());
1076        assert!(chebwin(16, 0.0).is_err());
1077        assert!(chebwin(16, f64::NAN).is_err());
1078    }
1079
1080    #[test]
1081    fn chebwin_matches_scipy_known_values() {
1082        // scipy.signal.windows.chebwin(11, 50) reference values.
1083        let want = [
1084            0.06228583, 0.20113857, 0.42847525, 0.69573494, 0.91497506, 1.0, 0.91497506,
1085            0.69573494, 0.42847525, 0.20113857, 0.06228583,
1086        ];
1087        let got = chebwin(11, 50.0).unwrap();
1088        let s = got.as_slice().unwrap();
1089        for (i, (&g, &w)) in s.iter().zip(want.iter()).enumerate() {
1090            assert!((g - w).abs() < 1e-6, "i={i}: got={g}, want={w}");
1091        }
1092    }
1093
1094    #[test]
1095    fn chebwin_matches_scipy_even_length() {
1096        // scipy.signal.windows.chebwin(8, 40) reference values.
1097        let want = [
1098            0.14609713, 0.41790422, 0.75944595, 1.0, 1.0, 0.75944595, 0.41790422, 0.14609713,
1099        ];
1100        let got = chebwin(8, 40.0).unwrap();
1101        let s = got.as_slice().unwrap();
1102        for (i, (&g, &w)) in s.iter().zip(want.iter()).enumerate() {
1103            assert!((g - w).abs() < 1e-6, "i={i}: got={g}, want={w}");
1104        }
1105    }
1106
1107    #[test]
1108    fn chebwin_m0_m1_edge_cases() {
1109        let z = chebwin(0, 50.0).unwrap();
1110        assert_eq!(z.shape(), &[0]);
1111        let one = chebwin(1, 50.0).unwrap();
1112        assert_eq!(one.as_slice().unwrap(), &[1.0]);
1113    }
1114
1115    // ---- additional scipy.signal windows (#738) ------------------------
1116
1117    #[test]
1118    fn boxcar_all_ones() {
1119        let w = boxcar(8).unwrap();
1120        let s = w.as_slice().unwrap();
1121        for &v in s {
1122            assert!((v - 1.0).abs() < 1e-15);
1123        }
1124    }
1125
1126    #[test]
1127    fn triang_centre_is_max_and_endpoints_smaller() {
1128        let w = triang(7).unwrap();
1129        let s = w.as_slice().unwrap();
1130        // Symmetry: w[n] == w[m-1-n].
1131        for i in 0..s.len() / 2 {
1132            assert!((s[i] - s[s.len() - 1 - i]).abs() < 1e-12);
1133        }
1134        // Centre is the peak.
1135        let mid = s.len() / 2;
1136        for i in 0..s.len() {
1137            if i != mid {
1138                assert!(s[i] <= s[mid] + 1e-12);
1139            }
1140        }
1141    }
1142
1143    #[test]
1144    fn bohman_zero_at_endpoints_and_one_at_centre() {
1145        let w = bohman(5).unwrap();
1146        let s = w.as_slice().unwrap();
1147        assert!(s[0].abs() < 1e-12);
1148        assert!(s[s.len() - 1].abs() < 1e-12);
1149        // Centre value is the maximum and very close to 1.
1150        let mid = s.len() / 2;
1151        assert!((s[mid] - 1.0).abs() < 1e-6);
1152    }
1153
1154    #[test]
1155    fn flattop_centre_around_one_and_negative_lobes_present() {
1156        // flattop is designed for amplitude flatness — centre ≈ 1
1157        // with the canonical scaling, and side regions go negative.
1158        let w = flattop(11).unwrap();
1159        let s = w.as_slice().unwrap();
1160        let mid = s.len() / 2;
1161        assert!((s[mid] - 1.0).abs() < 0.01);
1162        // Some side samples must go below zero (window has negative
1163        // side lobes by design).
1164        let any_negative = s.iter().any(|&v| v < -0.001);
1165        assert!(any_negative, "flattop should have negative side lobes");
1166    }
1167
1168    #[test]
1169    fn lanczos_centre_is_one_and_endpoints_zero() {
1170        let w = lanczos(9).unwrap();
1171        let s = w.as_slice().unwrap();
1172        let mid = s.len() / 2;
1173        assert!((s[mid] - 1.0).abs() < 1e-12);
1174        assert!(s[0].abs() < 1e-12);
1175        assert!(s[s.len() - 1].abs() < 1e-12);
1176    }
1177
1178    #[test]
1179    fn small_m_edge_cases_handled() {
1180        // m = 0 returns empty, m = 1 returns [1.0] for symmetric windows.
1181        for f in [triang as fn(usize) -> _, bohman, lanczos] {
1182            let z = f(0).unwrap();
1183            assert_eq!(z.shape(), &[0]);
1184            let one = f(1).unwrap();
1185            assert_eq!(one.shape(), &[1]);
1186            assert_eq!(one.as_slice().unwrap(), &[1.0]);
1187        }
1188    }
1189
1190    // ---- f32 sibling functions (#529) ----------------------------------
1191
1192    #[test]
1193    fn f32_windows_match_f64_within_f32_eps() {
1194        for m in [1usize, 5, 16, 64] {
1195            for (name, got, want) in [
1196                ("bartlett", bartlett_f32(m).unwrap(), bartlett(m).unwrap()),
1197                ("blackman", blackman_f32(m).unwrap(), blackman(m).unwrap()),
1198                ("hamming", hamming_f32(m).unwrap(), hamming(m).unwrap()),
1199                ("hanning", hanning_f32(m).unwrap(), hanning(m).unwrap()),
1200            ] {
1201                assert_eq!(got.shape(), want.shape(), "{name} m={m} shape");
1202                let g = got.as_slice().unwrap();
1203                let w = want.as_slice().unwrap();
1204                for (i, (&a, &b)) in g.iter().zip(w.iter()).enumerate() {
1205                    let bf = b as f32;
1206                    assert!(
1207                        (a - bf).abs() < 1e-6,
1208                        "{name} m={m} idx={i}: f32 {a} vs f64 {b}"
1209                    );
1210                }
1211            }
1212        }
1213    }
1214
1215    #[test]
1216    fn kaiser_f32_matches_kaiser() {
1217        let m = 20;
1218        let beta = 8.6;
1219        let f32_arr = kaiser_f32(m, beta).unwrap();
1220        let f64_arr = kaiser(m, beta).unwrap();
1221        for (a, b) in f32_arr.iter().zip(f64_arr.iter()) {
1222            assert!((a - *b as f32).abs() < 1e-5);
1223        }
1224    }
1225
1226    // Helper: compare two slices within tolerance
1227    fn assert_close(actual: &[f64], expected: &[f64], tol: f64) {
1228        assert_eq!(
1229            actual.len(),
1230            expected.len(),
1231            "length mismatch: {} vs {}",
1232            actual.len(),
1233            expected.len()
1234        );
1235        for (i, (a, e)) in actual.iter().zip(expected.iter()).enumerate() {
1236            assert!(
1237                (a - e).abs() <= tol,
1238                "element {i}: {a} vs {e} (diff = {})",
1239                (a - e).abs()
1240            );
1241        }
1242    }
1243
1244    // -----------------------------------------------------------------------
1245    // Edge cases: M=0 and M=1
1246    // -----------------------------------------------------------------------
1247
1248    #[test]
1249    fn bartlett_m0() {
1250        let w = bartlett(0).unwrap();
1251        assert_eq!(w.shape(), &[0]);
1252    }
1253
1254    #[test]
1255    fn bartlett_m1() {
1256        let w = bartlett(1).unwrap();
1257        assert_eq!(w.as_slice().unwrap(), &[1.0]);
1258    }
1259
1260    #[test]
1261    fn blackman_m0() {
1262        let w = blackman(0).unwrap();
1263        assert_eq!(w.shape(), &[0]);
1264    }
1265
1266    #[test]
1267    fn blackman_m1() {
1268        let w = blackman(1).unwrap();
1269        assert_eq!(w.as_slice().unwrap(), &[1.0]);
1270    }
1271
1272    #[test]
1273    fn hamming_m0() {
1274        let w = hamming(0).unwrap();
1275        assert_eq!(w.shape(), &[0]);
1276    }
1277
1278    #[test]
1279    fn hamming_m1() {
1280        let w = hamming(1).unwrap();
1281        assert_eq!(w.as_slice().unwrap(), &[1.0]);
1282    }
1283
1284    #[test]
1285    fn hanning_m0() {
1286        let w = hanning(0).unwrap();
1287        assert_eq!(w.shape(), &[0]);
1288    }
1289
1290    #[test]
1291    fn hanning_m1() {
1292        let w = hanning(1).unwrap();
1293        assert_eq!(w.as_slice().unwrap(), &[1.0]);
1294    }
1295
1296    #[test]
1297    fn kaiser_m0() {
1298        let w = kaiser(0, 14.0).unwrap();
1299        assert_eq!(w.shape(), &[0]);
1300    }
1301
1302    #[test]
1303    fn kaiser_m1() {
1304        let w = kaiser(1, 14.0).unwrap();
1305        assert_eq!(w.as_slice().unwrap(), &[1.0]);
1306    }
1307
1308    #[test]
1309    fn kaiser_negative_beta() {
1310        // Negative beta should produce the same result as positive beta (I0 is even)
1311        let pos = kaiser(5, 1.0).unwrap();
1312        let neg = kaiser(5, -1.0).unwrap();
1313        assert_eq!(pos.as_slice().unwrap(), neg.as_slice().unwrap());
1314    }
1315
1316    #[test]
1317    fn kaiser_nan_beta() {
1318        assert!(kaiser(5, f64::NAN).is_err());
1319    }
1320
1321    // -----------------------------------------------------------------------
1322    // AC-1: bartlett(5) matches np.bartlett(5) to within 4 ULPs
1323    // -----------------------------------------------------------------------
1324    // np.bartlett(5) = [0.0, 0.5, 1.0, 0.5, 0.0]
1325    #[test]
1326    fn bartlett_5_ac1() {
1327        let w = bartlett(5).unwrap();
1328        let expected = [0.0, 0.5, 1.0, 0.5, 0.0];
1329        assert_close(w.as_slice().unwrap(), &expected, 1e-14);
1330    }
1331
1332    // -----------------------------------------------------------------------
1333    // AC-2: kaiser(5, 14.0) matches np.kaiser(5, 14.0) to within 4 ULPs
1334    // -----------------------------------------------------------------------
1335    // np.kaiser(5, 14.0) = [7.72686684e-06, 1.64932188e-01, 1.0, 1.64932188e-01, 7.72686684e-06]
1336    #[test]
1337    fn kaiser_5_14_ac2() {
1338        let w = kaiser(5, 14.0).unwrap();
1339        let s = w.as_slice().unwrap();
1340        assert_eq!(s.len(), 5);
1341        // NumPy reference values (verified via np.kaiser(5, 14.0))
1342        let expected: [f64; 5] = [
1343            7.72686684e-06,
1344            1.64932188e-01,
1345            1.0,
1346            1.64932188e-01,
1347            7.72686684e-06,
1348        ];
1349        // Bessel polynomial approximation has limited precision (~6 digits)
1350        for (i, (&a, &e)) in s.iter().zip(expected.iter()).enumerate() {
1351            let tol = if e.abs() < 1e-4 { 1e-8 } else { 1e-6 };
1352            assert!(
1353                (a - e).abs() <= tol,
1354                "kaiser element {i}: {a} vs {e} (diff = {})",
1355                (a - e).abs()
1356            );
1357        }
1358    }
1359
1360    // -----------------------------------------------------------------------
1361    // AC-3: All 5 window functions return correct length and match NumPy fixtures
1362    // -----------------------------------------------------------------------
1363
1364    // np.blackman(5) = [-1.38777878e-17, 3.40000000e-01, 1.00000000e+00, 3.40000000e-01, -1.38777878e-17]
1365    #[test]
1366    fn blackman_5() {
1367        let w = blackman(5).unwrap();
1368        assert_eq!(w.shape(), &[5]);
1369        let s = w.as_slice().unwrap();
1370        let expected = [
1371            -1.38777878e-17,
1372            3.40000000e-01,
1373            1.00000000e+00,
1374            3.40000000e-01,
1375            -1.38777878e-17,
1376        ];
1377        assert_close(s, &expected, 1e-10);
1378    }
1379
1380    // np.hamming(5) = [0.08, 0.54, 1.0, 0.54, 0.08]
1381    #[test]
1382    fn hamming_5() {
1383        let w = hamming(5).unwrap();
1384        assert_eq!(w.shape(), &[5]);
1385        let s = w.as_slice().unwrap();
1386        let expected = [0.08, 0.54, 1.0, 0.54, 0.08];
1387        assert_close(s, &expected, 1e-14);
1388    }
1389
1390    // np.hanning(5) = [0.0, 0.5, 1.0, 0.5, 0.0]
1391    #[test]
1392    fn hanning_5() {
1393        let w = hanning(5).unwrap();
1394        assert_eq!(w.shape(), &[5]);
1395        let s = w.as_slice().unwrap();
1396        let expected = [0.0, 0.5, 1.0, 0.5, 0.0];
1397        assert_close(s, &expected, 1e-14);
1398    }
1399
1400    // Larger window: np.bartlett(12)
1401    #[test]
1402    fn bartlett_12() {
1403        let w = bartlett(12).unwrap();
1404        assert_eq!(w.shape(), &[12]);
1405        let s = w.as_slice().unwrap();
1406        // First and last should be 0, peak near center
1407        assert!((s[0] - 0.0).abs() < 1e-14);
1408        assert!((s[11] - 0.0).abs() < 1e-14);
1409        // Symmetric
1410        for i in 0..6 {
1411            assert!((s[i] - s[11 - i]).abs() < 1e-14, "symmetry at {i}");
1412        }
1413    }
1414
1415    // Symmetry test for all windows
1416    #[test]
1417    fn all_windows_symmetric() {
1418        let m = 7;
1419        let windows: Vec<Array<f64, Ix1>> = vec![
1420            bartlett(m).unwrap(),
1421            blackman(m).unwrap(),
1422            hamming(m).unwrap(),
1423            hanning(m).unwrap(),
1424            kaiser(m, 5.0).unwrap(),
1425        ];
1426        for (idx, w) in windows.iter().enumerate() {
1427            let s = w.as_slice().unwrap();
1428            for i in 0..m / 2 {
1429                assert!(
1430                    (s[i] - s[m - 1 - i]).abs() < 1e-12,
1431                    "window {idx} not symmetric at {i}"
1432                );
1433            }
1434        }
1435    }
1436
1437    // All windows peak at center
1438    #[test]
1439    fn all_windows_peak_at_center() {
1440        let m = 9;
1441        let windows: Vec<Array<f64, Ix1>> = vec![
1442            bartlett(m).unwrap(),
1443            blackman(m).unwrap(),
1444            hamming(m).unwrap(),
1445            hanning(m).unwrap(),
1446            kaiser(m, 5.0).unwrap(),
1447        ];
1448        for (idx, w) in windows.iter().enumerate() {
1449            let s = w.as_slice().unwrap();
1450            let center = s[m / 2];
1451            assert!(
1452                (center - 1.0).abs() < 1e-10,
1453                "window {idx} center = {center}, expected 1.0"
1454            );
1455        }
1456    }
1457
1458    // Kaiser with beta=0 should be a rectangular window (all ones)
1459    #[test]
1460    fn kaiser_beta_zero() {
1461        let w = kaiser(5, 0.0).unwrap();
1462        let s = w.as_slice().unwrap();
1463        for &v in s {
1464            assert!((v - 1.0).abs() < 1e-10, "expected 1.0, got {v}");
1465        }
1466    }
1467
1468    #[test]
1469    fn bessel_i0_scalar_zero() {
1470        // Cephes Chebyshev expansion (issue #750): bit-near-exact at origin.
1471        assert!((bessel_i0_scalar::<f64>(0.0) - 1.0).abs() < 1e-13);
1472    }
1473
1474    #[test]
1475    fn bessel_i0_scalar_known() {
1476        // Issue #750: bessel_i0_scalar replaced A&S 9.8.1/9.8.2 with the
1477        // Cephes Chebyshev expansion (~1e-15 relative precision). Reference
1478        // values are from scipy.special.i0 (also Cephes).
1479        // I0(1) ~ 1.2660658777520082
1480        assert!((bessel_i0_scalar::<f64>(1.0) - 1.266_065_877_752_008_2).abs() < 1e-13);
1481        // I0(5) ~ 27.239871823604442 (now in the |x| <= 8 branch)
1482        assert!(
1483            (bessel_i0_scalar::<f64>(5.0) - 27.239_871_823_604_442).abs() / 27.239_871_823_604_442
1484                < 1e-13
1485        );
1486        // I0(10) ~ 2815.7166284662544 (asymptotic |x| > 8 branch)
1487        let expected_i0_10 = 2_815.716_628_466_254;
1488        assert!((bessel_i0_scalar::<f64>(10.0) - expected_i0_10).abs() / expected_i0_10 < 1e-13);
1489    }
1490
1491    // ----- Edge-length window coverage (#295) -----
1492
1493    #[test]
1494    fn bartlett_m2_is_zeros() {
1495        // Bartlett of length 2: w(n) = 1 - |2n/(M-1) - 1|
1496        // n=0: 1 - |0 - 1| = 0
1497        // n=1: 1 - |2 - 1| = 0
1498        let w = bartlett(2).unwrap();
1499        assert_eq!(w.as_slice().unwrap(), &[0.0, 0.0]);
1500    }
1501
1502    #[test]
1503    fn hanning_m2_is_zeros() {
1504        // Hann of length 2: w(n) = 0.5 * (1 - cos(2*pi*n/(M-1)))
1505        // n=0: 0.5 * (1 - cos(0)) = 0
1506        // n=1: 0.5 * (1 - cos(2*pi)) = 0
1507        let w = hanning(2).unwrap();
1508        let s = w.as_slice().unwrap();
1509        assert!(s[0].abs() < 1e-14);
1510        assert!(s[1].abs() < 1e-14);
1511    }
1512
1513    #[test]
1514    fn bartlett_even_length_is_symmetric() {
1515        // Even-length windows have no single "center"; symmetry holds
1516        // around the midpoint between indices (M/2 - 1) and (M/2).
1517        for &m in &[4usize, 6, 8, 10] {
1518            let w = bartlett(m).unwrap();
1519            let s = w.as_slice().unwrap();
1520            for i in 0..m / 2 {
1521                assert!(
1522                    (s[i] - s[m - 1 - i]).abs() < 1e-14,
1523                    "bartlett({m}) not symmetric at {i}"
1524                );
1525            }
1526        }
1527    }
1528
1529    #[test]
1530    fn blackman_even_length_is_symmetric() {
1531        for &m in &[4usize, 6, 8, 10] {
1532            let w = blackman(m).unwrap();
1533            let s = w.as_slice().unwrap();
1534            for i in 0..m / 2 {
1535                assert!(
1536                    (s[i] - s[m - 1 - i]).abs() < 1e-14,
1537                    "blackman({m}) not symmetric at {i}"
1538                );
1539            }
1540        }
1541    }
1542
1543    #[test]
1544    fn kaiser_large_beta_errors() {
1545        // Issue #294: beta above the f64 safe range should error,
1546        // not silently produce NaN.
1547        assert!(kaiser(8, 800.0).is_err());
1548        assert!(kaiser(8, 1000.0).is_err());
1549        // 700 is still safe.
1550        assert!(kaiser(8, 700.0).is_ok());
1551    }
1552
1553    // =======================================================================
1554    // SciPy-extended windows (cosine, exponential, gaussian, general_cosine,
1555    // general_hamming, nuttall, parzen, taylor, tukey).
1556    // =======================================================================
1557
1558    fn close(a: f64, b: f64, tol: f64) -> bool {
1559        (a - b).abs() < tol
1560    }
1561
1562    // ----- cosine -----
1563
1564    #[test]
1565    fn cosine_length_and_endpoints() {
1566        // cosine(M) = sin(pi(n + 0.5)/M); endpoints are sin(pi/2M) and
1567        // sin(pi(M-0.5)/M) — both small but nonzero.
1568        let w = cosine(8).unwrap();
1569        let s = w.as_slice().unwrap();
1570        assert_eq!(s.len(), 8);
1571        // Symmetric.
1572        for i in 0..4 {
1573            assert!(close(s[i], s[7 - i], 1e-14));
1574        }
1575    }
1576
1577    #[test]
1578    fn cosine_m1_and_m0() {
1579        assert_eq!(cosine(0).unwrap().shape(), &[0]);
1580        assert_eq!(cosine(1).unwrap().as_slice().unwrap(), &[1.0]);
1581    }
1582
1583    // ----- exponential -----
1584
1585    #[test]
1586    fn exponential_centred_default() {
1587        // tau=1, default centre = (M-1)/2 = 3.5 for M=8 — midpoint
1588        // between samples 3 and 4.
1589        let w = exponential(8, None, 1.0).unwrap();
1590        let s = w.as_slice().unwrap();
1591        // Symmetric.
1592        for i in 0..4 {
1593            assert!(close(s[i], s[7 - i], 1e-14));
1594        }
1595        // Centre samples have largest values.
1596        let centre_max = s[3].max(s[4]);
1597        for &v in s {
1598            assert!(v <= centre_max + 1e-14);
1599        }
1600    }
1601
1602    #[test]
1603    fn exponential_rejects_nonpositive_tau() {
1604        assert!(exponential(8, None, 0.0).is_err());
1605        assert!(exponential(8, None, -1.0).is_err());
1606        assert!(exponential(8, None, f64::NAN).is_err());
1607    }
1608
1609    // ----- gaussian -----
1610
1611    #[test]
1612    fn gaussian_centre_is_one() {
1613        // For odd M, centre sample is exactly 1.
1614        let w = gaussian(11, 2.0).unwrap();
1615        let s = w.as_slice().unwrap();
1616        assert!(close(s[5], 1.0, 1e-14));
1617    }
1618
1619    #[test]
1620    fn gaussian_known_value() {
1621        // gaussian(7, 1) at n=4: z = (4 - 3) / 1 = 1, exp(-0.5) = 0.6065...
1622        let w = gaussian(7, 1.0).unwrap();
1623        let s = w.as_slice().unwrap();
1624        assert!(close(s[4], (-0.5_f64).exp(), 1e-14));
1625    }
1626
1627    #[test]
1628    fn gaussian_rejects_nonpositive_std() {
1629        assert!(gaussian(8, 0.0).is_err());
1630        assert!(gaussian(8, -1.0).is_err());
1631    }
1632
1633    // ----- general_cosine -----
1634
1635    #[test]
1636    fn general_cosine_with_hann_coeffs_matches_hann() {
1637        // [0.5, 0.5] → 0.5 - 0.5 cos(2πn/(M-1)) = Hann.
1638        let m = 9;
1639        let gc = general_cosine(m, &[0.5, 0.5]).unwrap();
1640        let hn = hanning(m).unwrap();
1641        for i in 0..m {
1642            assert!(close(
1643                gc.as_slice().unwrap()[i],
1644                hn.as_slice().unwrap()[i],
1645                1e-14,
1646            ));
1647        }
1648    }
1649
1650    #[test]
1651    fn general_cosine_with_blackman_coeffs_matches_blackman() {
1652        // [0.42, 0.5, 0.08] → classical Blackman.
1653        let m = 9;
1654        let gc = general_cosine(m, &[0.42, 0.5, 0.08]).unwrap();
1655        let bk = blackman(m).unwrap();
1656        for i in 0..m {
1657            assert!(close(
1658                gc.as_slice().unwrap()[i],
1659                bk.as_slice().unwrap()[i],
1660                1e-12,
1661            ));
1662        }
1663    }
1664
1665    #[test]
1666    fn general_cosine_rejects_empty_coeffs() {
1667        assert!(general_cosine(8, &[]).is_err());
1668    }
1669
1670    // ----- general_hamming -----
1671
1672    #[test]
1673    fn general_hamming_alpha_half_matches_hann() {
1674        let m = 9;
1675        let gh = general_hamming(m, 0.5).unwrap();
1676        let hn = hanning(m).unwrap();
1677        for i in 0..m {
1678            assert!(close(
1679                gh.as_slice().unwrap()[i],
1680                hn.as_slice().unwrap()[i],
1681                1e-14,
1682            ));
1683        }
1684    }
1685
1686    #[test]
1687    fn general_hamming_alpha_054_matches_hamming() {
1688        let m = 9;
1689        let gh = general_hamming(m, 0.54).unwrap();
1690        let hm = hamming(m).unwrap();
1691        for i in 0..m {
1692            assert!(close(
1693                gh.as_slice().unwrap()[i],
1694                hm.as_slice().unwrap()[i],
1695                1e-14,
1696            ));
1697        }
1698    }
1699
1700    // ----- nuttall -----
1701
1702    #[test]
1703    fn nuttall_endpoints_are_small() {
1704        // Nuttall is engineered to have very low sidelobes; endpoints
1705        // are O(1e-3) for typical M.
1706        let w = nuttall(64).unwrap();
1707        let s = w.as_slice().unwrap();
1708        assert!(s[0].abs() < 1e-2);
1709        assert!(s[s.len() - 1].abs() < 1e-2);
1710    }
1711
1712    #[test]
1713    fn nuttall_is_symmetric() {
1714        let m = 33;
1715        let w = nuttall(m).unwrap();
1716        let s = w.as_slice().unwrap();
1717        for i in 0..m / 2 {
1718            assert!(close(s[i], s[m - 1 - i], 1e-14));
1719        }
1720    }
1721
1722    // ----- parzen -----
1723
1724    #[test]
1725    fn parzen_endpoints_are_zero() {
1726        let w = parzen(16).unwrap();
1727        let s = w.as_slice().unwrap();
1728        // Outer edge of Parzen has very small value (cubic falloff).
1729        assert!(s[0].abs() < 1e-2);
1730        assert!(s[s.len() - 1].abs() < 1e-2);
1731    }
1732
1733    #[test]
1734    fn parzen_centre_is_one() {
1735        // For odd M, centre is at sample (M-1)/2; r = 0 → w = 1.
1736        let w = parzen(13).unwrap();
1737        let s = w.as_slice().unwrap();
1738        assert!(close(s[6], 1.0, 1e-14));
1739    }
1740
1741    #[test]
1742    fn parzen_is_symmetric() {
1743        let m = 21;
1744        let w = parzen(m).unwrap();
1745        let s = w.as_slice().unwrap();
1746        for i in 0..m / 2 {
1747            assert!(close(s[i], s[m - 1 - i], 1e-14));
1748        }
1749    }
1750
1751    // ----- taylor -----
1752
1753    #[test]
1754    fn taylor_default_normalised_centre_is_one() {
1755        // With norm=true the centre sample is 1.0.
1756        let w = taylor(33, 4, 30.0, true).unwrap();
1757        let s = w.as_slice().unwrap();
1758        assert!(close(s[16], 1.0, 1e-12));
1759    }
1760
1761    #[test]
1762    fn taylor_is_symmetric() {
1763        let m = 33;
1764        let w = taylor(m, 4, 30.0, true).unwrap();
1765        let s = w.as_slice().unwrap();
1766        for i in 0..m / 2 {
1767            assert!(close(s[i], s[m - 1 - i], 1e-12));
1768        }
1769    }
1770
1771    #[test]
1772    fn taylor_rejects_nbar_zero() {
1773        assert!(taylor(8, 0, 30.0, true).is_err());
1774        assert!(taylor(8, 4, f64::NAN, true).is_err());
1775    }
1776
1777    #[test]
1778    fn taylor_even_m_matches_scipy_analytic_peak() {
1779        // Regression for #810 (upstream from ferrotorch). Pre-fix, even-M
1780        // taylor windows diverged from scipy by ~1.5e-3 because the centre
1781        // value was estimated by averaging the two adjacent samples
1782        // (the secant midpoint) rather than evaluating the analytic peak
1783        // W((M-1)/2) = 1 + 2·Σ F_m at the fractional midpoint.
1784        //
1785        // scipy.signal.windows.taylor(16, nbar=4, sll=30, norm=True, sym=True)[0]
1786        //   = 0.252321041674507
1787        let w = taylor(16, 4, 30.0, true).unwrap();
1788        let s = w.as_slice().unwrap();
1789        let expected_first = 0.252_321_041_674_507_f64;
1790        assert!(
1791            (s[0] - expected_first).abs() < 1e-13,
1792            "taylor(16,4,30,true)[0] = {} expected {}",
1793            s[0],
1794            expected_first,
1795        );
1796        // The window is symmetric and even-length, so no sample sits
1797        // exactly at the centre — but the analytic peak (the value the
1798        // window would take at xn = 0) is, by construction, 1.0. The
1799        // closest pair of samples should lie just below 1.
1800        assert!(s[7] < 1.0 && s[8] < 1.0);
1801        assert!(s[7] > 0.99 && s[8] > 0.99);
1802    }
1803
1804    // ----- tukey -----
1805
1806    #[test]
1807    fn tukey_alpha_zero_is_rectangular() {
1808        let w = tukey(8, 0.0).unwrap();
1809        for &v in w.as_slice().unwrap() {
1810            assert!(close(v, 1.0, 1e-14));
1811        }
1812    }
1813
1814    #[test]
1815    fn tukey_alpha_one_matches_hann() {
1816        // alpha=1 → fully Hann-shaped.
1817        let m = 9;
1818        let tk = tukey(m, 1.0).unwrap();
1819        let hn = hanning(m).unwrap();
1820        for i in 0..m {
1821            assert!(close(
1822                tk.as_slice().unwrap()[i],
1823                hn.as_slice().unwrap()[i],
1824                1e-12,
1825            ));
1826        }
1827    }
1828
1829    #[test]
1830    fn tukey_centre_is_one() {
1831        let m = 21;
1832        let w = tukey(m, 0.5).unwrap();
1833        // Wide flat-top region: middle must be 1.
1834        assert!(close(w.as_slice().unwrap()[m / 2], 1.0, 1e-14));
1835    }
1836
1837    #[test]
1838    fn tukey_rejects_invalid_alpha() {
1839        assert!(tukey(8, -0.1).is_err());
1840        assert!(tukey(8, 1.1).is_err());
1841        assert!(tukey(8, f64::NAN).is_err());
1842    }
1843}
1844
1845#[cfg(test)]
1846mod debug_test_removed {
1847    // Debug test scaffolding for #740 chebwin development;
1848    // removed once the closed-form output matched scipy.
1849}