LASP 1.0
Library for Acoustic Signal Processing
Loading...
Searching...
No Matches
lasp_fft.cpp
Go to the documentation of this file.
1/* #define DEBUGTRACE_ENABLED */
2#include <memory>
3#include <cassert>
4#include "debugtrace.hpp"
5#include "lasp_config.h"
6#include "lasp_fft.h"
7using rte = std::runtime_error;
8
9#if LASP_FFT_BACKEND == Armadillo
10class Fft_impl {
11public:
13 Fft_impl(const us nfft) : nfft(nfft) { DEBUGTRACE_ENTER; }
14
15 vc fft(const vd &time) {
16 DEBUGTRACE_ENTER;
17 vc res = arma::fft(time);
18 return res.rows(0, nfft / 2);
19 }
20 vd ifft(const vc &f) {
21 DEBUGTRACE_ENTER;
22 vc f_full(nfft);
23 vc f_above_nyq = arma::conj(arma::reverse(f));
24 f_full.rows(0, nfft / 2) = f;
25 f_full.rows(nfft / 2, nfft - 1) = f_above_nyq.rows(0, nfft/2-1);
26 return arma::real(arma::ifft(f_full));
27 }
28};
29#elif LASP_FFT_BACKEND == FFTW
30#include <fftw3.h>
31
35class Fft_impl {
36public:
37 us nfft;
38 fftw_plan forward_plan = nullptr;
39 fftw_plan reverse_plan = nullptr;
40 fftw_complex *frequencyDomain = nullptr;
41 d *timeDomain = nullptr;
42
43 Fft_impl(const us nfft) : nfft(nfft) {
44
45 timeDomain = (d *)fftw_malloc(sizeof(d) * nfft);
46
47 frequencyDomain =
48 (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * (nfft / 2 + 1));
49
50 forward_plan =
51 fftw_plan_dft_r2c_1d(nfft, timeDomain, frequencyDomain, FFTW_MEASURE);
52
53 reverse_plan =
54 fftw_plan_dft_c2r_1d(nfft, frequencyDomain, timeDomain, FFTW_MEASURE);
55 if (!forward_plan || !reverse_plan || !timeDomain || !frequencyDomain) {
56 throw rte("Error allocating FFT");
57 }
58 };
59
60 ~Fft_impl() {
61 fftw_destroy_plan(forward_plan);
62 fftw_destroy_plan(reverse_plan);
63 fftw_free(frequencyDomain);
64 fftw_free(timeDomain);
65 }
66
67 vc fft(const vd &time) {
68 if (time.n_rows != nfft) {
69 throw rte("Invalid size of input vector, should be equal to nfft");
70 }
71
72 vc res(nfft / 2 + 1);
73 memcpy(timeDomain, time.memptr(), sizeof(d) * nfft);
74
75 fftw_execute(forward_plan);
76
77 memcpy(reinterpret_cast<void *>(res.memptr()), frequencyDomain,
78 sizeof(c) * (nfft / 2 + 1));
79
80 return res;
81 }
82 vd ifft(const vc &f) {
83 DEBUGTRACE_ENTER;
84 if (f.n_rows != nfft / 2 + 1) {
85 throw rte("Invalid size of input vector, should be equal to nfft/2+1");
86 }
87 memcpy(frequencyDomain,
88 reinterpret_cast<void *>(const_cast<c *>(f.memptr())),
89 (nfft / 2 + 1) * sizeof(c));
90
91 fftw_execute(reverse_plan);
92
93 vd res(nfft);
94 memcpy(res.memptr(), timeDomain, nfft * sizeof(d));
95
96 /* Scale by dividing by nfft. Checked with numpy implementation
97 * that this indeed needs to be done for FFTW. */
98 res *= 1.0 / nfft;
99
100 return res;
101 }
102};
103#else
104#error \
105 "Cannot compile lasp_ffc.c, no FFT backend specified. Should either be FFTPack, or FFTW"
106#endif
107
108Fft::Fft(const us nfft) {
109 DEBUGTRACE_ENTER;
110 if (nfft == 0) {
111 throw rte("Invalid nfft: 0");
112 }
113 if (nfft >= LASP_MAX_NFFT) {
114 throw rte("Invalid nfft, should be smaller than: " +
115 std::to_string(LASP_MAX_NFFT));
116 }
117 _impl = std::make_unique<Fft_impl>(nfft);
118}
120
121us Fft::nfft() const {
122 assert(_impl);
123 return _impl->nfft;
124}
125
126vc Fft::fft(const vd &timedata) {
127 assert(_impl);
128 return _impl->fft(timedata);
129}
130
131cmat Fft::fft(const dmat &freqdata) {
132 DEBUGTRACE_ENTER;
133 assert(_impl);
134 cmat res(_impl->nfft / 2 + 1, freqdata.n_cols);
138 // #pragma omp parallel for
139 for (us colno = 0; colno < freqdata.n_cols; colno++) {
140 res.col(colno) = _impl->fft(freqdata.col(colno));
141 }
142 return res;
143}
144
145vd Fft::ifft(const vc &freqdata) {
146 DEBUGTRACE_ENTER;
147 assert(_impl);
148 return _impl->ifft(freqdata);
149}
150dmat Fft::ifft(const cmat &freqdata) {
151 dmat res(_impl->nfft, freqdata.n_cols);
155 // #pragma omp parallel for
156
157 for (us colno = 0; colno < freqdata.n_cols; colno++) {
158 res.col(colno) = _impl->ifft(freqdata.col(colno));
159 }
160 return res;
161}
162
163void Fft::load_fft_wisdom(const std::string &wisdom) {
164#if LASP_FFT_BACKEND == Armadillo
165#elif LASP_FFT_BACKEND == FFTW
166 if (wisdom.length() > 0) {
167 int rv = fftw_import_wisdom_from_string(wisdom.c_str());
168 if (rv != 1) {
169 throw rte("Error loading FFTW wisdom");
170 }
171 }
172#endif
173}
174
176#if LASP_FFT_BACKEND == Armadillo
177 return "";
178#elif LASP_FFT_BACKEND == FFTW
179 // It is not possible to let FFTW directly return a C++ string. We have to
180 // put it in this container by copying in. Not a good solution if this has to
181 // happen often.
182 // Fortunately, this function is only called at the end of the program.
183 char *wis = fftw_export_wisdom_to_string();
184 std::string res{wis};
185 free(wis);
186 return res;
187#endif
188}
Fft_impl(const us nfft)
Definition lasp_fft.cpp:13
vd ifft(const vc &f)
Definition lasp_fft.cpp:20
vc fft(const vd &time)
Definition lasp_fft.cpp:15
~Fft()
Definition lasp_fft.cpp:119
vd ifft(const vc &freqdata)
Definition lasp_fft.cpp:145
static std::string store_fft_wisdom()
Return a string containing FFT wisdom storage. String is empty for backend != FFTW.
Definition lasp_fft.cpp:175
vc fft(const vd &timedata)
Definition lasp_fft.cpp:126
us nfft() const
Return nfft.
Definition lasp_fft.cpp:121
static void load_fft_wisdom(const std::string &wisdom)
Load FFT wisdom from a wisdom string. Function does nothing if FFT backend is not FFTW.
Definition lasp_fft.cpp:163
Fft(const us nfft)
Initialize FFT.
Definition lasp_fft.cpp:108
std::runtime_error rte
Definition lasp_daq.cpp:16
std::runtime_error rte
Definition lasp_fft.cpp:7
arma::Col< d > vd
arma::Mat< d > dmat
arma::Col< c > vc
arma::Mat< c > cmat
double complex c
Definition lasp_types.h:49
size_t us
We often use boolean values.
Definition lasp_types.h:29