Mercurial > hg > qm-dsp
comparison dsp/transforms/FFT.cpp @ 358:7a225d665ed2
Merge from branch kissfft
author | Chris Cannam <c.cannam@qmul.ac.uk> |
---|---|
date | Wed, 16 Oct 2013 12:52:44 +0100 |
parents | 650bbacf8288 |
children | 8c86761a5533 |
comparison
equal
deleted
inserted
replaced
352:5ff562f5b521 | 358:7a225d665ed2 |
---|---|
9 | 9 |
10 #include "FFT.h" | 10 #include "FFT.h" |
11 | 11 |
12 #include "maths/MathUtilities.h" | 12 #include "maths/MathUtilities.h" |
13 | 13 |
14 #include "kiss_fft.h" | |
15 #include "kiss_fftr.h" | |
16 | |
14 #include <cmath> | 17 #include <cmath> |
15 | 18 |
16 #include <iostream> | 19 #include <iostream> |
17 | 20 |
21 #include <stdexcept> | |
22 | |
23 class FFT::D | |
24 { | |
25 public: | |
26 D(int n) : m_n(n) { | |
27 m_planf = kiss_fft_alloc(m_n, 0, NULL, NULL); | |
28 m_plani = kiss_fft_alloc(m_n, 1, NULL, NULL); | |
29 m_kin = new kiss_fft_cpx[m_n]; | |
30 m_kout = new kiss_fft_cpx[m_n]; | |
31 } | |
32 | |
33 ~D() { | |
34 kiss_fft_free(m_planf); | |
35 kiss_fft_free(m_plani); | |
36 delete[] m_kin; | |
37 delete[] m_kout; | |
38 } | |
39 | |
40 void process(bool inverse, | |
41 const double *ri, | |
42 const double *ii, | |
43 double *ro, | |
44 double *io) { | |
45 | |
46 for (int i = 0; i < m_n; ++i) { | |
47 m_kin[i].r = ri[i]; | |
48 m_kin[i].i = (ii ? ii[i] : 0.0); | |
49 } | |
50 | |
51 if (!inverse) { | |
52 | |
53 kiss_fft(m_planf, m_kin, m_kout); | |
54 | |
55 for (int i = 0; i < m_n; ++i) { | |
56 ro[i] = m_kout[i].r; | |
57 io[i] = m_kout[i].i; | |
58 } | |
59 | |
60 } else { | |
61 | |
62 kiss_fft(m_plani, m_kin, m_kout); | |
63 | |
64 double scale = 1.0 / m_n; | |
65 | |
66 for (int i = 0; i < m_n; ++i) { | |
67 ro[i] = m_kout[i].r * scale; | |
68 io[i] = m_kout[i].i * scale; | |
69 } | |
70 } | |
71 } | |
72 | |
73 private: | |
74 int m_n; | |
75 kiss_fft_cfg m_planf; | |
76 kiss_fft_cfg m_plani; | |
77 kiss_fft_cpx *m_kin; | |
78 kiss_fft_cpx *m_kout; | |
79 }; | |
80 | |
18 FFT::FFT(int n) : | 81 FFT::FFT(int n) : |
19 m_n(n) | 82 m_d(new D(n)) |
20 { | 83 { |
21 if( !MathUtilities::isPowerOfTwo(m_n) ) | |
22 { | |
23 std::cerr << "ERROR: FFT: Non-power-of-two FFT size " | |
24 << m_n << " not supported in this implementation" | |
25 << std::endl; | |
26 return; | |
27 } | |
28 } | 84 } |
29 | 85 |
30 FFT::~FFT() | 86 FFT::~FFT() |
31 { | 87 { |
32 | 88 delete m_d; |
33 } | 89 } |
34 | 90 |
35 FFTReal::FFTReal(int n) : | 91 void |
36 m_n(n), | 92 FFT::process(bool inverse, |
37 m_fft(new FFT(n)), | |
38 m_r(new double[n]), | |
39 m_i(new double[n]), | |
40 m_discard(new double[n]) | |
41 { | |
42 for (int i = 0; i < n; ++i) { | |
43 m_r[i] = 0; | |
44 m_i[i] = 0; | |
45 m_discard[i] = 0; | |
46 } | |
47 } | |
48 | |
49 FFTReal::~FFTReal() | |
50 { | |
51 delete m_fft; | |
52 delete[] m_discard; | |
53 delete[] m_r; | |
54 delete[] m_i; | |
55 } | |
56 | |
57 void | |
58 FFTReal::forward(const double *realIn, | |
59 double *realOut, double *imagOut) | |
60 { | |
61 m_fft->process(false, realIn, 0, realOut, imagOut); | |
62 } | |
63 | |
64 void | |
65 FFTReal::inverse(const double *realIn, const double *imagIn, | |
66 double *realOut) | |
67 { | |
68 for (int i = 0; i < m_n/2 + 1; ++i) { | |
69 m_r[i] = realIn[i]; | |
70 m_i[i] = imagIn[i]; | |
71 if (i > 0 && i < m_n/2) { | |
72 m_r[m_n - i] = realIn[i]; | |
73 m_i[m_n - i] = -imagIn[i]; | |
74 } | |
75 } | |
76 m_fft->process(true, m_r, m_i, realOut, m_discard); | |
77 } | |
78 | |
79 static int numberOfBitsNeeded(int p_nSamples) | |
80 { | |
81 int i; | |
82 | |
83 if( p_nSamples < 2 ) | |
84 { | |
85 return 0; | |
86 } | |
87 | |
88 for ( i=0; ; i++ ) | |
89 { | |
90 if( p_nSamples & (1 << i) ) return i; | |
91 } | |
92 } | |
93 | |
94 static int reverseBits(int p_nIndex, int p_nBits) | |
95 { | |
96 int i, rev; | |
97 | |
98 for(i=rev=0; i < p_nBits; i++) | |
99 { | |
100 rev = (rev << 1) | (p_nIndex & 1); | |
101 p_nIndex >>= 1; | |
102 } | |
103 | |
104 return rev; | |
105 } | |
106 | |
107 void | |
108 FFT::process(bool p_bInverseTransform, | |
109 const double *p_lpRealIn, const double *p_lpImagIn, | 93 const double *p_lpRealIn, const double *p_lpImagIn, |
110 double *p_lpRealOut, double *p_lpImagOut) | 94 double *p_lpRealOut, double *p_lpImagOut) |
111 { | 95 { |
112 if (!p_lpRealIn || !p_lpRealOut || !p_lpImagOut) return; | 96 m_d->process(inverse, |
113 | 97 p_lpRealIn, p_lpImagIn, |
114 // std::cerr << "FFT::process(" << m_n << "," << p_bInverseTransform << ")" << std::endl; | 98 p_lpRealOut, p_lpImagOut); |
115 | 99 } |
116 int NumBits; | 100 |
117 int i, j, k, n; | 101 class FFTReal::D |
118 int BlockSize, BlockEnd; | 102 { |
119 | 103 public: |
120 double angle_numerator = 2.0 * M_PI; | 104 D(int n) : m_n(n) { |
121 double tr, ti; | 105 if (n % 2) { |
122 | 106 throw std::invalid_argument |
123 if( !MathUtilities::isPowerOfTwo(m_n) ) | 107 ("nsamples must be even in FFTReal constructor"); |
124 { | 108 } |
125 std::cerr << "ERROR: FFT::process: Non-power-of-two FFT size " | 109 m_planf = kiss_fftr_alloc(m_n, 0, NULL, NULL); |
126 << m_n << " not supported in this implementation" | 110 m_plani = kiss_fftr_alloc(m_n, 1, NULL, NULL); |
127 << std::endl; | 111 m_c = new kiss_fft_cpx[m_n]; |
128 return; | 112 } |
129 } | 113 |
130 | 114 ~D() { |
131 if( p_bInverseTransform ) angle_numerator = -angle_numerator; | 115 kiss_fftr_free(m_planf); |
132 | 116 kiss_fftr_free(m_plani); |
133 NumBits = numberOfBitsNeeded ( m_n ); | 117 delete[] m_c; |
134 | 118 } |
135 | 119 |
136 | 120 void forward(const double *ri, double *ro, double *io) { |
137 for( i=0; i < m_n; i++ ) | 121 |
138 { | 122 kiss_fftr(m_planf, ri, m_c); |
139 j = reverseBits ( i, NumBits ); | 123 |
140 p_lpRealOut[j] = p_lpRealIn[i]; | 124 for (int i = 0; i <= m_n/2; ++i) { |
141 p_lpImagOut[j] = (p_lpImagIn == 0) ? 0.0 : p_lpImagIn[i]; | 125 ro[i] = m_c[i].r; |
142 } | 126 io[i] = m_c[i].i; |
143 | 127 } |
144 | 128 |
145 BlockEnd = 1; | 129 for (int i = 0; i + 1 < m_n/2; ++i) { |
146 for( BlockSize = 2; BlockSize <= m_n; BlockSize <<= 1 ) | 130 ro[m_n - i - 1] = ro[i + 1]; |
147 { | 131 io[m_n - i - 1] = -io[i + 1]; |
148 double delta_angle = angle_numerator / (double)BlockSize; | 132 } |
149 double sm2 = -sin ( -2 * delta_angle ); | 133 } |
150 double sm1 = -sin ( -delta_angle ); | 134 |
151 double cm2 = cos ( -2 * delta_angle ); | 135 void forwardMagnitude(const double *ri, double *mo) { |
152 double cm1 = cos ( -delta_angle ); | 136 |
153 double w = 2 * cm1; | 137 double *io = new double[m_n]; |
154 double ar[3], ai[3]; | 138 |
155 | 139 forward(ri, mo, io); |
156 for( i=0; i < m_n; i += BlockSize ) | 140 |
157 { | 141 for (int i = 0; i < m_n; ++i) { |
158 | 142 mo[i] = sqrt(mo[i] * mo[i] + io[i] * io[i]); |
159 ar[2] = cm2; | 143 } |
160 ar[1] = cm1; | 144 |
161 | 145 delete[] io; |
162 ai[2] = sm2; | 146 } |
163 ai[1] = sm1; | 147 |
164 | 148 void inverse(const double *ri, const double *ii, double *ro) { |
165 for ( j=i, n=0; n < BlockEnd; j++, n++ ) | 149 |
166 { | 150 for (int i = 0; i < m_n; ++i) { |
167 | 151 m_c[i].r = ri[i]; |
168 ar[0] = w*ar[1] - ar[2]; | 152 m_c[i].i = ii[i]; |
169 ar[2] = ar[1]; | 153 } |
170 ar[1] = ar[0]; | 154 |
171 | 155 kiss_fftri(m_plani, m_c, ro); |
172 ai[0] = w*ai[1] - ai[2]; | 156 |
173 ai[2] = ai[1]; | 157 double scale = 1.0 / m_n; |
174 ai[1] = ai[0]; | 158 |
175 | 159 for (int i = 0; i < m_n; ++i) { |
176 k = j + BlockEnd; | 160 ro[i] *= scale; |
177 tr = ar[0]*p_lpRealOut[k] - ai[0]*p_lpImagOut[k]; | 161 } |
178 ti = ar[0]*p_lpImagOut[k] + ai[0]*p_lpRealOut[k]; | 162 } |
179 | 163 |
180 p_lpRealOut[k] = p_lpRealOut[j] - tr; | 164 private: |
181 p_lpImagOut[k] = p_lpImagOut[j] - ti; | 165 int m_n; |
182 | 166 kiss_fftr_cfg m_planf; |
183 p_lpRealOut[j] += tr; | 167 kiss_fftr_cfg m_plani; |
184 p_lpImagOut[j] += ti; | 168 kiss_fft_cpx *m_c; |
185 | 169 }; |
186 } | 170 |
187 } | 171 FFTReal::FFTReal(int n) : |
188 | 172 m_d(new D(n)) |
189 BlockEnd = BlockSize; | 173 { |
190 | 174 } |
191 } | 175 |
192 | 176 FFTReal::~FFTReal() |
193 | 177 { |
194 if( p_bInverseTransform ) | 178 delete m_d; |
195 { | 179 } |
196 double denom = (double)m_n; | 180 |
197 | 181 void |
198 for ( i=0; i < m_n; i++ ) | 182 FFTReal::forward(const double *ri, double *ro, double *io) |
199 { | 183 { |
200 p_lpRealOut[i] /= denom; | 184 m_d->forward(ri, ro, io); |
201 p_lpImagOut[i] /= denom; | 185 } |
202 } | 186 |
203 } | 187 void |
204 } | 188 FFTReal::forwardMagnitude(const double *ri, double *mo) |
205 | 189 { |
190 m_d->forwardMagnitude(ri, mo); | |
191 } | |
192 | |
193 void | |
194 FFTReal::inverse(const double *ri, const double *ii, double *ro) | |
195 { | |
196 m_d->inverse(ri, ii, ro); | |
197 } | |
198 | |
199 | |
200 |