numerics
Loading...
Searching...
No Matches
fft.cpp
Go to the documentation of this file.
1#include "spectral/fft.hpp"
2#include "backends/seq/impl.hpp"
5#ifdef NUMERICS_HAS_STD_SIMD
7#endif
8#include <stdexcept>
9
10namespace num {
11namespace spectral {
12
13// -- One-shot dispatch --------------------------------------------------------
14
15void fft(const CVector& in, CVector& out, FFTBackend b) {
16 if (out.size() != in.size())
17 throw std::invalid_argument("fft: in and out must have the same size");
18#ifdef NUMERICS_HAS_FFTW
19 if (b == FFTBackend::fftw) { backends::fftw::fft(in, out); return; }
20#endif
21#if defined(NUMERICS_HAS_AVX2) || defined(NUMERICS_HAS_NEON)
22 if (b == FFTBackend::simd) { backends::opt::fft(in, out); return; }
23#endif
24#ifdef NUMERICS_HAS_STD_SIMD
25 if (b == FFTBackend::stdsimd) { backends::stdsimd::fft(in, out); return; }
26#endif
27 // seq is the fallback for simd/stdsimd on unsupported platforms
29}
30
32 if (out.size() != in.size())
33 throw std::invalid_argument("ifft: in and out must have the same size");
34#ifdef NUMERICS_HAS_FFTW
35 if (b == FFTBackend::fftw) { backends::fftw::ifft(in, out); return; }
36#endif
37#if defined(NUMERICS_HAS_AVX2) || defined(NUMERICS_HAS_NEON)
38 if (b == FFTBackend::simd) { backends::opt::ifft(in, out); return; }
39#endif
40#ifdef NUMERICS_HAS_STD_SIMD
41 if (b == FFTBackend::stdsimd) { backends::stdsimd::ifft(in, out); return; }
42#endif
44}
45
46void rfft(const Vector& in, CVector& out, FFTBackend b) {
47 if (static_cast<int>(out.size()) != static_cast<int>(in.size()) / 2 + 1)
48 throw std::invalid_argument("rfft: out must have size n/2+1");
49#ifdef NUMERICS_HAS_FFTW
50 if (b == FFTBackend::fftw) { backends::fftw::rfft(in, out); return; }
51#endif
52#if defined(NUMERICS_HAS_AVX2) || defined(NUMERICS_HAS_NEON)
53 if (b == FFTBackend::simd) { backends::opt::rfft(in, out); return; }
54#endif
55#ifdef NUMERICS_HAS_STD_SIMD
56 if (b == FFTBackend::stdsimd) { backends::stdsimd::rfft(in, out); return; }
57#endif
59}
60
61void irfft(const CVector& in, int n, Vector& out, FFTBackend b) {
62 if (static_cast<int>(in.size()) != n / 2 + 1)
63 throw std::invalid_argument("irfft: in must have size n/2+1");
64 if (static_cast<int>(out.size()) != n)
65 throw std::invalid_argument("irfft: out must have size n");
66#ifdef NUMERICS_HAS_FFTW
67 if (b == FFTBackend::fftw) { backends::fftw::irfft(in, n, out); return; }
68#endif
69#if defined(NUMERICS_HAS_AVX2) || defined(NUMERICS_HAS_NEON)
70 if (b == FFTBackend::simd) { backends::opt::irfft(in, n, out); return; }
71#endif
72#ifdef NUMERICS_HAS_STD_SIMD
73 if (b == FFTBackend::stdsimd) { backends::stdsimd::irfft(in, n, out); return; }
74#endif
76}
77
78// -- FFTPlan ------------------------------------------------------------------
79
80FFTPlan::FFTPlan(int n, bool forward, FFTBackend b) : n_(n), backend_(b) {
81#ifdef NUMERICS_HAS_FFTW
82 if (b == FFTBackend::fftw) {
83 impl_ = new backends::fftw::FFTPlanImpl(n, forward);
84 return;
85 }
86#endif
87#if defined(NUMERICS_HAS_AVX2) || defined(NUMERICS_HAS_NEON)
88 if (b == FFTBackend::simd) {
89 impl_ = new backends::opt::FFTPlanImpl(n, !forward);
90 return;
91 }
92#endif
93#ifdef NUMERICS_HAS_STD_SIMD
94 if (b == FFTBackend::stdsimd) {
95 impl_ = new backends::stdsimd::FFTPlanImpl(n, !forward);
96 return;
97 }
98#endif
99 impl_ = new backends::seq::FFTPlanImpl(n, !forward);
100}
101
103#ifdef NUMERICS_HAS_FFTW
104 if (backend_ == FFTBackend::fftw) {
105 delete static_cast<backends::fftw::FFTPlanImpl*>(impl_);
106 return;
107 }
108#endif
109#if defined(NUMERICS_HAS_AVX2) || defined(NUMERICS_HAS_NEON)
110 if (backend_ == FFTBackend::simd) {
111 delete static_cast<backends::opt::FFTPlanImpl*>(impl_);
112 return;
113 }
114#endif
115#ifdef NUMERICS_HAS_STD_SIMD
116 if (backend_ == FFTBackend::stdsimd) {
117 delete static_cast<backends::stdsimd::FFTPlanImpl*>(impl_);
118 return;
119 }
120#endif
121 delete static_cast<backends::seq::FFTPlanImpl*>(impl_);
122}
123
124void FFTPlan::execute(const CVector& in, CVector& out) const {
125 if (static_cast<int>(in.size()) != n_ || static_cast<int>(out.size()) != n_)
126 throw std::invalid_argument("FFTPlan::execute: size mismatch");
127#ifdef NUMERICS_HAS_FFTW
128 if (backend_ == FFTBackend::fftw) {
129 static_cast<backends::fftw::FFTPlanImpl*>(impl_)->execute(in, out);
130 return;
131 }
132#endif
133#if defined(NUMERICS_HAS_AVX2) || defined(NUMERICS_HAS_NEON)
134 if (backend_ == FFTBackend::simd) {
135 for (idx i = 0; i < static_cast<idx>(n_); ++i) out[i] = in[i];
136 static_cast<backends::opt::FFTPlanImpl*>(impl_)->execute(out);
137 return;
138 }
139#endif
140#ifdef NUMERICS_HAS_STD_SIMD
141 if (backend_ == FFTBackend::stdsimd) {
142 for (idx i = 0; i < static_cast<idx>(n_); ++i) out[i] = in[i];
143 static_cast<backends::stdsimd::FFTPlanImpl*>(impl_)->execute(out);
144 return;
145 }
146#endif
147 for (idx i = 0; i < static_cast<idx>(n_); ++i) out[i] = in[i];
148 static_cast<backends::seq::FFTPlanImpl*>(impl_)->execute(out);
149}
150
151} // namespace spectral
152} // namespace num
Dense vector with optional GPU storage, templated over scalar type T.
Definition vector.hpp:23
FFTPlan(int n, bool forward=true, FFTBackend b=default_fft_backend)
Definition fft.cpp:80
void execute(const CVector &in, CVector &out) const
Definition fft.cpp:124
FFT interface with backend dispatch.
void irfft(const num::CVector &in, int n, num::Vector &out)
Definition impl.hpp:198
void ifft(const num::CVector &in, num::CVector &out)
Definition impl.hpp:182
void fft(const num::CVector &in, num::CVector &out)
Definition impl.hpp:175
void rfft(const num::Vector &in, num::CVector &out)
Definition impl.hpp:189
void ifft(const num::CVector &in, num::CVector &out)
Definition impl.hpp:51
void rfft(const num::Vector &in, num::CVector &out)
Definition impl.hpp:56
void irfft(const num::CVector &in, int n, num::Vector &out)
Definition impl.hpp:64
void fft(const num::CVector &in, num::CVector &out)
Definition impl.hpp:46
void ifft(const CVector &in, CVector &out, FFTBackend b=default_fft_backend)
Inverse complex DFT (unnormalised: result = n * true_inverse).
Definition fft.cpp:31
void fft(const CVector &in, CVector &out, FFTBackend b=default_fft_backend)
Forward complex DFT. out must be pre-allocated to in.size().
Definition fft.cpp:15
FFTBackend
Selects which backend handles an FFT operation.
Definition fft.hpp:26
@ simd
Handwritten AVX2 / NEON butterfly – falls back to seq if unavailable.
@ stdsimd
std::experimental::simd butterfly – requires NUMERICS_HAS_STD_SIMD
@ fftw
FFTW3 (mixed-radix, AVX / NEON optimised) – requires NUMERICS_HAS_FFTW.
void irfft(const CVector &in, int n, Vector &out, FFTBackend b=default_fft_backend)
Complex-to-real inverse DFT (unnormalised).
Definition fft.cpp:61
void rfft(const Vector &in, CVector &out, FFTBackend b=default_fft_backend)
Real-to-complex forward DFT. out must be pre-allocated to n/2+1.
Definition fft.cpp:46
constexpr T ipow(T x) noexcept
Compute x^N at compile time via repeated squaring.
std::size_t idx
Definition types.hpp:11
FFTW3 backend wrappers. Only included by src/spectral/fft.cpp.
Handwritten AVX2 / NEON butterfly for the FFT.
std::experimental::simd butterfly for FFT.
Precomputed twiddle factors + SIMD butterfly execution.
Definition impl.hpp:142