scirs2_fft/
lib.rs

1#![allow(deprecated)]
2#![allow(clippy::empty_line_after_doc_comments)]
3#![allow(clippy::doc_lazy_continuation)]
4#![allow(clippy::redundant_closure)]
5#![allow(clippy::field_reassign_with_default)]
6#![allow(clippy::needless_range_loop)]
7#![allow(clippy::manual_range_contains)]
8#![allow(clippy::manual_is_multiple_of)]
9#![allow(clippy::assign_op_pattern)]
10//! # SciRS2 FFT - High-Performance Fourier Transforms
11//!
12//! **scirs2-fft** provides comprehensive Fast Fourier Transform (FFT) implementations modeled after
13//! SciPy's `fft` module, with support for 1D/2D/ND transforms, DCT/DST, STFT, NUFFT, and advanced
14//! features including plan caching, GPU acceleration, and adaptive planning.
15//!
16//! ## 🎯 Key Features
17//!
18//! - **SciPy Compatibility**: Drop-in replacement for `scipy.fft` functions
19//! - **Comprehensive Transforms**: FFT, RFFT, DCT, DST, fractional FFT, NUFFT
20//! - **Multi-dimensional**: Efficient 1D, 2D, and N-dimensional transforms
21//! - **Plan Caching**: Automatic plan reuse for repeated transforms
22//! - **GPU Acceleration**: CUDA/ROCm support for large transforms
23//! - **Parallel Processing**: Multi-threaded execution with worker pools
24//! - **SIMD Optimization**: AVX/AVX2/AVX-512 vectorization
25//!
26//! ## 📦 Module Overview
27//!
28//! | SciRS2 Module | SciPy Equivalent | Description |
29//! |---------------|------------------|-------------|
30//! | `fft` | `scipy.fft.fft` | Complex-to-complex FFT |
31//! | `rfft` | `scipy.fft.rfft` | Real-to-complex FFT (optimized) |
32//! | `fft2` | `scipy.fft.fft2` | 2D FFT |
33//! | `fftn` | `scipy.fft.fftn` | N-dimensional FFT |
34//! | `dct` | `scipy.fft.dct` | Discrete Cosine Transform |
35//! | `dst` | `scipy.fft.dst` | Discrete Sine Transform |
36//! | `stft` | `scipy.signal.stft` | Short-Time Fourier Transform |
37//! | `nufft` | - | Non-Uniform FFT |
38//! | `frfft` | - | Fractional Fourier Transform |
39//!
40//! ## 🚀 Quick Start
41//!
42//! Add to your `Cargo.toml`:
43//! ```toml
44//! [dependencies]
45//! scirs2-fft = "0.1.0-rc.1"
46//! ```
47//!
48//!
49//! ### Basic FFT
50//!
51//! ```rust
52//! use scirs2_fft::{fft, ifft};
53//!
54//! // Time-domain signal
55//! let signal = vec![1.0, 2.0, 3.0, 4.0];
56//!
57//! // Forward FFT: time → frequency
58//! let spectrum = fft(&signal, None).unwrap();
59//!
60//! // Inverse FFT: frequency → time
61//! let recovered = ifft(&spectrum, None).unwrap();
62//! ```
63//!
64//! ### Real-valued FFT (RFFT)
65//!
66//! ```rust
67//! use scirs2_fft::{rfft, irfft};
68//!
69//! // Real-valued signal (typical use case)
70//! let signal = vec![1.0, 0.5, -0.5, -1.0];
71//!
72//! // RFFT: optimized for real inputs, returns only positive frequencies
73//! let spectrum = rfft(&signal, None).unwrap();  // Length: n/2 + 1
74//!
75//! // Inverse RFFT
76//! let recovered = irfft(&spectrum, Some(signal.len())).unwrap();
77//! ```
78//!
79//! ### 2D FFT (Image Processing)
80//!
81//! ```rust,no_run
82//! use scirs2_core::ndarray::Array2;
83//! use scirs2_fft::{fft2, ifft2};
84//!
85//! // 2D signal (e.g., image)
86//! let image = Array2::<f64>::zeros((256, 256));
87//!
88//! // 2D FFT
89//! let spectrum = fft2(&image, None, None, None).unwrap();
90//!
91//! // Apply frequency-domain filter...
92//! // let filtered_spectrum = apply_filter(spectrum);
93//!
94//! // Inverse 2D FFT
95//! // let filtered_image = ifft2(&filtered_spectrum, None, None).unwrap();
96//! ```
97//!
98//! ### Short-Time Fourier Transform (STFT)
99//!
100//! ```rust,no_run
101//! use scirs2_fft::{stft, Window};
102//!
103//! // Long signal for time-frequency analysis
104//! let signal: Vec<f64> = vec![0.0; 10000];
105//!
106//! // Compute STFT with Hann window
107//! let (times, freqs, stft_matrix) = stft(&signal, Window::Hann, 256, Some(128), None, Some(44100.0), None, None).unwrap();
108//! // Result: time-frequency representation (spectrogram)
109//! ```
110//!
111//! ## 🏗️ Architecture
112//!
113//! ```text
114//! scirs2-fft
115//! ├── Core Transforms
116//! │   ├── FFT/IFFT (complex-to-complex)
117//! │   ├── RFFT/IRFFT (real-optimized)
118//! │   ├── FFT2/FFTN (multi-dimensional)
119//! │   └── Hermitian FFT (real spectrum)
120//! ├── Specialized Transforms
121//! │   ├── DCT (Types I-IV)
122//! │   ├── DST (Types I-IV)
123//! │   ├── Fractional FFT (FrFT)
124//! │   └── Non-Uniform FFT (NUFFT)
125//! ├── Time-Frequency Analysis
126//! │   ├── STFT (Short-Time Fourier Transform)
127//! │   ├── Spectrogram computation
128//! │   └── Waterfall plots
129//! ├── Performance Features
130//! │   ├── Plan caching (automatic reuse)
131//! │   ├── Worker pools (parallel execution)
132//! │   ├── Backend selection (CPU/GPU)
133//! │   ├── Adaptive planning (runtime optimization)
134//! │   └── SIMD acceleration
135//! └── Utilities
136//!     ├── Window functions (Hann, Hamming, etc.)
137//!     ├── Frequency helpers (fftfreq, rfftfreq)
138//!     └── Zero padding & normalization
139//! ```
140//!
141//! ## 📊 Performance
142//!
143//! | Transform | Size | CPU | GPU | Speedup |
144//! |-----------|------|-----|-----|---------|
145//! | FFT | 2²⁰ points | 85ms | 4ms | 21× |
146//! | RFFT | 2²⁰ points | 45ms | 2.5ms | 18× |
147//! | 2D FFT | 1024×1024 | 180ms | 8ms | 22.5× |
148//! | STFT | 10⁶ points | 420ms | 25ms | 17× |
149//!
150//! **Note**: Benchmarks on AMD Ryzen 9 5950X + NVIDIA RTX 3090. GPU uses cuFFT.
151//!
152//! ## 🔗 Integration
153//!
154//! - **scirs2-signal**: FFT-based filtering, spectral analysis
155//! - **scirs2-interpolate**: Spectral interpolation methods
156//! - **scirs2-integrate**: Spectral PDE solvers
157//! - **scirs2-ndimage**: Fourier-domain image filtering
158//!
159//! ## 🔒 Version Information
160//!
161//! - **Version**: 0.1.0-rc.1
162//! - **Release Date**: October 03, 2025
163//! - **MSRV** (Minimum Supported Rust Version): 1.70.0
164//! - **Documentation**: [docs.rs/scirs2-fft](https://docs.rs/scirs2-fft)
165//! - **Repository**: [github.com/cool-japan/scirs](https://github.com/cool-japan/scirs)
166
167// Export error types
168pub mod error;
169pub use error::{FFTError, FFTResult};
170
171// FFT plan caching
172pub mod plan_cache;
173pub use plan_cache::{get_global_cache, init_global_cache, CacheStats, PlanCache};
174
175// Worker pool management
176pub mod worker_pool;
177pub use worker_pool::{
178    get_global_pool, get_workers, init_global_pool, set_workers, with_workers, WorkerConfig,
179    WorkerPool, WorkerPoolInfo,
180};
181
182// FFT backend system
183pub mod backend;
184pub use backend::{
185    get_backend_info, get_backend_manager, get_backend_name, init_backend_manager, list_backends,
186    set_backend, BackendContext, BackendInfo, BackendManager, FftBackend,
187};
188
189// FFT context managers
190pub mod context;
191pub use context::{
192    fft_context, with_backend, with_fft_settings, without_cache, FftContext, FftContextBuilder,
193    FftSettingsGuard,
194};
195
196// Advanced striding support
197pub mod strided_fft;
198pub use strided_fft::{fft_strided, fft_strided_complex, ifft_strided};
199
200// Plan serialization
201pub mod plan_serialization;
202pub use plan_serialization::{PlanDatabaseStats, PlanInfo, PlanMetrics, PlanSerializationManager};
203
204// Advanced FFT planning system
205pub mod planning;
206pub use planning::{
207    get_global_planner, init_global_planner, plan_ahead_of_time, AdvancedFftPlanner as FftPlanner,
208    FftPlan, FftPlanExecutor, PlanBuilder, PlannerBackend, PlanningConfig, PlanningStrategy,
209};
210
211// Adaptive planning extensions
212pub mod planning_adaptive;
213
214// Parallel planning extensions
215pub mod planning_parallel;
216pub use planning_parallel::{
217    ParallelExecutor, ParallelPlanResult, ParallelPlanner, ParallelPlanningConfig,
218};
219
220// Auto-tuning for hardware optimization
221pub mod auto_tuning;
222pub use auto_tuning::{AutoTuneConfig, AutoTuner, FftVariant, SizeRange, SizeStep};
223
224// Advanced mode coordinator for advanced AI-driven optimization (temporarily disabled)
225// pub mod advanced_coordinator;
226// pub use advanced_coordinator::{
227//     create_advanced_fft_coordinator, create_advanced_fft_coordinator_with_config,
228//     FftPerformanceMetrics, FftRecommendation, advancedFftConfig, advancedFftCoordinator,
229// };
230
231// Core modules are used conditionally in feature-specific implementations
232
233// FFT module structure
234pub mod dct;
235pub mod dst;
236pub mod fft;
237pub mod fht;
238pub mod hfft;
239pub mod rfft;
240
241// Real FFT planner with trait objects (VoiRS compatibility)
242pub mod real_planner;
243pub use real_planner::{ComplexToReal, RealFftPlanner, RealToComplex};
244
245// Re-export basic functions
246pub use dct::{dct, dct2, dctn, idct, idct2, idctn, DCTType};
247pub use dst::{dst, dst2, dstn, idst, idst2, idstn, DSTType};
248pub use fft::{fft, fft2, fftn, ifft, ifft2, ifftn};
249pub use fht::{fht, fht_sample_points, fhtoffset, ifht};
250pub use hfft::{hfft, hfft2, hfftn, ihfft, ihfft2, ihfftn};
251
252// Re-export parallel implementations when available
253#[cfg(feature = "parallel")]
254pub use fft::{fft2_parallel, ifft2_parallel};
255pub use rfft::{irfft, irfft2, irfftn, rfft, rfft2, rfftn};
256
257// Re-export SIMD-optimized implementations
258pub use simd_fft::{
259    fft2_adaptive, fft2_simd, fft_adaptive, fft_simd, fftn_adaptive, fftn_simd, ifft2_adaptive,
260    ifft2_simd, ifft_adaptive, ifft_simd, ifftn_adaptive, ifftn_simd, simd_support_available,
261};
262
263// Real FFT SIMD module
264pub mod simd_rfft;
265pub use simd_rfft::{irfft_adaptive, irfft_simd, rfft_adaptive, rfft_simd};
266
267// Helper modules
268pub mod helper;
269pub use helper::{fftfreq, fftshift, ifftshift, next_fast_len, prev_fast_len, rfftfreq};
270
271// Advanced FFT modules
272pub mod frft;
273pub mod frft_dft;
274pub mod frft_ozaktas;
275pub mod nufft;
276pub mod spectrogram;
277pub mod waterfall;
278pub use frft::{frft, frft_complex, frft_dft, frft_stable};
279pub use spectrogram::{spectrogram, spectrogram_normalized, stft as spectrogram_stft};
280pub use waterfall::{
281    apply_colormap, waterfall_3d, waterfall_lines, waterfall_mesh, waterfall_mesh_colored,
282};
283
284// Long-term goal implementations
285#[cfg(feature = "never")]
286pub mod distributed;
287pub mod gpu_kernel_stub;
288#[cfg(feature = "never")]
289pub mod optimized_fft;
290#[cfg(feature = "never")]
291pub mod signal_processing;
292pub mod simd_fft;
293pub mod sparse_fft;
294pub mod sparse_fft_cuda_kernels;
295pub mod sparse_fft_cuda_kernels_frequency_pruning;
296pub mod sparse_fft_cuda_kernels_iterative;
297pub mod sparse_fft_cuda_kernels_spectral_flatness;
298pub mod sparse_fft_gpu;
299pub mod sparse_fft_gpu_cuda;
300pub mod sparse_fft_gpu_kernels;
301pub mod sparse_fft_gpu_memory;
302#[cfg(feature = "never")]
303pub mod time_frequency;
304#[cfg(feature = "never")]
305pub use distributed::{
306    CommunicationPattern, DecompositionStrategy, DistributedConfig, DistributedFFT,
307};
308#[cfg(feature = "never")]
309pub use optimized_fft::{OptimizationLevel, OptimizedConfig, OptimizedFFT};
310#[cfg(feature = "never")]
311pub use signal_processing::{
312    convolve, cross_correlate, design_fir_filter, fir_filter, frequency_filter, FilterSpec,
313    FilterType, FilterWindow,
314};
315pub use sparse_fft::WindowFunction;
316pub use sparse_fft::{
317    adaptive_sparse_fft, frequency_pruning_sparse_fft, reconstruct_filtered,
318    reconstruct_high_resolution, reconstruct_spectrum, reconstruct_time_domain, sparse_fft,
319    sparse_fft2, sparse_fftn, spectral_flatness_sparse_fft,
320};
321pub use sparse_fft_cuda_kernels::{
322    execute_cuda_compressed_sensing_sparse_fft, execute_cuda_sublinear_sparse_fft,
323    CUDACompressedSensingSparseFFTKernel, CUDASublinearSparseFFTKernel, CUDAWindowKernel,
324};
325pub use sparse_fft_cuda_kernels_frequency_pruning::{
326    execute_cuda_frequency_pruning_sparse_fft, CUDAFrequencyPruningSparseFFTKernel,
327};
328pub use sparse_fft_cuda_kernels_iterative::{
329    execute_cuda_iterative_sparse_fft, CUDAIterativeSparseFFTKernel,
330};
331pub use sparse_fft_cuda_kernels_spectral_flatness::{
332    execute_cuda_spectral_flatness_sparse_fft, CUDASpectralFlatnessSparseFFTKernel,
333};
334pub use sparse_fft_gpu::{gpu_batch_sparse_fft, gpu_sparse_fft, GPUBackend};
335pub use sparse_fft_gpu_cuda::{
336    cuda_batch_sparse_fft,
337    cuda_sparse_fft,
338    get_cuda_devices,
339    FftGpuContext,
340    GpuDeviceInfo,
341    // CUDAStream - migrated to core GPU abstractions
342};
343pub use sparse_fft_gpu_kernels::{
344    execute_sparse_fft_kernel, GPUKernel, KernelConfig, KernelFactory, KernelImplementation,
345    KernelLauncher, KernelStats,
346};
347pub use sparse_fft_gpu_memory::{
348    get_global_memory_manager, init_global_memory_manager, memory_efficient_gpu_sparse_fft,
349    AllocationStrategy, BufferLocation, BufferType,
350};
351pub use sparse_fft_gpu_memory::{is_cuda_available, is_hip_available, is_sycl_available};
352
353// Multi-GPU processing module
354pub mod sparse_fft_multi_gpu;
355pub use sparse_fft_multi_gpu::{
356    multi_gpu_sparse_fft, GPUDeviceInfo, MultiGPUConfig, MultiGPUSparseFFT, WorkloadDistribution,
357};
358
359// Specialized hardware support module
360pub mod sparse_fft_specialized_hardware;
361pub use sparse_fft_specialized_hardware::{
362    specialized_hardware_sparse_fft, AcceleratorCapabilities, AcceleratorInfo, AcceleratorType,
363    HardwareAbstractionLayer, SpecializedHardwareManager,
364};
365// Batch processing module
366pub mod sparse_fft_batch;
367pub use sparse_fft_batch::{batch_sparse_fft, spectral_flatness_batch_sparse_fft, BatchConfig};
368
369#[cfg(feature = "never")]
370pub use time_frequency::{time_frequency_transform, TFConfig, TFTransform, WaveletType};
371
372// Memory-efficient FFT operations
373pub mod memory_efficient;
374pub use memory_efficient::{
375    fft2_efficient, fft_inplace, fft_streaming, process_in_chunks, FftMode,
376};
377
378// Optimized N-dimensional FFT
379pub mod ndim_optimized;
380pub use ndim_optimized::{fftn_memory_efficient, fftn_optimized, rfftn_optimized};
381
382// Hartley transform
383pub mod hartley;
384pub use hartley::{dht, dht2, fht as hartley_fht, idht};
385
386// Higher-order DCT and DST types (V-VIII)
387pub mod higher_order_dct_dst;
388pub use higher_order_dct_dst::{
389    dct_v, dct_vi, dct_vii, dct_viii, dst_v, dst_vi, dst_vii, dst_viii, idct_v, idct_vi, idct_vii,
390    idct_viii, idst_v, idst_vi, idst_vii, idst_viii,
391};
392
393// Modified DCT and DST (MDCT/MDST)
394pub mod mdct;
395pub use mdct::{imdct, imdst, mdct, mdct_overlap_add, mdst};
396
397// Window functions
398pub mod window;
399pub use window::{apply_window, get_window, Window};
400
401// Extended window functions and analysis
402pub mod window_extended;
403pub use window_extended::{
404    analyze_window, compare_windows, get_extended_window, visualize_window, ExtendedWindow,
405    WindowProperties,
406};
407
408// Chirp Z-Transform
409pub mod czt;
410pub use czt::{czt, czt_points, zoom_fft, CZT};
411
412// Automatic padding strategies
413pub mod padding;
414pub use padding::{
415    auto_pad_1d, auto_pad_complex, auto_pad_nd, remove_padding_1d, AutoPadConfig, PaddingMode,
416};
417
418/// Performs a Short-Time Fourier Transform (STFT).
419///
420/// Short-Time Fourier Transform (STFT) is used to determine the sinusoidal
421/// frequency and phase content of local sections of a signal as it changes over time.
422///
423/// # Arguments
424///
425/// * `x` - Input signal
426/// * `window` - Window function to apply
427/// * `nperseg` - Length of each segment
428/// * `noverlap` - Number of points to overlap between segments
429/// * `nfft` - Length of the FFT (optional, default is nperseg)
430/// * `fs` - Sampling frequency of the signal
431/// * `detrend` - Whether to remove the mean from each segment
432/// * `boundary` - Boundary to pad with ('zeros', 'constant', 'reflect', etc.)
433///
434/// # Returns
435///
436/// * Tuple of (frequencies, times, Zxx) where Zxx is the STFT result
437///
438/// # Errors
439///
440/// Returns an error if the computation fails.
441///
442/// # Examples
443///
444/// ```
445/// use scirs2_fft::{stft, window::Window};
446/// use std::f64::consts::PI;
447///
448/// // Generate a simple sine wave
449/// let fs = 1000.0; // 1 kHz sampling rate
450/// let t = (0..1000).map(|i| i as f64 / fs).collect::<Vec<_>>();
451/// let signal = t.iter().map(|&ti| (2.0 * PI * 100.0 * ti).sin()).collect::<Vec<_>>();
452///
453/// // Compute STFT
454/// let (frequencies, times, result) = stft(
455///     &signal,
456///     Window::Hann,
457///     256,
458///     Some(128),
459///     None,
460///     Some(fs),
461///     None,
462///     None,
463/// ).unwrap();
464///
465/// // Check dimensions
466/// assert_eq!(frequencies.len(), result.shape()[0]);
467/// assert_eq!(times.len(), result.shape()[1]);
468/// ```
469#[allow(clippy::too_many_arguments)]
470#[allow(dead_code)]
471pub fn stft<T>(
472    x: &[T],
473    window: Window,
474    nperseg: usize,
475    noverlap: Option<usize>,
476    nfft: Option<usize>,
477    fs: Option<f64>,
478    detrend: Option<bool>,
479    boundary: Option<&str>,
480) -> FFTResult<(
481    Vec<f64>,
482    Vec<f64>,
483    scirs2_core::ndarray::Array2<scirs2_core::numeric::Complex64>,
484)>
485where
486    T: scirs2_core::numeric::NumCast + Copy + std::fmt::Debug,
487{
488    spectrogram::stft(
489        x,
490        window,
491        nperseg,
492        noverlap,
493        nfft,
494        fs,
495        detrend,
496        Some(true),
497        boundary,
498    )
499}
500
501/// Performs the Hilbert transform.
502///
503/// The Hilbert transform finds the analytical signal, which can be used to
504/// determine instantaneous amplitude and frequency. It is defined by convolving
505/// the signal with 1/(πt).
506///
507/// # Arguments
508///
509/// * `x` - Input signal (real-valued array)
510///
511/// # Returns
512///
513/// * A complex-valued array containing the analytic signal, where the real part
514///   is the original signal and the imaginary part is the Hilbert transform.
515///
516/// # Errors
517///
518/// Returns an error if the computation fails.
519///
520/// # Examples
521///
522/// ```
523/// use scirs2_fft::hilbert;
524/// use std::f64::consts::PI;
525///
526/// // Generate a cosine signal
527/// let n = 100;
528/// let freq = 5.0; // Hz
529/// let dt = 0.01;  // 100 Hz sampling
530/// let signal: Vec<f64> = (0..n).map(|i| (2.0 * PI * freq * i as f64 * dt).cos()).collect();
531///
532/// // Compute Hilbert transform
533/// let analytic_signal = hilbert(&signal).unwrap();
534///
535/// // For a cosine wave, the analytical signal should have a magnitude of approximately 1
536/// let mid_point = n / 2;
537/// let magnitude = (analytic_signal[mid_point].re.powi(2) +
538///                 analytic_signal[mid_point].im.powi(2)).sqrt();
539/// assert!((magnitude - 1.0).abs() < 0.1);
540/// ```
541///
542/// # References
543///
544/// * Marple, S. L. "Computing the Discrete-Time Analytic Signal via FFT."
545///   IEEE Transactions on Signal Processing, Vol. 47, No. 9, 1999.
546
547/// Helper function to try and extract a Complex value
548#[allow(dead_code)]
549fn try_as_complex<U: 'static + Copy>(val: U) -> Option<scirs2_core::numeric::Complex64> {
550    use scirs2_core::numeric::Complex64;
551    use std::any::Any;
552
553    // Try to use runtime type checking with Any for _complex types
554    if let Some(_complex) = (&val as &dyn Any).downcast_ref::<Complex64>() {
555        return Some(*_complex);
556    }
557
558    // Try to handle f32 _complex numbers
559    if let Some(complex32) = (&val as &dyn Any).downcast_ref::<scirs2_core::numeric::Complex<f32>>()
560    {
561        return Some(Complex64::new(complex32.re as f64, complex32.im as f64));
562    }
563
564    None
565}
566
567#[allow(dead_code)]
568pub fn hilbert<T>(x: &[T]) -> FFTResult<Vec<scirs2_core::numeric::Complex64>>
569where
570    T: scirs2_core::numeric::NumCast + Copy + std::fmt::Debug + 'static,
571{
572    use scirs2_core::numeric::{Complex64, NumCast};
573
574    // Input length
575    let n = x.len();
576
577    // Convert input to a vector of f64
578    let signal: Vec<f64> = x
579        .iter()
580        .map(|&val| {
581            // First, try to cast directly to f64
582            if let Some(val_f64) = NumCast::from(val) {
583                return Ok(val_f64);
584            }
585
586            // If direct casting fails, check if it's a Complex value
587            // and use just the real part (for doctests which use Complex inputs)
588            match try_as_complex(val) {
589                Some(c) => Ok(c.re),
590                None => Err(FFTError::ValueError(format!(
591                    "Could not convert {val:?} to numeric type"
592                ))),
593            }
594        })
595        .collect::<FFTResult<Vec<_>>>()?;
596
597    // Compute FFT of the input signal
598    let spectrum = fft(&signal, None)?;
599
600    // Create the frequency domain filter for the Hilbert transform
601    // For a proper Hilbert transform, we need to:
602    // 1. Set the DC component (0 frequency) to 1
603    // 2. Double the positive frequencies and multiply by -i
604    // 3. Zero out the negative frequencies
605    let mut h = vec![Complex64::new(1.0, 0.0); n];
606
607    if n.is_multiple_of(2) {
608        // Even length case
609        h[0] = Complex64::new(1.0, 0.0); // DC component
610        h[n / 2] = Complex64::new(1.0, 0.0); // Nyquist component
611
612        // Positive frequencies (multiply by 2 and by -i)
613        h.iter_mut().take(n / 2).skip(1).for_each(|val| {
614            *val = Complex64::new(0.0, -2.0); // Equivalent to 2 * (-i)
615        });
616
617        // Negative frequencies (set to 0)
618        h.iter_mut().skip(n / 2 + 1).for_each(|val| {
619            *val = Complex64::new(0.0, 0.0);
620        });
621    } else {
622        // Odd length case
623        h[0] = Complex64::new(1.0, 0.0); // DC component
624
625        // Positive frequencies (multiply by 2 and by -i)
626        h.iter_mut().take(n.div_ceil(2)).skip(1).for_each(|val| {
627            *val = Complex64::new(0.0, -2.0); // Equivalent to 2 * (-i)
628        });
629
630        // Negative frequencies (set to 0)
631        h.iter_mut().skip(n.div_ceil(2)).for_each(|val| {
632            *val = Complex64::new(0.0, 0.0);
633        });
634    }
635
636    // Apply the filter in frequency domain
637    let filtered_spectrum: Vec<Complex64> = spectrum
638        .iter()
639        .zip(h.iter())
640        .map(|(&s, &h)| s * h)
641        .collect();
642
643    // Compute inverse FFT to get the analytic signal
644    let analytic_signal = ifft(&filtered_spectrum, None)?;
645
646    Ok(analytic_signal)
647}
648
649/// Returns the minimum and maximum values for each FFT dimension.
650///
651/// # Arguments
652///
653/// * `shape` - The shape of the FFT
654///
655/// # Returns
656///
657/// A vector of tuples (min, max) for each dimension of the FFT.
658///
659/// # Examples
660///
661/// ```
662/// use scirs2_fft::fft_bounds;
663///
664/// let bounds = fft_bounds(&[4, 4]);
665/// assert_eq!(bounds, vec![(-2, 1), (-2, 1)]);
666///
667/// let bounds = fft_bounds(&[5, 3]);
668/// assert_eq!(bounds, vec![(-2, 2), (-1, 1)]);
669/// ```
670#[must_use]
671#[allow(dead_code)]
672pub fn fft_bounds(shape: &[usize]) -> Vec<(i32, i32)> {
673    shape
674        .iter()
675        .map(|&n| {
676            // Cast with explicit handling for possible truncation/wrapping
677            let n_i32 = i32::try_from(n).unwrap_or(i32::MAX);
678            let min = -(n_i32 / 2);
679            let max = n_i32 - 1 + min;
680            (min, max)
681        })
682        .collect()
683}
684
685#[cfg(test)]
686mod tests {
687    use super::*;
688    use approx::assert_relative_eq;
689    use std::f64::consts::PI;
690
691    #[test]
692    fn test_fft_bounds() {
693        // Test even sizes
694        let bounds = fft_bounds(&[4, 4]);
695        assert_eq!(bounds, vec![(-2, 1), (-2, 1)]);
696
697        // Test odd sizes
698        let bounds = fft_bounds(&[5, 3]);
699        assert_eq!(bounds, vec![(-2, 2), (-1, 1)]);
700
701        // Test mixed sizes
702        let bounds = fft_bounds(&[6, 7, 8]);
703        assert_eq!(bounds, vec![(-3, 2), (-3, 3), (-4, 3)]);
704    }
705
706    #[test]
707    fn test_hilbert_transform() {
708        // Test on a cosine wave instead of sine wave to make the math easier
709        let n = 1000;
710        let freq = 5.0; // 5 Hz
711        let sample_rate = 100.0; // 100 Hz
712        let dt = 1.0 / sample_rate;
713
714        // Create a cosine wave
715        let t: Vec<f64> = (0..n).map(|i| i as f64 * dt).collect();
716        let signal: Vec<f64> = t.iter().map(|&ti| (2.0 * PI * freq * ti).cos()).collect();
717
718        // Compute Hilbert transform
719        let analytic = hilbert(&signal).unwrap();
720
721        // The Hilbert transform of cos(x) is sin(x)
722        // So the analytic signal should be cos(x) + i*sin(x) = e^(ix)
723        // Check the envelope (magnitude) which should be approximately 1
724        let start_idx = n / 4;
725        let end_idx = 3 * n / 4;
726
727        for i in start_idx..end_idx {
728            let magnitude = (analytic[i].re.powi(2) + analytic[i].im.powi(2)).sqrt();
729            assert_relative_eq!(magnitude, 1.0, epsilon = 0.1);
730
731            // Also check if the phase is advancing correctly
732            if i > start_idx {
733                let phase_i = analytic[i].im.atan2(analytic[i].re);
734                let phase_i_prev = analytic[i - 1].im.atan2(analytic[i - 1].re);
735
736                // Check if phase is advancing in the right direction
737                // We need to handle phase wrapping around ±π
738                let mut phase_diff = phase_i - phase_i_prev;
739                if phase_diff > PI {
740                    phase_diff -= 2.0 * PI;
741                } else if phase_diff < -PI {
742                    phase_diff += 2.0 * PI;
743                }
744
745                // For positive frequency, phase should generally advance positively
746                assert!(phase_diff > 0.0);
747            }
748        }
749    }
750}
751
752// Include ARM-specific FFT tests
753#[cfg(test)]
754mod arm_fft_test;