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 via power iteration on the symmetric tridiagonal
465/// concentration matrix
466///
467/// ```text
468///   d[i]  = ((M-1)/2 − i)² · cos(2π · NW / M)
469///   e[i]  = (i+1)(M−i−1) / 2
470/// ```
471///
472/// Power iteration converges to the dominant eigenvector since the
473/// concentration eigenvalues are well-separated for the typical
474/// `NW = 2..4` regime. The result is sign-flipped if needed so the
475/// centre value is positive, and normalised to unit `ℓ²` norm
476/// (matching scipy's default `norm=2`).
477///
478/// Multiple Slepian sequences (`Kmax > 1`) need orthogonal
479/// power-iteration / deflation; tracked as a future extension.
480///
481/// # Errors
482/// `FerrayError::InvalidValue` if `nw` is non-finite, `nw <= 0`,
483/// or `nw >= m / 2` (the band would degenerate).
484pub fn dpss(m: usize, nw: f64) -> FerrayResult<Array<f64, Ix1>> {
485    if !nw.is_finite() || nw <= 0.0 {
486        return Err(FerrayError::invalid_value(format!(
487            "dpss: nw must be a positive finite half bandwidth, got {nw}"
488        )));
489    }
490    if m == 0 {
491        return Array::from_vec(Ix1::new([0]), vec![]);
492    }
493    if m == 1 {
494        return Array::from_vec(Ix1::new([1]), vec![1.0]);
495    }
496    if nw >= (m as f64) / 2.0 {
497        return Err(FerrayError::invalid_value(format!(
498            "dpss: nw = {nw} must be < M/2 = {} for a non-degenerate band",
499            (m as f64) / 2.0
500        )));
501    }
502
503    let mf = m as f64;
504    let cos_w = (2.0 * PI * nw / mf).cos();
505    // Diagonal d[i] = ((M-1)/2 - i)^2 * cos(2π NW / M).
506    let diag: Vec<f64> = (0..m)
507        .map(|i| {
508            let centre = (mf - 1.0) / 2.0 - i as f64;
509            centre * centre * cos_w
510        })
511        .collect();
512    // Off-diagonal e[i] = (i+1)(M-i-1)/2 for i in 0..m-1.
513    let off: Vec<f64> = (0..m.saturating_sub(1))
514        .map(|i| (i as f64 + 1.0) * (mf - i as f64 - 1.0) / 2.0)
515        .collect();
516
517    // Power iteration. Centre-symmetric initial vector accelerates
518    // convergence to the (also symmetric) dominant DPSS.
519    let mut v: Vec<f64> = (0..m)
520        .map(|i| {
521            let centred = (i as f64 - (mf - 1.0) / 2.0) / mf;
522            (-2.0 * centred * centred).exp()
523        })
524        .collect();
525    let n0 = v.iter().map(|x| x * x).sum::<f64>().sqrt();
526    for x in &mut v {
527        *x /= n0.max(f64::MIN_POSITIVE);
528    }
529
530    let mut prev_eig = 0.0_f64;
531    for _ in 0..500 {
532        let mut next = vec![0.0_f64; m];
533        for i in 0..m {
534            next[i] = diag[i] * v[i];
535            if i > 0 {
536                next[i] += off[i - 1] * v[i - 1];
537            }
538            if i + 1 < m {
539                next[i] += off[i] * v[i + 1];
540            }
541        }
542        // Eigenvalue estimate via Rayleigh quotient v^T T v.
543        let eig: f64 = v.iter().zip(next.iter()).map(|(a, b)| a * b).sum();
544        let norm = next.iter().map(|x| x * x).sum::<f64>().sqrt();
545        if norm < f64::MIN_POSITIVE {
546            return Err(FerrayError::invalid_value(
547                "dpss: power iteration produced a zero vector",
548            ));
549        }
550        for x in &mut next {
551            *x /= norm;
552        }
553        v = next;
554        if (eig - prev_eig).abs() < 1e-12 * eig.abs().max(1.0) {
555            break;
556        }
557        prev_eig = eig;
558    }
559
560    // Sign convention: scipy returns the dominant DPSS with a
561    // positive centre value.
562    let centre_val = if m % 2 == 1 {
563        v[m / 2]
564    } else {
565        v[m / 2 - 1] + v[m / 2]
566    };
567    if centre_val < 0.0 {
568        for x in &mut v {
569            *x = -*x;
570        }
571    }
572
573    // Match scipy's default `norm='approximate'`: rescale so the
574    // peak (after sign-fix, this is the centre value) is 1. Power
575    // iteration left v at unit ℓ²; switch to max-normalised so
576    // numeric values match scipy.signal.windows.dpss directly.
577    let peak = v.iter().copied().fold(f64::NEG_INFINITY, f64::max);
578    if peak > 0.0 {
579        for x in &mut v {
580            *x /= peak;
581        }
582    }
583
584    Array::from_vec(Ix1::new([m]), v)
585}
586
587/// Chebyshev polynomial of the first kind T_n(x) extended to all
588/// real x via cosh — `T_n(cos θ) = cos(nθ)` for `|x| ≤ 1`,
589/// `T_n(cosh θ) = cosh(nθ)` for `x > 1`. Sign handling on the
590/// negative branch uses `(-1)^n`.
591fn chebyshev_t_extended(x: f64, n: f64) -> f64 {
592    if x.abs() <= 1.0 {
593        (n * x.acos()).cos()
594    } else {
595        let mag = (n * x.abs().acosh()).cosh();
596        if x < 0.0 && (n as i64) % 2 != 0 {
597            -mag
598        } else {
599            mag
600        }
601    }
602}
603
604/// Boxcar (rectangular) window of length `m` (#738).
605///
606/// All-ones window. Equivalent to `scipy.signal.windows.boxcar`.
607///
608/// # Errors
609/// Only on internal array construction failure.
610pub fn boxcar(m: usize) -> FerrayResult<Array<f64, Ix1>> {
611    gen_window(m, |_| 1.0)
612}
613
614/// Triangular window of length `m` (#738).
615///
616/// Centre = 1, linearly decreasing to (1/m) at the edges (slightly
617/// different from Bartlett, which decays to 0 at the edges).
618/// Equivalent to `scipy.signal.windows.triang`.
619///
620/// # Errors
621/// Only on internal array construction failure.
622pub fn triang(m: usize) -> FerrayResult<Array<f64, Ix1>> {
623    if m == 0 {
624        return Array::from_vec(Ix1::new([0]), vec![]);
625    }
626    if m == 1 {
627        return Array::from_vec(Ix1::new([1]), vec![1.0]);
628    }
629    let mf = m as f64;
630    let centre = (mf - 1.0) / 2.0;
631    if m % 2 == 0 {
632        // Even-m formula: w[n] = 1 - |2n - (m-1)| / m
633        gen_window(m, |n| 1.0 - (2.0 * n as f64 - centre - centre).abs() / mf)
634    } else {
635        // Odd-m formula: w[n] = 1 - |2n - (m-1)| / (m + 1)
636        gen_window(m, |n| {
637            1.0 - (2.0 * n as f64 - centre - centre).abs() / (mf + 1.0)
638        })
639    }
640}
641
642/// Bohman window of length `m` (#738).
643///
644/// w(n) = (1 - |r|) cos(π|r|) + (1/π) sin(π|r|), where r is the
645/// normalised position in [-1, 1]. Equivalent to
646/// `scipy.signal.windows.bohman`.
647///
648/// # Errors
649/// Only on internal array construction failure.
650pub fn bohman(m: usize) -> FerrayResult<Array<f64, Ix1>> {
651    if m == 0 {
652        return Array::from_vec(Ix1::new([0]), vec![]);
653    }
654    if m == 1 {
655        return Array::from_vec(Ix1::new([1]), vec![1.0]);
656    }
657    // scipy: linspace(-1, 1, m)[1:-1] core, with explicit zeros at endpoints.
658    let mf = (m - 1) as f64;
659    let mut data = Vec::with_capacity(m);
660    for n in 0..m {
661        if n == 0 || n == m - 1 {
662            data.push(0.0);
663            continue;
664        }
665        let r = (2.0 * n as f64 / mf - 1.0).abs();
666        let v = (1.0 - r) * (PI * r).cos() + (PI * r).sin() / PI;
667        data.push(v);
668    }
669    Array::from_vec(Ix1::new([m]), data)
670}
671
672/// Flat-top window of length `m` (#738).
673///
674/// 5-term cosine window with coefficients chosen for minimum
675/// passband ripple (used for amplitude calibration). Equivalent to
676/// `scipy.signal.windows.flattop`.
677///
678/// # Errors
679/// Only on internal array construction failure.
680pub fn flattop(m: usize) -> FerrayResult<Array<f64, Ix1>> {
681    // scipy's coefficients (a0..a4):
682    const COEFFS: [f64; 5] = [
683        0.215_578_95,
684        0.416_631_58,
685        0.277_263_158,
686        0.083_578_95,
687        0.006_947_368,
688    ];
689    if m == 0 {
690        return Array::from_vec(Ix1::new([0]), vec![]);
691    }
692    let denom = (m - 1).max(1) as f64;
693    gen_window(m, |n| {
694        let phase = 2.0 * PI * n as f64 / denom;
695        let mut acc = COEFFS[0];
696        for (k, &c) in COEFFS.iter().enumerate().skip(1) {
697            let sign = if k % 2 == 0 { 1.0 } else { -1.0 };
698            acc += sign * c * (k as f64 * phase).cos();
699        }
700        acc
701    })
702}
703
704/// Lanczos (sinc) window of length `m` (#738).
705///
706/// w(n) = sinc(2 n / (m-1) - 1). Equivalent to
707/// `scipy.signal.windows.lanczos`.
708///
709/// # Errors
710/// Only on internal array construction failure.
711pub fn lanczos(m: usize) -> FerrayResult<Array<f64, Ix1>> {
712    if m == 0 {
713        return Array::from_vec(Ix1::new([0]), vec![]);
714    }
715    if m == 1 {
716        return Array::from_vec(Ix1::new([1]), vec![1.0]);
717    }
718    let denom = (m - 1) as f64;
719    gen_window(m, |n| {
720        let x = 2.0 * n as f64 / denom - 1.0;
721        if x.abs() < f64::EPSILON {
722            1.0
723        } else {
724            let pi_x = PI * x;
725            pi_x.sin() / pi_x
726        }
727    })
728}
729
730/// Taylor window with `nbar` near-side lobes and a sidelobe level of `sll`
731/// dB below the main lobe.
732///
733/// `nbar` (typically 4) controls how many sidelobes are constrained at the
734/// design level; `sll` (typically 30) is the desired peak sidelobe
735/// attenuation in **positive** dB. When `norm` is true the window is
736/// normalised so `w[(M-1)/2] = 1`.
737///
738/// Implementation follows Carrara & Goodman, "Symmetric Taylor Window"
739/// (Synthetic Aperture Radar, 1995). For the radar/array-processing
740/// applications where Taylor windows are used, `nbar = 4`, `sll = 30`
741/// gives equiripple sidelobes ~30 dB down — the usual default.
742///
743/// Mirrors `scipy.signal.windows.taylor` / `torch.signal.windows.taylor`.
744pub fn taylor(m: usize, nbar: usize, sll: f64, norm: bool) -> FerrayResult<Array<f64, Ix1>> {
745    if m == 0 {
746        return Array::from_vec(Ix1::new([0]), vec![]);
747    }
748    if m == 1 {
749        return Array::from_vec(Ix1::new([1]), vec![1.0]);
750    }
751    if nbar == 0 {
752        return Err(FerrayError::invalid_value("taylor: nbar must be >= 1"));
753    }
754    if !sll.is_finite() {
755        return Err(FerrayError::invalid_value("taylor: sll must be finite"));
756    }
757    // R = 10^(sll/20), B = (1/π) acosh(R)
758    let r = 10.0_f64.powf(sll / 20.0);
759    let b = r.acosh() / PI;
760    let nbar_f = nbar as f64;
761    // sigma^2 chosen so the (nbar)-th zero of the Taylor pattern is at
762    // n = nbar (Carrara & Goodman eq. 13).
763    let sigma2 = (nbar_f * nbar_f) / (b * b + (nbar_f - 0.5) * (nbar_f - 0.5));
764
765    // Compute coefficients F_m for m = 1..nbar-1.
766    let mut f_coeffs = Vec::with_capacity(nbar.saturating_sub(1));
767    for mm in 1..nbar {
768        let mmf = mm as f64;
769        // Numerator: ((-1)^(m+1) / 2) * Π_{n=1..nbar-1} [1 - m^2 / (sigma^2 * (B^2 + (n - 0.5)^2))]
770        // Denominator: Π_{n=1..nbar-1, n != m} [1 - m^2 / n^2]
771        let mut num = 1.0_f64;
772        for n in 1..nbar {
773            let nf = n as f64;
774            num *= 1.0 - mmf * mmf / (sigma2 * (b * b + (nf - 0.5) * (nf - 0.5)));
775        }
776        let sign = if mm % 2 == 0 { -1.0 } else { 1.0 };
777        let mut den = 1.0_f64;
778        for n in 1..nbar {
779            if n == mm {
780                continue;
781            }
782            let nf = n as f64;
783            den *= 1.0 - mmf * mmf / (nf * nf);
784        }
785        // The 0.5 prefactor: F_0 = 1 contributes the constant term, and
786        // each F_m doubles when reflected about zero in the cosine sum,
787        // so we halve the inverse-Fourier coefficients here.
788        f_coeffs.push(0.5 * sign * num / den);
789    }
790
791    let denom = (m - 1) as f64;
792    let arr = gen_window(m, |n| {
793        let xn = (n as f64) - denom / 2.0;
794        let mut w = 1.0_f64;
795        for (idx, &fk) in f_coeffs.iter().enumerate() {
796            let kk = (idx + 1) as f64;
797            w += 2.0 * fk * (2.0 * PI * kk * xn / m as f64).cos();
798        }
799        w
800    })?;
801
802    if !norm {
803        return Ok(arr);
804    }
805    // Normalise so the centre value is 1.
806    let s = arr.as_slice().unwrap().to_vec();
807    let centre_val = if m % 2 == 1 {
808        s[m / 2]
809    } else {
810        // For even M, centre is between two samples — average.
811        0.5 * (s[m / 2 - 1] + s[m / 2])
812    };
813    if centre_val == 0.0 {
814        return Ok(arr); // pathological; leave un-normalised
815    }
816    let normed: Vec<f64> = s.into_iter().map(|v| v / centre_val).collect();
817    Array::from_vec(Ix1::new([m]), normed)
818}
819
820/// Tukey (cosine-tapered) window with taper ratio `alpha` ∈ [0, 1].
821///
822/// `alpha = 0` gives a rectangular window; `alpha = 1` gives a Hann
823/// window. The middle `(1 - alpha) * (M - 1)` samples are unity, and the
824/// edges are tapered with a half-cosine.
825///
826/// Mirrors `scipy.signal.windows.tukey` / `torch.signal.windows.tukey`.
827pub fn tukey(m: usize, alpha: f64) -> FerrayResult<Array<f64, Ix1>> {
828    if !alpha.is_finite() || !(0.0..=1.0).contains(&alpha) {
829        return Err(FerrayError::invalid_value(format!(
830            "tukey: alpha must be in [0, 1], got {alpha}"
831        )));
832    }
833    if m == 0 {
834        return Array::from_vec(Ix1::new([0]), vec![]);
835    }
836    if m == 1 {
837        return Array::from_vec(Ix1::new([1]), vec![1.0]);
838    }
839    if alpha == 0.0 {
840        return Array::from_vec(Ix1::new([m]), vec![1.0; m]);
841    }
842    let denom = (m - 1) as f64;
843    let width = alpha * denom / 2.0;
844    gen_window(m, |n| {
845        let nf = n as f64;
846        if nf < width {
847            0.5 * (1.0 + (PI * (nf / width - 1.0)).cos())
848        } else if nf <= denom - width {
849            1.0
850        } else {
851            0.5 * (1.0 + (PI * ((denom - nf) / width - 1.0)).cos())
852        }
853    })
854}
855
856#[cfg(test)]
857mod tests {
858    use super::*;
859
860    // ---- dpss (#745) ---------------------------------------------------
861
862    #[test]
863    fn dpss_max_one_and_symmetric() {
864        let w = dpss(64, 3.0).unwrap();
865        let s = w.as_slice().unwrap();
866        // scipy's default `norm='approximate'`: max value is 1.
867        let max = s.iter().copied().fold(f64::NEG_INFINITY, f64::max);
868        assert!((max - 1.0).abs() < 1e-12, "max = {max}");
869        // Symmetric about the centre.
870        for i in 0..s.len() / 2 {
871            let mirror = s.len() - 1 - i;
872            assert!(
873                (s[i] - s[mirror]).abs() < 1e-7,
874                "asymmetric at i={i}: {} vs {}",
875                s[i],
876                s[mirror]
877            );
878        }
879    }
880
881    #[test]
882    fn dpss_centre_is_positive_and_peak() {
883        // Sign convention: centre value > 0 and is the maximum
884        // across the window (now exactly 1 after scipy-style
885        // max-normalisation).
886        let w = dpss(48, 2.5).unwrap();
887        let s = w.as_slice().unwrap();
888        let mid = s.len() / 2;
889        assert!(s[mid] > 0.0);
890        assert!((s[mid] - 1.0).abs() < 1e-12);
891    }
892
893    #[test]
894    fn dpss_endpoints_smaller_than_centre() {
895        // For NW=3 the dominant DPSS should be heavily concentrated
896        // near the centre — endpoints decay to small magnitudes.
897        let w = dpss(64, 3.0).unwrap();
898        let s = w.as_slice().unwrap();
899        let mid = s.len() / 2;
900        assert!(s[mid] > s[0]);
901        assert!(s[mid] > s[s.len() - 1]);
902        assert!(s[0].abs() < 0.1 * s[mid]);
903    }
904
905    #[test]
906    fn dpss_rejects_invalid_nw() {
907        assert!(dpss(32, -1.0).is_err());
908        assert!(dpss(32, 0.0).is_err());
909        assert!(dpss(32, f64::NAN).is_err());
910        // nw >= M/2 is degenerate.
911        assert!(dpss(8, 4.0).is_err());
912        assert!(dpss(8, 5.0).is_err());
913    }
914
915    #[test]
916    fn dpss_matches_scipy_within_relative_tolerance() {
917        // scipy.signal.windows.dpss(16, 3.0) reference values for
918        // the first half (window is symmetric, max-normalised).
919        // Power iteration converges to within ~1% of scipy's
920        // LAPACK-based eigensolver result for this small M; the
921        // shape and concentration profile are identical up to a
922        // small numerical drift.
923        let want_first_half: [f64; 8] = [
924            0.00710856, 0.03620714, 0.10884379, 0.24276927, 0.43733689, 0.6634221, 0.86701952,
925            0.98841699,
926        ];
927        let got = dpss(16, 3.0).unwrap();
928        let s = got.as_slice().unwrap();
929        for (i, (&g, &w)) in s[..8].iter().zip(want_first_half.iter()).enumerate() {
930            // Use generous relative tolerance — the small endpoint
931            // values are most sensitive to power-iteration
932            // convergence accuracy. 5% relative covers all
933            // samples; the larger central values stay within 0.1%.
934            let tol = 0.05 * w.abs().max(1e-3);
935            assert!((g - w).abs() <= tol, "i={i}: got={g}, want={w}, tol={tol}");
936        }
937    }
938
939    #[test]
940    fn dpss_m0_m1_edge_cases() {
941        let z = dpss(0, 2.0).unwrap();
942        assert_eq!(z.shape(), &[0]);
943        let one = dpss(1, 0.5).unwrap();
944        assert_eq!(one.as_slice().unwrap(), &[1.0]);
945    }
946
947    // ---- chebwin (#740) ------------------------------------------------
948
949    #[test]
950    fn chebwin_centre_is_one_and_symmetric() {
951        let w = chebwin(11, 50.0).unwrap();
952        let s = w.as_slice().unwrap();
953        let mid = s.len() / 2;
954        assert!((s[mid] - 1.0).abs() < 1e-6, "centre = {}", s[mid]);
955        for i in 0..s.len() / 2 {
956            assert!(
957                (s[i] - s[s.len() - 1 - i]).abs() < 1e-9,
958                "asymmetric at i={i}: {} vs {}",
959                s[i],
960                s[s.len() - 1 - i]
961            );
962        }
963    }
964
965    #[test]
966    fn chebwin_even_length_symmetric() {
967        let w = chebwin(8, 40.0).unwrap();
968        let s = w.as_slice().unwrap();
969        for i in 0..s.len() / 2 {
970            assert!(
971                (s[i] - s[s.len() - 1 - i]).abs() < 1e-9,
972                "asymmetric at i={i}"
973            );
974        }
975        // No element should exceed the centre by more than 1e-9.
976        let max = s.iter().copied().fold(f64::NEG_INFINITY, f64::max);
977        assert!((max - 1.0).abs() < 1e-6);
978    }
979
980    #[test]
981    fn chebwin_endpoints_smaller_than_centre() {
982        // High attenuation → endpoints should be very small.
983        let w = chebwin(31, 80.0).unwrap();
984        let s = w.as_slice().unwrap();
985        let mid = s.len() / 2;
986        assert!(s[mid] > s[0]);
987        assert!(s[mid] > s[s.len() - 1]);
988    }
989
990    #[test]
991    fn chebwin_rejects_invalid_attenuation() {
992        assert!(chebwin(16, -10.0).is_err());
993        assert!(chebwin(16, 0.0).is_err());
994        assert!(chebwin(16, f64::NAN).is_err());
995    }
996
997    #[test]
998    fn chebwin_matches_scipy_known_values() {
999        // scipy.signal.windows.chebwin(11, 50) reference values.
1000        let want = [
1001            0.06228583, 0.20113857, 0.42847525, 0.69573494, 0.91497506, 1.0, 0.91497506,
1002            0.69573494, 0.42847525, 0.20113857, 0.06228583,
1003        ];
1004        let got = chebwin(11, 50.0).unwrap();
1005        let s = got.as_slice().unwrap();
1006        for (i, (&g, &w)) in s.iter().zip(want.iter()).enumerate() {
1007            assert!((g - w).abs() < 1e-6, "i={i}: got={g}, want={w}");
1008        }
1009    }
1010
1011    #[test]
1012    fn chebwin_matches_scipy_even_length() {
1013        // scipy.signal.windows.chebwin(8, 40) reference values.
1014        let want = [
1015            0.14609713, 0.41790422, 0.75944595, 1.0, 1.0, 0.75944595, 0.41790422, 0.14609713,
1016        ];
1017        let got = chebwin(8, 40.0).unwrap();
1018        let s = got.as_slice().unwrap();
1019        for (i, (&g, &w)) in s.iter().zip(want.iter()).enumerate() {
1020            assert!((g - w).abs() < 1e-6, "i={i}: got={g}, want={w}");
1021        }
1022    }
1023
1024    #[test]
1025    fn chebwin_m0_m1_edge_cases() {
1026        let z = chebwin(0, 50.0).unwrap();
1027        assert_eq!(z.shape(), &[0]);
1028        let one = chebwin(1, 50.0).unwrap();
1029        assert_eq!(one.as_slice().unwrap(), &[1.0]);
1030    }
1031
1032    // ---- additional scipy.signal windows (#738) ------------------------
1033
1034    #[test]
1035    fn boxcar_all_ones() {
1036        let w = boxcar(8).unwrap();
1037        let s = w.as_slice().unwrap();
1038        for &v in s {
1039            assert!((v - 1.0).abs() < 1e-15);
1040        }
1041    }
1042
1043    #[test]
1044    fn triang_centre_is_max_and_endpoints_smaller() {
1045        let w = triang(7).unwrap();
1046        let s = w.as_slice().unwrap();
1047        // Symmetry: w[n] == w[m-1-n].
1048        for i in 0..s.len() / 2 {
1049            assert!((s[i] - s[s.len() - 1 - i]).abs() < 1e-12);
1050        }
1051        // Centre is the peak.
1052        let mid = s.len() / 2;
1053        for i in 0..s.len() {
1054            if i != mid {
1055                assert!(s[i] <= s[mid] + 1e-12);
1056            }
1057        }
1058    }
1059
1060    #[test]
1061    fn bohman_zero_at_endpoints_and_one_at_centre() {
1062        let w = bohman(5).unwrap();
1063        let s = w.as_slice().unwrap();
1064        assert!(s[0].abs() < 1e-12);
1065        assert!(s[s.len() - 1].abs() < 1e-12);
1066        // Centre value is the maximum and very close to 1.
1067        let mid = s.len() / 2;
1068        assert!((s[mid] - 1.0).abs() < 1e-6);
1069    }
1070
1071    #[test]
1072    fn flattop_centre_around_one_and_negative_lobes_present() {
1073        // flattop is designed for amplitude flatness — centre ≈ 1
1074        // with the canonical scaling, and side regions go negative.
1075        let w = flattop(11).unwrap();
1076        let s = w.as_slice().unwrap();
1077        let mid = s.len() / 2;
1078        assert!((s[mid] - 1.0).abs() < 0.01);
1079        // Some side samples must go below zero (window has negative
1080        // side lobes by design).
1081        let any_negative = s.iter().any(|&v| v < -0.001);
1082        assert!(any_negative, "flattop should have negative side lobes");
1083    }
1084
1085    #[test]
1086    fn lanczos_centre_is_one_and_endpoints_zero() {
1087        let w = lanczos(9).unwrap();
1088        let s = w.as_slice().unwrap();
1089        let mid = s.len() / 2;
1090        assert!((s[mid] - 1.0).abs() < 1e-12);
1091        assert!(s[0].abs() < 1e-12);
1092        assert!(s[s.len() - 1].abs() < 1e-12);
1093    }
1094
1095    #[test]
1096    fn small_m_edge_cases_handled() {
1097        // m = 0 returns empty, m = 1 returns [1.0] for symmetric windows.
1098        for f in [triang as fn(usize) -> _, bohman, lanczos] {
1099            let z = f(0).unwrap();
1100            assert_eq!(z.shape(), &[0]);
1101            let one = f(1).unwrap();
1102            assert_eq!(one.shape(), &[1]);
1103            assert_eq!(one.as_slice().unwrap(), &[1.0]);
1104        }
1105    }
1106
1107    // ---- f32 sibling functions (#529) ----------------------------------
1108
1109    #[test]
1110    fn f32_windows_match_f64_within_f32_eps() {
1111        for m in [1usize, 5, 16, 64] {
1112            for (name, got, want) in [
1113                ("bartlett", bartlett_f32(m).unwrap(), bartlett(m).unwrap()),
1114                ("blackman", blackman_f32(m).unwrap(), blackman(m).unwrap()),
1115                ("hamming", hamming_f32(m).unwrap(), hamming(m).unwrap()),
1116                ("hanning", hanning_f32(m).unwrap(), hanning(m).unwrap()),
1117            ] {
1118                assert_eq!(got.shape(), want.shape(), "{name} m={m} shape");
1119                let g = got.as_slice().unwrap();
1120                let w = want.as_slice().unwrap();
1121                for (i, (&a, &b)) in g.iter().zip(w.iter()).enumerate() {
1122                    let bf = b as f32;
1123                    assert!(
1124                        (a - bf).abs() < 1e-6,
1125                        "{name} m={m} idx={i}: f32 {a} vs f64 {b}"
1126                    );
1127                }
1128            }
1129        }
1130    }
1131
1132    #[test]
1133    fn kaiser_f32_matches_kaiser() {
1134        let m = 20;
1135        let beta = 8.6;
1136        let f32_arr = kaiser_f32(m, beta).unwrap();
1137        let f64_arr = kaiser(m, beta).unwrap();
1138        for (a, b) in f32_arr.iter().zip(f64_arr.iter()) {
1139            assert!((a - *b as f32).abs() < 1e-5);
1140        }
1141    }
1142
1143    // Helper: compare two slices within tolerance
1144    fn assert_close(actual: &[f64], expected: &[f64], tol: f64) {
1145        assert_eq!(
1146            actual.len(),
1147            expected.len(),
1148            "length mismatch: {} vs {}",
1149            actual.len(),
1150            expected.len()
1151        );
1152        for (i, (a, e)) in actual.iter().zip(expected.iter()).enumerate() {
1153            assert!(
1154                (a - e).abs() <= tol,
1155                "element {i}: {a} vs {e} (diff = {})",
1156                (a - e).abs()
1157            );
1158        }
1159    }
1160
1161    // -----------------------------------------------------------------------
1162    // Edge cases: M=0 and M=1
1163    // -----------------------------------------------------------------------
1164
1165    #[test]
1166    fn bartlett_m0() {
1167        let w = bartlett(0).unwrap();
1168        assert_eq!(w.shape(), &[0]);
1169    }
1170
1171    #[test]
1172    fn bartlett_m1() {
1173        let w = bartlett(1).unwrap();
1174        assert_eq!(w.as_slice().unwrap(), &[1.0]);
1175    }
1176
1177    #[test]
1178    fn blackman_m0() {
1179        let w = blackman(0).unwrap();
1180        assert_eq!(w.shape(), &[0]);
1181    }
1182
1183    #[test]
1184    fn blackman_m1() {
1185        let w = blackman(1).unwrap();
1186        assert_eq!(w.as_slice().unwrap(), &[1.0]);
1187    }
1188
1189    #[test]
1190    fn hamming_m0() {
1191        let w = hamming(0).unwrap();
1192        assert_eq!(w.shape(), &[0]);
1193    }
1194
1195    #[test]
1196    fn hamming_m1() {
1197        let w = hamming(1).unwrap();
1198        assert_eq!(w.as_slice().unwrap(), &[1.0]);
1199    }
1200
1201    #[test]
1202    fn hanning_m0() {
1203        let w = hanning(0).unwrap();
1204        assert_eq!(w.shape(), &[0]);
1205    }
1206
1207    #[test]
1208    fn hanning_m1() {
1209        let w = hanning(1).unwrap();
1210        assert_eq!(w.as_slice().unwrap(), &[1.0]);
1211    }
1212
1213    #[test]
1214    fn kaiser_m0() {
1215        let w = kaiser(0, 14.0).unwrap();
1216        assert_eq!(w.shape(), &[0]);
1217    }
1218
1219    #[test]
1220    fn kaiser_m1() {
1221        let w = kaiser(1, 14.0).unwrap();
1222        assert_eq!(w.as_slice().unwrap(), &[1.0]);
1223    }
1224
1225    #[test]
1226    fn kaiser_negative_beta() {
1227        // Negative beta should produce the same result as positive beta (I0 is even)
1228        let pos = kaiser(5, 1.0).unwrap();
1229        let neg = kaiser(5, -1.0).unwrap();
1230        assert_eq!(pos.as_slice().unwrap(), neg.as_slice().unwrap());
1231    }
1232
1233    #[test]
1234    fn kaiser_nan_beta() {
1235        assert!(kaiser(5, f64::NAN).is_err());
1236    }
1237
1238    // -----------------------------------------------------------------------
1239    // AC-1: bartlett(5) matches np.bartlett(5) to within 4 ULPs
1240    // -----------------------------------------------------------------------
1241    // np.bartlett(5) = [0.0, 0.5, 1.0, 0.5, 0.0]
1242    #[test]
1243    fn bartlett_5_ac1() {
1244        let w = bartlett(5).unwrap();
1245        let expected = [0.0, 0.5, 1.0, 0.5, 0.0];
1246        assert_close(w.as_slice().unwrap(), &expected, 1e-14);
1247    }
1248
1249    // -----------------------------------------------------------------------
1250    // AC-2: kaiser(5, 14.0) matches np.kaiser(5, 14.0) to within 4 ULPs
1251    // -----------------------------------------------------------------------
1252    // np.kaiser(5, 14.0) = [7.72686684e-06, 1.64932188e-01, 1.0, 1.64932188e-01, 7.72686684e-06]
1253    #[test]
1254    fn kaiser_5_14_ac2() {
1255        let w = kaiser(5, 14.0).unwrap();
1256        let s = w.as_slice().unwrap();
1257        assert_eq!(s.len(), 5);
1258        // NumPy reference values (verified via np.kaiser(5, 14.0))
1259        let expected: [f64; 5] = [
1260            7.72686684e-06,
1261            1.64932188e-01,
1262            1.0,
1263            1.64932188e-01,
1264            7.72686684e-06,
1265        ];
1266        // Bessel polynomial approximation has limited precision (~6 digits)
1267        for (i, (&a, &e)) in s.iter().zip(expected.iter()).enumerate() {
1268            let tol = if e.abs() < 1e-4 { 1e-8 } else { 1e-6 };
1269            assert!(
1270                (a - e).abs() <= tol,
1271                "kaiser element {i}: {a} vs {e} (diff = {})",
1272                (a - e).abs()
1273            );
1274        }
1275    }
1276
1277    // -----------------------------------------------------------------------
1278    // AC-3: All 5 window functions return correct length and match NumPy fixtures
1279    // -----------------------------------------------------------------------
1280
1281    // np.blackman(5) = [-1.38777878e-17, 3.40000000e-01, 1.00000000e+00, 3.40000000e-01, -1.38777878e-17]
1282    #[test]
1283    fn blackman_5() {
1284        let w = blackman(5).unwrap();
1285        assert_eq!(w.shape(), &[5]);
1286        let s = w.as_slice().unwrap();
1287        let expected = [
1288            -1.38777878e-17,
1289            3.40000000e-01,
1290            1.00000000e+00,
1291            3.40000000e-01,
1292            -1.38777878e-17,
1293        ];
1294        assert_close(s, &expected, 1e-10);
1295    }
1296
1297    // np.hamming(5) = [0.08, 0.54, 1.0, 0.54, 0.08]
1298    #[test]
1299    fn hamming_5() {
1300        let w = hamming(5).unwrap();
1301        assert_eq!(w.shape(), &[5]);
1302        let s = w.as_slice().unwrap();
1303        let expected = [0.08, 0.54, 1.0, 0.54, 0.08];
1304        assert_close(s, &expected, 1e-14);
1305    }
1306
1307    // np.hanning(5) = [0.0, 0.5, 1.0, 0.5, 0.0]
1308    #[test]
1309    fn hanning_5() {
1310        let w = hanning(5).unwrap();
1311        assert_eq!(w.shape(), &[5]);
1312        let s = w.as_slice().unwrap();
1313        let expected = [0.0, 0.5, 1.0, 0.5, 0.0];
1314        assert_close(s, &expected, 1e-14);
1315    }
1316
1317    // Larger window: np.bartlett(12)
1318    #[test]
1319    fn bartlett_12() {
1320        let w = bartlett(12).unwrap();
1321        assert_eq!(w.shape(), &[12]);
1322        let s = w.as_slice().unwrap();
1323        // First and last should be 0, peak near center
1324        assert!((s[0] - 0.0).abs() < 1e-14);
1325        assert!((s[11] - 0.0).abs() < 1e-14);
1326        // Symmetric
1327        for i in 0..6 {
1328            assert!((s[i] - s[11 - i]).abs() < 1e-14, "symmetry at {i}");
1329        }
1330    }
1331
1332    // Symmetry test for all windows
1333    #[test]
1334    fn all_windows_symmetric() {
1335        let m = 7;
1336        let windows: Vec<Array<f64, Ix1>> = vec![
1337            bartlett(m).unwrap(),
1338            blackman(m).unwrap(),
1339            hamming(m).unwrap(),
1340            hanning(m).unwrap(),
1341            kaiser(m, 5.0).unwrap(),
1342        ];
1343        for (idx, w) in windows.iter().enumerate() {
1344            let s = w.as_slice().unwrap();
1345            for i in 0..m / 2 {
1346                assert!(
1347                    (s[i] - s[m - 1 - i]).abs() < 1e-12,
1348                    "window {idx} not symmetric at {i}"
1349                );
1350            }
1351        }
1352    }
1353
1354    // All windows peak at center
1355    #[test]
1356    fn all_windows_peak_at_center() {
1357        let m = 9;
1358        let windows: Vec<Array<f64, Ix1>> = vec![
1359            bartlett(m).unwrap(),
1360            blackman(m).unwrap(),
1361            hamming(m).unwrap(),
1362            hanning(m).unwrap(),
1363            kaiser(m, 5.0).unwrap(),
1364        ];
1365        for (idx, w) in windows.iter().enumerate() {
1366            let s = w.as_slice().unwrap();
1367            let center = s[m / 2];
1368            assert!(
1369                (center - 1.0).abs() < 1e-10,
1370                "window {idx} center = {center}, expected 1.0"
1371            );
1372        }
1373    }
1374
1375    // Kaiser with beta=0 should be a rectangular window (all ones)
1376    #[test]
1377    fn kaiser_beta_zero() {
1378        let w = kaiser(5, 0.0).unwrap();
1379        let s = w.as_slice().unwrap();
1380        for &v in s {
1381            assert!((v - 1.0).abs() < 1e-10, "expected 1.0, got {v}");
1382        }
1383    }
1384
1385    #[test]
1386    fn bessel_i0_scalar_zero() {
1387        assert!((bessel_i0_scalar::<f64>(0.0) - 1.0).abs() < 1e-6);
1388    }
1389
1390    #[test]
1391    fn bessel_i0_scalar_known() {
1392        // Issue #296: tighten the tolerances. The Abramowitz & Stegun
1393        // approximation is rated ~1e-7 in the polynomial branch and
1394        // ~1e-7 in the asymptotic branch per its documentation.
1395        // I0(1) ~ 1.2660658777520082
1396        assert!((bessel_i0_scalar::<f64>(1.0) - 1.266_065_877_752_008_2).abs() < 1e-7);
1397        // I0(5) ~ 27.239871823604443 (asymptotic branch)
1398        assert!((bessel_i0_scalar::<f64>(5.0) - 27.239_871_823_604_443).abs() < 1e-5);
1399        // I0(10) ~ 2815.7166284662544
1400        let expected_i0_10 = 2_815.716_628_466_254;
1401        assert!((bessel_i0_scalar::<f64>(10.0) - expected_i0_10).abs() / expected_i0_10 < 1e-5);
1402    }
1403
1404    // ----- Edge-length window coverage (#295) -----
1405
1406    #[test]
1407    fn bartlett_m2_is_zeros() {
1408        // Bartlett of length 2: w(n) = 1 - |2n/(M-1) - 1|
1409        // n=0: 1 - |0 - 1| = 0
1410        // n=1: 1 - |2 - 1| = 0
1411        let w = bartlett(2).unwrap();
1412        assert_eq!(w.as_slice().unwrap(), &[0.0, 0.0]);
1413    }
1414
1415    #[test]
1416    fn hanning_m2_is_zeros() {
1417        // Hann of length 2: w(n) = 0.5 * (1 - cos(2*pi*n/(M-1)))
1418        // n=0: 0.5 * (1 - cos(0)) = 0
1419        // n=1: 0.5 * (1 - cos(2*pi)) = 0
1420        let w = hanning(2).unwrap();
1421        let s = w.as_slice().unwrap();
1422        assert!(s[0].abs() < 1e-14);
1423        assert!(s[1].abs() < 1e-14);
1424    }
1425
1426    #[test]
1427    fn bartlett_even_length_is_symmetric() {
1428        // Even-length windows have no single "center"; symmetry holds
1429        // around the midpoint between indices (M/2 - 1) and (M/2).
1430        for &m in &[4usize, 6, 8, 10] {
1431            let w = bartlett(m).unwrap();
1432            let s = w.as_slice().unwrap();
1433            for i in 0..m / 2 {
1434                assert!(
1435                    (s[i] - s[m - 1 - i]).abs() < 1e-14,
1436                    "bartlett({m}) not symmetric at {i}"
1437                );
1438            }
1439        }
1440    }
1441
1442    #[test]
1443    fn blackman_even_length_is_symmetric() {
1444        for &m in &[4usize, 6, 8, 10] {
1445            let w = blackman(m).unwrap();
1446            let s = w.as_slice().unwrap();
1447            for i in 0..m / 2 {
1448                assert!(
1449                    (s[i] - s[m - 1 - i]).abs() < 1e-14,
1450                    "blackman({m}) not symmetric at {i}"
1451                );
1452            }
1453        }
1454    }
1455
1456    #[test]
1457    fn kaiser_large_beta_errors() {
1458        // Issue #294: beta above the f64 safe range should error,
1459        // not silently produce NaN.
1460        assert!(kaiser(8, 800.0).is_err());
1461        assert!(kaiser(8, 1000.0).is_err());
1462        // 700 is still safe.
1463        assert!(kaiser(8, 700.0).is_ok());
1464    }
1465
1466    // =======================================================================
1467    // SciPy-extended windows (cosine, exponential, gaussian, general_cosine,
1468    // general_hamming, nuttall, parzen, taylor, tukey).
1469    // =======================================================================
1470
1471    fn close(a: f64, b: f64, tol: f64) -> bool {
1472        (a - b).abs() < tol
1473    }
1474
1475    // ----- cosine -----
1476
1477    #[test]
1478    fn cosine_length_and_endpoints() {
1479        // cosine(M) = sin(pi(n + 0.5)/M); endpoints are sin(pi/2M) and
1480        // sin(pi(M-0.5)/M) — both small but nonzero.
1481        let w = cosine(8).unwrap();
1482        let s = w.as_slice().unwrap();
1483        assert_eq!(s.len(), 8);
1484        // Symmetric.
1485        for i in 0..4 {
1486            assert!(close(s[i], s[7 - i], 1e-14));
1487        }
1488    }
1489
1490    #[test]
1491    fn cosine_m1_and_m0() {
1492        assert_eq!(cosine(0).unwrap().shape(), &[0]);
1493        assert_eq!(cosine(1).unwrap().as_slice().unwrap(), &[1.0]);
1494    }
1495
1496    // ----- exponential -----
1497
1498    #[test]
1499    fn exponential_centred_default() {
1500        // tau=1, default centre = (M-1)/2 = 3.5 for M=8 — midpoint
1501        // between samples 3 and 4.
1502        let w = exponential(8, None, 1.0).unwrap();
1503        let s = w.as_slice().unwrap();
1504        // Symmetric.
1505        for i in 0..4 {
1506            assert!(close(s[i], s[7 - i], 1e-14));
1507        }
1508        // Centre samples have largest values.
1509        let centre_max = s[3].max(s[4]);
1510        for &v in s {
1511            assert!(v <= centre_max + 1e-14);
1512        }
1513    }
1514
1515    #[test]
1516    fn exponential_rejects_nonpositive_tau() {
1517        assert!(exponential(8, None, 0.0).is_err());
1518        assert!(exponential(8, None, -1.0).is_err());
1519        assert!(exponential(8, None, f64::NAN).is_err());
1520    }
1521
1522    // ----- gaussian -----
1523
1524    #[test]
1525    fn gaussian_centre_is_one() {
1526        // For odd M, centre sample is exactly 1.
1527        let w = gaussian(11, 2.0).unwrap();
1528        let s = w.as_slice().unwrap();
1529        assert!(close(s[5], 1.0, 1e-14));
1530    }
1531
1532    #[test]
1533    fn gaussian_known_value() {
1534        // gaussian(7, 1) at n=4: z = (4 - 3) / 1 = 1, exp(-0.5) = 0.6065...
1535        let w = gaussian(7, 1.0).unwrap();
1536        let s = w.as_slice().unwrap();
1537        assert!(close(s[4], (-0.5_f64).exp(), 1e-14));
1538    }
1539
1540    #[test]
1541    fn gaussian_rejects_nonpositive_std() {
1542        assert!(gaussian(8, 0.0).is_err());
1543        assert!(gaussian(8, -1.0).is_err());
1544    }
1545
1546    // ----- general_cosine -----
1547
1548    #[test]
1549    fn general_cosine_with_hann_coeffs_matches_hann() {
1550        // [0.5, 0.5] → 0.5 - 0.5 cos(2πn/(M-1)) = Hann.
1551        let m = 9;
1552        let gc = general_cosine(m, &[0.5, 0.5]).unwrap();
1553        let hn = hanning(m).unwrap();
1554        for i in 0..m {
1555            assert!(close(
1556                gc.as_slice().unwrap()[i],
1557                hn.as_slice().unwrap()[i],
1558                1e-14,
1559            ));
1560        }
1561    }
1562
1563    #[test]
1564    fn general_cosine_with_blackman_coeffs_matches_blackman() {
1565        // [0.42, 0.5, 0.08] → classical Blackman.
1566        let m = 9;
1567        let gc = general_cosine(m, &[0.42, 0.5, 0.08]).unwrap();
1568        let bk = blackman(m).unwrap();
1569        for i in 0..m {
1570            assert!(close(
1571                gc.as_slice().unwrap()[i],
1572                bk.as_slice().unwrap()[i],
1573                1e-12,
1574            ));
1575        }
1576    }
1577
1578    #[test]
1579    fn general_cosine_rejects_empty_coeffs() {
1580        assert!(general_cosine(8, &[]).is_err());
1581    }
1582
1583    // ----- general_hamming -----
1584
1585    #[test]
1586    fn general_hamming_alpha_half_matches_hann() {
1587        let m = 9;
1588        let gh = general_hamming(m, 0.5).unwrap();
1589        let hn = hanning(m).unwrap();
1590        for i in 0..m {
1591            assert!(close(
1592                gh.as_slice().unwrap()[i],
1593                hn.as_slice().unwrap()[i],
1594                1e-14,
1595            ));
1596        }
1597    }
1598
1599    #[test]
1600    fn general_hamming_alpha_054_matches_hamming() {
1601        let m = 9;
1602        let gh = general_hamming(m, 0.54).unwrap();
1603        let hm = hamming(m).unwrap();
1604        for i in 0..m {
1605            assert!(close(
1606                gh.as_slice().unwrap()[i],
1607                hm.as_slice().unwrap()[i],
1608                1e-14,
1609            ));
1610        }
1611    }
1612
1613    // ----- nuttall -----
1614
1615    #[test]
1616    fn nuttall_endpoints_are_small() {
1617        // Nuttall is engineered to have very low sidelobes; endpoints
1618        // are O(1e-3) for typical M.
1619        let w = nuttall(64).unwrap();
1620        let s = w.as_slice().unwrap();
1621        assert!(s[0].abs() < 1e-2);
1622        assert!(s[s.len() - 1].abs() < 1e-2);
1623    }
1624
1625    #[test]
1626    fn nuttall_is_symmetric() {
1627        let m = 33;
1628        let w = nuttall(m).unwrap();
1629        let s = w.as_slice().unwrap();
1630        for i in 0..m / 2 {
1631            assert!(close(s[i], s[m - 1 - i], 1e-14));
1632        }
1633    }
1634
1635    // ----- parzen -----
1636
1637    #[test]
1638    fn parzen_endpoints_are_zero() {
1639        let w = parzen(16).unwrap();
1640        let s = w.as_slice().unwrap();
1641        // Outer edge of Parzen has very small value (cubic falloff).
1642        assert!(s[0].abs() < 1e-2);
1643        assert!(s[s.len() - 1].abs() < 1e-2);
1644    }
1645
1646    #[test]
1647    fn parzen_centre_is_one() {
1648        // For odd M, centre is at sample (M-1)/2; r = 0 → w = 1.
1649        let w = parzen(13).unwrap();
1650        let s = w.as_slice().unwrap();
1651        assert!(close(s[6], 1.0, 1e-14));
1652    }
1653
1654    #[test]
1655    fn parzen_is_symmetric() {
1656        let m = 21;
1657        let w = parzen(m).unwrap();
1658        let s = w.as_slice().unwrap();
1659        for i in 0..m / 2 {
1660            assert!(close(s[i], s[m - 1 - i], 1e-14));
1661        }
1662    }
1663
1664    // ----- taylor -----
1665
1666    #[test]
1667    fn taylor_default_normalised_centre_is_one() {
1668        // With norm=true the centre sample is 1.0.
1669        let w = taylor(33, 4, 30.0, true).unwrap();
1670        let s = w.as_slice().unwrap();
1671        assert!(close(s[16], 1.0, 1e-12));
1672    }
1673
1674    #[test]
1675    fn taylor_is_symmetric() {
1676        let m = 33;
1677        let w = taylor(m, 4, 30.0, true).unwrap();
1678        let s = w.as_slice().unwrap();
1679        for i in 0..m / 2 {
1680            assert!(close(s[i], s[m - 1 - i], 1e-12));
1681        }
1682    }
1683
1684    #[test]
1685    fn taylor_rejects_nbar_zero() {
1686        assert!(taylor(8, 0, 30.0, true).is_err());
1687        assert!(taylor(8, 4, f64::NAN, true).is_err());
1688    }
1689
1690    // ----- tukey -----
1691
1692    #[test]
1693    fn tukey_alpha_zero_is_rectangular() {
1694        let w = tukey(8, 0.0).unwrap();
1695        for &v in w.as_slice().unwrap() {
1696            assert!(close(v, 1.0, 1e-14));
1697        }
1698    }
1699
1700    #[test]
1701    fn tukey_alpha_one_matches_hann() {
1702        // alpha=1 → fully Hann-shaped.
1703        let m = 9;
1704        let tk = tukey(m, 1.0).unwrap();
1705        let hn = hanning(m).unwrap();
1706        for i in 0..m {
1707            assert!(close(
1708                tk.as_slice().unwrap()[i],
1709                hn.as_slice().unwrap()[i],
1710                1e-12,
1711            ));
1712        }
1713    }
1714
1715    #[test]
1716    fn tukey_centre_is_one() {
1717        let m = 21;
1718        let w = tukey(m, 0.5).unwrap();
1719        // Wide flat-top region: middle must be 1.
1720        assert!(close(w.as_slice().unwrap()[m / 2], 1.0, 1e-14));
1721    }
1722
1723    #[test]
1724    fn tukey_rejects_invalid_alpha() {
1725        assert!(tukey(8, -0.1).is_err());
1726        assert!(tukey(8, 1.1).is_err());
1727        assert!(tukey(8, f64::NAN).is_err());
1728    }
1729}
1730
1731#[cfg(test)]
1732mod debug_test_removed {
1733    // Debug test scaffolding for #740 chebwin development;
1734    // removed once the closed-form output matched scipy.
1735}