eFFT
An efficient method for the calculation of the exact Fourier transform of an asynchronous event stream
 
Loading...
Searching...
No Matches
efft.hpp
Go to the documentation of this file.
1
9#ifndef EFFT_HPP
10#define EFFT_HPP
11
12#include <algorithm>
13#include <array>
14#include <cmath>
15#include <complex>
16#include <cstddef>
17#include <eigen3/Eigen/Core>
18#include <ostream>
19#include <stdint.h>
20#include <unordered_map>
21#include <utility>
22#include <vector>
23
24#ifdef EFFT_USE_FFTW3
25#include <fftw3.h>
26#include <set>
27#endif
28
29class Stimulus {
30public:
31 unsigned int row{0};
32 unsigned int col{0};
33 bool state{true};
34 Stimulus() = default;
35 Stimulus(const unsigned int row, const unsigned int col) : row{row}, col{col} {};
36 Stimulus(const unsigned int row, const unsigned int col, const bool state) : row{row}, col{col}, state{state} {};
37 bool operator==(const Stimulus &p) const { return (row == p.row && col == p.col); }
38 bool operator!=(const Stimulus &p) const { return (row != p.row || col != p.col); }
39 Stimulus &on() {
40 state = true;
41 return *this;
42 }
43 Stimulus &off() {
44 state = false;
45 return *this;
46 }
47 Stimulus &set(const bool s) {
48 state = s;
49 return *this;
50 }
51 Stimulus &toggle() {
52 state = !state;
53 return *this;
54 }
55 friend std::ostream &operator<<(std::ostream &os, const Stimulus &stimulus) {
56 os << "Stimulus(row: " << stimulus.row << ", col: " << stimulus.col << ", state: " << (stimulus.state ? "on" : "off") << ")";
57 return os;
58 }
59};
60
61class Stimuli : public std::vector<Stimulus> {
62 using std::vector<Stimulus>::vector;
63
64 static inline uint64_t pack(int row, int col) {
65 return (static_cast<uint64_t>(static_cast<uint32_t>(row)) << 32) | static_cast<uint32_t>(col);
66 }
67
68public:
69 void on() {
70 std::for_each(begin(), end(), [](Stimulus &stimulus) { stimulus.state = true; });
71 }
72 void off() {
73 std::for_each(begin(), end(), [](Stimulus &stimulus) { stimulus.state = false; });
74 }
75 void set(const bool state) {
76 std::for_each(begin(), end(), [state](Stimulus &stimulus) { stimulus.state = state; });
77 }
78 void toggle() {
79 std::for_each(begin(), end(), [](Stimulus &stimulus) { stimulus.state = !stimulus.state; });
80 }
81
86 void filter() {
87 std::vector<Stimulus> out;
88 std::unordered_map<uint64_t, size_t> pos;
89 out.reserve(size());
90 pos.reserve(size() * 2);
91
92 for(const auto &s : *this) {
93 uint64_t key = pack(s.row, s.col);
94 auto [it, inserted] = pos.emplace(key, out.size());
95 if(inserted) {
96 out.emplace_back(s);
97 } else {
98 auto &chosen = out[it->second];
99 if(s.state && !chosen.state) {
100 chosen = s; // prefer state true over false
101 }
102 }
103 }
104 this->assign(out.begin(), out.end());
105 }
106};
107
108using cfloat = std::complex<float>;
109using cfloatmat = Eigen::Matrix<cfloat, Eigen::Dynamic, Eigen::Dynamic>;
110
111constexpr unsigned int LOG2(const unsigned int n) { return ((n < 2) ? 0 : 1 + LOG2(n >> 1U)); }
112
113#ifdef __has_builtin
114#if __has_builtin(__builtin_ffs)
115#define EFFT_HAS_BUILTIN_FFS
116#endif
117#endif
118#ifdef EFFT_HAS_BUILTIN_FFS
119inline unsigned int log2i(const unsigned int n) { return __builtin_ffs(static_cast<int>(n)) - 1; }
120#else
121inline unsigned int log2i(const unsigned int n) { return static_cast<unsigned int>(std::log2f(n)); }
122#endif
123
124template <unsigned int N>
125class eFFT {
126private:
127 static constexpr unsigned int LOG2_N = LOG2(N);
128 std::array<std::vector<cfloatmat>, LOG2_N + 1> tree_;
129 std::array<cfloat, static_cast<std::size_t>(N) * static_cast<std::size_t>(N + 1)> twiddle_;
130#ifdef EFFT_USE_FFTW3
131 std::array<fftw_complex, static_cast<std::size_t>(N) * static_cast<std::size_t>(N)> fftwInput_;
132 std::array<fftw_complex, static_cast<std::size_t>(N) * static_cast<std::size_t>(N)> fftwOutput_;
133 fftw_plan plan_{nullptr};
134#endif
135
136public:
137 eFFT() {
138 constexpr float MINUS_TWO_PI = -2 * M_PI;
139 for(unsigned int i = 0; i < N; i++) {
140 for(unsigned int n = 0; n <= N; n++) {
141 twiddle_[i + N * n] = static_cast<cfloat>(std::polar(1.0F, MINUS_TWO_PI * static_cast<float>(i) / static_cast<float>(n)));
142 }
143 }
144 }
145
146 ~eFFT() {
147#ifdef EFFT_USE_FFTW3
148 if(plan_ != nullptr) {
149 fftw_destroy_plan(plan_);
150 }
151#endif
152 }
153
154 eFFT(const eFFT &) = default;
155 eFFT(eFFT &&) noexcept = default;
156 eFFT &operator=(const eFFT &) = default;
157 eFFT &operator=(eFFT &&) noexcept = default;
158
163 [[nodiscard]] constexpr unsigned int framesize() const {
164 return N;
165 }
166
170 void initialize() {
171 cfloatmat f0(cfloatmat::Zero(N, N));
172 initialize(f0);
173 }
174
181 void initialize(cfloatmat &x, const int offset = 0) {
182 const unsigned int n = x.rows();
183 if(n == 1) {
184 tree_[0].emplace_back(x);
185 return;
186 }
187 const unsigned int ndiv2 = n >> 1U;
188 const unsigned int idx = log2i(ndiv2);
189
190 cfloatmat s00(x(Eigen::seq(Eigen::fix<0>, Eigen::last, Eigen::fix<2>), Eigen::seq(Eigen::fix<0>, Eigen::last, Eigen::fix<2>)));
191 cfloatmat s01(x(Eigen::seq(Eigen::fix<0>, Eigen::last, Eigen::fix<2>), Eigen::seq(Eigen::fix<1>, Eigen::last, Eigen::fix<2>)));
192 cfloatmat s10(x(Eigen::seq(Eigen::fix<1>, Eigen::last, Eigen::fix<2>), Eigen::seq(Eigen::fix<0>, Eigen::last, Eigen::fix<2>)));
193 cfloatmat s11(x(Eigen::seq(Eigen::fix<1>, Eigen::last, Eigen::fix<2>), Eigen::seq(Eigen::fix<1>, Eigen::last, Eigen::fix<2>)));
194 initialize(s00, 4 * offset);
195 initialize(s01, 4 * offset + 4);
196 initialize(s10, 4 * offset + 8);
197 initialize(s11, 4 * offset + 12);
198
199 const cfloatmat &x00 = tree_[idx][offset];
200 const cfloatmat &x01 = tree_[idx][offset + 1];
201 const cfloatmat &x10 = tree_[idx][offset + 2];
202 const cfloatmat &x11 = tree_[idx][offset + 3];
203 const unsigned int Nn = N * n;
204 for(unsigned int i = 0; i < ndiv2; i++) {
205 for(unsigned int j = 0; j < ndiv2; j++) {
206 const cfloat tu = twiddle_[j + Nn] * x01(i, j);
207 const cfloat td = twiddle_[(i + j) + Nn] * x11(i, j);
208 const cfloat ts = twiddle_[i + Nn] * x10(i, j);
209
210 const cfloat a = x00(i, j) + tu;
211 const cfloat b = x00(i, j) - tu;
212 const cfloat c = ts + td;
213 const cfloat d = ts - td;
214
215 x(i, j) = a + c;
216 x(i, j + ndiv2) = b + d;
217 x(i + ndiv2, j) = a - c;
218 x(i + ndiv2, j + ndiv2) = b - d;
219 }
220 }
221 tree_[log2i(n)].emplace_back(x);
222 }
223
230 bool update(const Stimulus &p) { return update(tree_[LOG2_N][0], p); }
231
239 bool update(cfloatmat &x, const Stimulus &p, const unsigned int offset = 0) {
240 const unsigned int n = x.rows();
241 if(n == 1) {
242 return std::exchange(x(0, 0), p.state).real() != static_cast<float>(p.state);
243 }
244 const unsigned int ndiv2 = n >> 1U;
245 const unsigned int nndiv2 = n * ndiv2;
246 const unsigned int idx = log2i(ndiv2);
247
248 bool changed = false;
249 if(static_cast<bool>(p.row & 1U)) {
250 if(static_cast<bool>(p.col & 1U)) {
251 changed = update(tree_[idx][offset + 3], {p.row >> 1U, p.col >> 1U, p.state}, 4 * offset + 12); // odd-odd
252 } else {
253 changed = update(tree_[idx][offset + 2], {p.row >> 1U, p.col >> 1U, p.state}, 4 * offset + 8); // odd-even
254 }
255 } else {
256 if(static_cast<bool>(p.col & 1U)) {
257 changed = update(tree_[idx][offset + 1], {p.row >> 1U, p.col >> 1U, p.state}, 4 * offset + 4); // even-odd
258 } else {
259 changed = update(tree_[idx][offset], {p.row >> 1U, p.col >> 1U, p.state}, 4 * offset); // even-even
260 }
261 }
262
263 if(changed) {
264 const cfloat *x00 = tree_[idx][offset].data();
265 const cfloat *x01 = tree_[idx][offset + 1].data();
266 const cfloat *x10 = tree_[idx][offset + 2].data();
267 const cfloat *x11 = tree_[idx][offset + 3].data();
268 cfloat *xp = x.data();
269 const unsigned int Nn = N * n;
270
271 for(unsigned int j = 0; j < ndiv2; j++) {
272 const unsigned int ndiv2j = ndiv2 * j, nj = n * j;
273 for(unsigned int i = 0; i < ndiv2; i++) {
274 const unsigned int k = i + ndiv2j, k1 = i + nj, k2 = k1 + ndiv2;
275
276 const cfloat tu = twiddle_[j + Nn] * x01[k];
277 const cfloat td = twiddle_[i + j + Nn] * x11[k];
278 const cfloat ts = twiddle_[i + Nn] * x10[k];
279
280 const cfloat x00_k = x00[k];
281 const cfloat a = x00_k + tu;
282 const cfloat b = x00_k - tu;
283 const cfloat c = ts + td;
284 const cfloat d = ts - td;
285
286 xp[k1] = a + c;
287 xp[k1 + nndiv2] = b + d;
288 xp[k2] = a - c;
289 xp[k2 + nndiv2] = b - d;
290 }
291 }
292 }
293 return changed;
294 }
295
302 bool update(Stimuli &pv) { return update(tree_[LOG2_N][0], pv.begin(), pv.end()); }
303
312 bool update(cfloatmat &x, Stimuli::iterator b0, Stimuli::iterator e0, const unsigned int offset = 0) {
313 const unsigned int n = x.rows();
314 if(n == 1) {
315 if(e0 - b0 == 1) {
316 return std::exchange(x(0, 0), b0->state).real() != static_cast<float>(b0->state);
317 }
318 const bool state = std::any_of(b0, e0, [](const Stimulus &p) { return p.state; });
319 return std::exchange(x(0, 0), state).real() != static_cast<float>(state);
320 }
321 const unsigned int ndiv2 = n >> 1U;
322 const unsigned int nndiv2 = n * ndiv2;
323 const unsigned int idx = log2i(ndiv2);
324
325 Stimuli::iterator e1, e2, e3;
326 e2 = std::partition(b0, e0, [](const Stimulus &p) { return p.row & 1U; });
327 e1 = std::partition(b0, e2, [](const Stimulus &p) { return p.col & 1U; });
328 e3 = std::partition(e2, e0, [](const Stimulus &p) { return p.col & 1U; });
329
330 auto transformStimuli = [](Stimuli::iterator begin, Stimuli::iterator end) {
331 for(auto it = begin; it != end; ++it) {
332 it->row >>= 1U;
333 it->col >>= 1U;
334 }
335 };
336 transformStimuli(b0, e1);
337 transformStimuli(e1, e2);
338 transformStimuli(e2, e3);
339 transformStimuli(e3, e0);
340
341 bool changed = false;
342 if(b0 != e1) {
343 changed = update(tree_[idx][offset + 3], b0, e1, 4 * offset + 12) || changed; // odd-odd
344 }
345 if(e1 != e2) {
346 changed = update(tree_[idx][offset + 2], e1, e2, 4 * offset + 8) || changed; // odd-even
347 }
348 if(e2 != e3) {
349 changed = update(tree_[idx][offset + 1], e2, e3, 4 * offset + 4) || changed; // even-odd
350 }
351 if(e3 != e0) {
352 changed = update(tree_[idx][offset], e3, e0, 4 * offset) || changed; // even-even
353 }
354
355 if(changed) {
356 const cfloat *x00 = tree_[idx][offset].data();
357 const cfloat *x01 = tree_[idx][offset + 1].data();
358 const cfloat *x10 = tree_[idx][offset + 2].data();
359 const cfloat *x11 = tree_[idx][offset + 3].data();
360 cfloat *xp = x.data();
361 const unsigned int Nn = N * n;
362
363 for(unsigned int j = 0; j < ndiv2; j++) {
364 const unsigned int ndiv2j = ndiv2 * j, nj = n * j;
365 for(unsigned int i = 0; i < ndiv2; i++) {
366 const unsigned int k = i + ndiv2j, k1 = i + nj, k2 = k1 + ndiv2;
367
368 const cfloat tu = twiddle_[j + Nn] * x01[k];
369 const cfloat td = twiddle_[i + j + Nn] * x11[k];
370 const cfloat ts = twiddle_[i + Nn] * x10[k];
371
372 const cfloat x00_k = x00[k];
373 const cfloat a = x00_k + tu;
374 const cfloat b = x00_k - tu;
375 const cfloat c = ts + td;
376 const cfloat d = ts - td;
377
378 xp[k1] = a + c;
379 xp[k1 + nndiv2] = b + d;
380 xp[k2] = a - c;
381 xp[k2 + nndiv2] = b - d;
382 }
383 }
384 }
385 return changed;
386 }
387
388#ifdef EFFT_USE_FFTW3
394 void initializeGroundTruth(const cfloatmat &image = cfloatmat::Zero(N, N)) {
395 plan_ = fftw_plan_dft_2d(N, N, fftwInput_.data(), fftwOutput_.data(), FFTW_FORWARD, FFTW_ESTIMATE | FFTW_UNALIGNED | FFTW_NO_SIMD | FFTW_PRESERVE_INPUT);
396 for(Eigen::Index i = 0; i < image.rows(); i++) {
397 for(Eigen::Index j = 0; j < image.cols(); j++) {
398 fftwInput_[N * i + j][0] = image(i, j).real();
399 fftwInput_[N * i + j][1] = image(i, j).imag();
400 }
401 }
402 fftw_execute(plan_);
403 }
404
410 void updateGroundTruth(const Stimulus &p) {
411 fftwInput_[N * p.row + p.col][0] = static_cast<double>(p.state);
412 fftwInput_[N * p.row + p.col][1] = 0;
413 fftw_execute(plan_);
414 }
415
421 void updateGroundTruth(const Stimuli &pv) {
422 std::set<std::pair<unsigned int, unsigned int>> activated;
423 for(const Stimulus &p : pv) {
424 if(p.state) {
425 activated.insert({p.row, p.col});
426 } else if(activated.find({p.row, p.col}) != activated.end()) {
427 continue;
428 }
429 fftwInput_[N * p.row + p.col][0] = static_cast<double>(p.state);
430 fftwInput_[N * p.row + p.col][1] = 0;
431 }
432 fftw_execute(plan_);
433 }
434#endif
435
441 [[nodiscard]] inline Eigen::Matrix<cfloat, N, N> getFFT() const {
442 return tree_[LOG2_N][0];
443 }
444
445#ifdef EFFT_USE_FFTW3
451 [[nodiscard]] inline Eigen::Matrix<cfloat, N, N> getGroundTruthFFT() const {
452 return Eigen::Map<const Eigen::Matrix<std::complex<double>, N, N, Eigen::RowMajor>>(reinterpret_cast<const std::complex<double> *>(fftwOutput_.data())).template cast<cfloat>();
453 }
454
460 [[nodiscard]] inline double check() const {
461 return (getFFT() - getGroundTruthFFT()).norm();
462 }
463#endif
464};
465
466#endif // EFFT_HPP
Definition efft.hpp:61
void filter()
Filters out stimuli.
Definition efft.hpp:86
Definition efft.hpp:29
void initialize()
Initializes the FFT computation with zero matrix.
Definition efft.hpp:170
bool update(const Stimulus &p)
Updates the FFT with a single stimulus.
Definition efft.hpp:230
bool update(Stimuli &pv)
Updates the FFT with multiple stimuli.
Definition efft.hpp:302
Eigen::Matrix< cfloat, N, N > getFFT() const
Get the FFT result as an Eigen matrix of complex floats.
Definition efft.hpp:441
bool update(cfloatmat &x, Stimuli::iterator b0, Stimuli::iterator e0, const unsigned int offset=0)
Updates the FFT with multiple stimuli in the specified matrix.
Definition efft.hpp:312
constexpr unsigned int framesize() const
Get the frame size of the FFT.
Definition efft.hpp:163
void initialize(cfloatmat &x, const int offset=0)
Initializes the FFT computation with the provided matrix.
Definition efft.hpp:181
bool update(cfloatmat &x, const Stimulus &p, const unsigned int offset=0)
Updates the FFT with a single stimulus in the specified matrix.
Definition efft.hpp:239