3828
|
1 /* |
|
2 |
|
3 This file is part of Octave. |
|
4 |
|
5 Octave is free software; you can redistribute it and/or modify it |
|
6 under the terms of the GNU General Public License as published by the |
|
7 Free Software Foundation; either version 2, or (at your option) any |
|
8 later version. |
|
9 |
|
10 Octave is distributed in the hope that it will be useful, but WITHOUT |
|
11 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
|
12 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
|
13 for more details. |
|
14 |
|
15 You should have received a copy of the GNU General Public License |
|
16 along with Octave; see the file COPYING. If not, write to the Free |
|
17 Software Foundation, 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. |
|
18 |
|
19 */ |
|
20 |
|
21 #ifdef HAVE_CONFIG_H |
|
22 #include <config.h> |
|
23 #endif |
|
24 |
4773
|
25 #if defined (HAVE_FFTW3) |
3828
|
26 |
4775
|
27 #include <iostream> |
|
28 #include <vector> |
|
29 |
3828
|
30 #include "oct-fftw.h" |
|
31 #include "lo-error.h" |
|
32 |
|
33 // Helper class to create and cache fftw plans for both 1d and 2d. This |
|
34 // implementation uses FFTW_ESTIMATE to create the plans, which in theory |
4773
|
35 // is suboptimal, but provides quite reasonable performance. |
|
36 |
|
37 // Also note that if FFTW_ESTIMATE is not used the planner in FFTW3 |
|
38 // destroys the input and output arrays. So with the form of the |
|
39 // current code we definitely want FFTW_ESTIMATE!! However, we use |
|
40 // any wsidom that is available, either in a FFTW3 system wide file |
|
41 // or as supplied by the user. |
|
42 |
|
43 // XXX FIXME XXX If we can ensure 16 byte alignment in Array<T> (<T> *data) |
|
44 // the FFTW3 can use SIMD instructions for further acceleration. |
|
45 |
|
46 // Note that it is profitable to store the FFTW3 plans, for small ffts |
3828
|
47 |
|
48 class |
|
49 octave_fftw_planner |
|
50 { |
|
51 public: |
|
52 octave_fftw_planner (); |
|
53 |
4773
|
54 fftw_plan create_plan (int dir, const int rank, const dim_vector dims, |
|
55 int howmany, int stride, int dist, |
|
56 const Complex *in, Complex *out); |
|
57 fftw_plan create_plan (const int rank, const dim_vector dims, |
|
58 int howmany, int stride, int dist, |
|
59 const double *in, Complex *out); |
3828
|
60 |
|
61 private: |
|
62 int plan_flags; |
|
63 |
4773
|
64 // Plan for fft and ifft of complex values |
3828
|
65 fftw_plan plan[2]; |
4773
|
66 int d[2]; // dist |
|
67 int s[2]; // stride |
|
68 int r[2]; // rank |
|
69 int h[2]; // howmany |
|
70 dim_vector n[2]; // dims |
|
71 char ialign[2]; |
|
72 char oalign[2]; |
3828
|
73 |
4773
|
74 // Plan for fft of real values |
|
75 fftw_plan rplan; |
|
76 int rd; // dist |
|
77 int rs; // stride |
|
78 int rr; // rank |
|
79 int rh; // howmany |
|
80 dim_vector rn; // dims |
|
81 char rialign; |
|
82 char roalign; |
3828
|
83 }; |
|
84 |
|
85 octave_fftw_planner::octave_fftw_planner () |
|
86 { |
|
87 plan_flags = FFTW_ESTIMATE; |
|
88 |
|
89 plan[0] = plan[1] = 0; |
4773
|
90 d[0] = d[1] = s[0] = s[1] = r[0] = r[1] = h[0] = h[1] = 0; |
|
91 ialign[0] = ialign[1] = oalign[0] = oalign[1] = 0; |
|
92 n[0] = n[1] = dim_vector(); |
|
93 |
|
94 rplan = 0; |
|
95 rd = rs = rr = rh = 0; |
|
96 rialign = roalign = 0; |
|
97 rn = dim_vector (); |
|
98 |
|
99 // If we have a system wide wisdom file, import it |
|
100 fftw_import_system_wisdom ( ); |
3828
|
101 } |
|
102 |
|
103 fftw_plan |
4773
|
104 octave_fftw_planner::create_plan (int dir, const int rank, |
|
105 const dim_vector dims, int howmany, |
|
106 int stride, int dist, |
|
107 const Complex *in, Complex *out) |
3828
|
108 { |
4773
|
109 int which = (dir == FFTW_FORWARD) ? 0 : 1; |
3828
|
110 fftw_plan *cur_plan_p = &plan[which]; |
|
111 bool create_new_plan = false; |
4783
|
112 char in_align = (reinterpret_cast<long> (in)) & 0xF; |
|
113 char out_align = (reinterpret_cast<long> (out)) & 0xF; |
3828
|
114 |
4783
|
115 if (plan[which] == 0 || d[which] != dist || s[which] != stride |
|
116 || r[which] != rank || h[which] != howmany |
|
117 || ialign[which] != in_align || oalign[which] != out_align) |
4773
|
118 create_new_plan = true; |
|
119 else |
|
120 // We still might not have the same shape of array |
|
121 for (int i = 0; i < rank; i++) |
|
122 if (dims(i) != n[which](i)) |
|
123 { |
|
124 create_new_plan = true; |
|
125 break; |
|
126 } |
3828
|
127 |
|
128 if (create_new_plan) |
|
129 { |
4773
|
130 d[which] = dist; |
|
131 s[which] = stride; |
|
132 r[which] = rank; |
|
133 h[which] = howmany; |
|
134 ialign[which] = in_align; |
|
135 oalign[which] = out_align; |
|
136 n[which] = dims; |
|
137 |
3828
|
138 if (*cur_plan_p) |
|
139 fftw_destroy_plan (*cur_plan_p); |
|
140 |
4773
|
141 // Note reversal of dimensions for column major storage in FFTW |
|
142 OCTAVE_LOCAL_BUFFER (int, tmp, rank); |
|
143 for (int i = 0, j = rank-1; i < rank; i++, j--) |
|
144 tmp[i] = dims(j); |
|
145 |
|
146 *cur_plan_p = |
|
147 fftw_plan_many_dft (rank, tmp, howmany, |
|
148 reinterpret_cast<fftw_complex *> (const_cast<Complex *> (in)), |
4774
|
149 0, stride, dist, reinterpret_cast<fftw_complex *> (out), |
|
150 0, stride, dist, dir, plan_flags); |
3828
|
151 |
|
152 if (*cur_plan_p == 0) |
|
153 (*current_liboctave_error_handler) ("Error creating fftw plan"); |
|
154 } |
|
155 |
|
156 return *cur_plan_p; |
|
157 } |
|
158 |
4773
|
159 fftw_plan |
|
160 octave_fftw_planner::create_plan (const int rank, const dim_vector dims, |
|
161 int howmany, int stride, int dist, |
|
162 const double *in, Complex *out) |
3828
|
163 { |
4773
|
164 fftw_plan *cur_plan_p = &rplan; |
3828
|
165 bool create_new_plan = false; |
4783
|
166 char in_align = (reinterpret_cast<long> (in)) & 0xF; |
|
167 char out_align = (reinterpret_cast<long> (out)) & 0xF; |
3828
|
168 |
4783
|
169 if (rplan == 0 || rd != dist || rs != stride || rr != rank |
|
170 || rh != howmany || rialign != in_align || roalign != out_align) |
4773
|
171 create_new_plan = true; |
|
172 else |
|
173 // We still might not have the same shape of array |
|
174 for (int i = 0; i < rank; i++) |
|
175 if (dims(i) != rn(i)) |
|
176 { |
|
177 create_new_plan = true; |
|
178 break; |
|
179 } |
3828
|
180 |
|
181 if (create_new_plan) |
|
182 { |
4773
|
183 rd = dist; |
|
184 rs = stride; |
|
185 rr = rank; |
|
186 rh = howmany; |
|
187 rialign = in_align; |
|
188 roalign = out_align; |
|
189 rn = dims; |
|
190 |
3828
|
191 if (*cur_plan_p) |
4773
|
192 fftw_destroy_plan (*cur_plan_p); |
3828
|
193 |
4773
|
194 // Note reversal of dimensions for column major storage in FFTW |
|
195 OCTAVE_LOCAL_BUFFER (int, tmp, rank); |
|
196 for (int i = 0, j = rank-1; i < rank; i++, j--) |
|
197 tmp[i] = dims(j); |
|
198 |
|
199 *cur_plan_p = |
|
200 fftw_plan_many_dft_r2c (rank, tmp, howmany, |
|
201 (const_cast<double *> (in)), |
4774
|
202 0, stride, dist, reinterpret_cast<fftw_complex *> (out), |
|
203 0, stride, dist, plan_flags); |
3828
|
204 |
|
205 if (*cur_plan_p == 0) |
4773
|
206 (*current_liboctave_error_handler) ("Error creating fftw plan"); |
3828
|
207 } |
|
208 |
|
209 return *cur_plan_p; |
|
210 } |
|
211 |
|
212 static octave_fftw_planner fftw_planner; |
|
213 |
4775
|
214 static inline void |
|
215 convert_packcomplex_1d (Complex *out, size_t nr, size_t nc, |
|
216 int stride, int dist) |
4773
|
217 { |
|
218 // Fill in the missing data |
|
219 for (size_t i = 0; i < nr; i++) |
|
220 for (size_t j = nc/2+1; j < nc; j++) |
|
221 out[j*stride + i*dist] = conj(out[(nc - j)*stride + i*dist]); |
|
222 } |
|
223 |
4775
|
224 static inline void |
|
225 convert_packcomplex_Nd (Complex *out, const dim_vector &dv) |
3828
|
226 { |
4773
|
227 size_t nc = dv(0); |
|
228 size_t nr = dv(1); |
|
229 size_t np = (dv.length() > 2 ? dv.numel () / nc / nr : 1); |
|
230 size_t nrp = nr * np; |
|
231 Complex *ptr1, *ptr2; |
|
232 |
|
233 // Create space for the missing elements |
|
234 for (size_t i = 0; i < nrp; i++) |
|
235 { |
|
236 ptr1 = out + i * (nc/2 + 1) + nrp*((nc-1)/2); |
|
237 ptr2 = out + i * nc; |
|
238 for (size_t j = 0; j < nc/2+1; j++) |
|
239 *ptr2++ = *ptr1++; |
|
240 } |
|
241 |
|
242 // Fill in the missing data for the rank = 2 case directly for speed |
|
243 for (size_t i = 0; i < np; i++) |
|
244 { |
|
245 for (size_t j = 1; j < nr; j++) |
|
246 for (size_t k = nc/2+1; k < nc; k++) |
|
247 out[k + (j + i*nr)*nc] = conj(out[nc - k + ((i+1)*nr - j)*nc]); |
|
248 |
|
249 for (size_t j = nc/2+1; j < nc; j++) |
|
250 out[j + i*nr*nc] = conj(out[(i*nr+1)*nc - j]); |
|
251 } |
|
252 |
|
253 // Now do the permutations needed for rank > 2 cases |
|
254 size_t jstart = dv(0) * dv(1); |
|
255 size_t kstep = dv(0); |
|
256 size_t nel = dv.numel (); |
|
257 for (int inner = 2; inner < dv.length(); inner++) |
|
258 { |
|
259 size_t jmax = jstart * dv(inner); |
|
260 for (size_t i = 0; i < nel; i+=jmax) |
|
261 for (size_t j = jstart, jj = jmax-jstart; j < jj; |
|
262 j+=jstart, jj-=jstart) |
|
263 for (size_t k = 0; k < jstart; k+= kstep) |
|
264 for (size_t l = nc/2+1; l < nc; l++) |
|
265 { |
|
266 Complex tmp = out[i+ j + k + l]; |
|
267 out[i + j + k + l] = out[i + jj + k + l]; |
|
268 out[i + jj + k + l] = tmp; |
|
269 } |
|
270 jstart = jmax; |
|
271 } |
|
272 } |
|
273 |
|
274 int |
|
275 octave_fftw::fft (const double *in, Complex *out, size_t npts, |
|
276 size_t nsamples, int stride, int dist) |
|
277 { |
|
278 dist = (dist < 0 ? npts : dist); |
|
279 |
|
280 dim_vector dv (npts); |
|
281 fftw_plan plan = fftw_planner.create_plan (1, dv, nsamples, stride, dist, |
|
282 in, out); |
|
283 |
|
284 fftw_execute_dft_r2c (plan, (const_cast<double *>(in)), |
|
285 reinterpret_cast<fftw_complex *> (out)); |
|
286 |
|
287 // Need to create other half of the transform |
|
288 convert_packcomplex_1d (out, nsamples, npts, stride, dist); |
3828
|
289 |
|
290 return 0; |
|
291 } |
|
292 |
|
293 int |
4773
|
294 octave_fftw::fft (const Complex *in, Complex *out, size_t npts, |
|
295 size_t nsamples, int stride, int dist) |
3828
|
296 { |
4773
|
297 dist = (dist < 0 ? npts : dist); |
|
298 |
|
299 dim_vector dv (npts); |
|
300 fftw_plan plan = fftw_planner.create_plan (FFTW_FORWARD, 1, dv, nsamples, |
|
301 stride, dist, in, out); |
|
302 |
|
303 fftw_execute_dft (plan, |
|
304 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
|
305 reinterpret_cast<fftw_complex *> (out)); |
|
306 |
|
307 return 0; |
|
308 } |
|
309 |
|
310 int |
|
311 octave_fftw::ifft (const Complex *in, Complex *out, size_t npts, |
|
312 size_t nsamples, int stride, int dist) |
|
313 { |
|
314 dist = (dist < 0 ? npts : dist); |
|
315 |
|
316 dim_vector dv (npts); |
|
317 fftw_plan plan = fftw_planner.create_plan (FFTW_BACKWARD, 1, dv, nsamples, |
|
318 stride, dist, in, out); |
|
319 |
|
320 fftw_execute_dft (plan, |
|
321 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
|
322 reinterpret_cast<fftw_complex *> (out)); |
3828
|
323 |
|
324 const Complex scale = npts; |
4773
|
325 for (size_t j = 0; j < nsamples; j++) |
|
326 for (size_t i = 0; i < npts; i++) |
|
327 out[i*stride + j*dist] /= scale; |
3828
|
328 |
|
329 return 0; |
|
330 } |
|
331 |
|
332 int |
4773
|
333 octave_fftw::fftNd (const double *in, Complex *out, const int rank, |
|
334 const dim_vector &dv) |
3828
|
335 { |
4773
|
336 int dist = 1; |
|
337 for (int i = 0; i < rank; i++) |
|
338 dist *= dv(i); |
|
339 |
|
340 // Fool with the position of the start of the output matrix, so that |
|
341 // creating other half of the matrix won't cause cache problems |
|
342 int offset = (dv.numel () / dv(0)) * ((dv(0) - 1) / 2); |
|
343 |
|
344 fftw_plan plan = fftw_planner.create_plan (rank, dv, 1, 1, dist, |
|
345 in, out + offset); |
|
346 |
|
347 fftw_execute_dft_r2c (plan, (const_cast<double *>(in)), |
|
348 reinterpret_cast<fftw_complex *> (out+ offset)); |
|
349 |
|
350 // Need to create other half of the transform |
|
351 convert_packcomplex_Nd (out, dv); |
3828
|
352 |
|
353 return 0; |
|
354 } |
|
355 |
|
356 int |
4773
|
357 octave_fftw::fftNd (const Complex *in, Complex *out, const int rank, |
|
358 const dim_vector &dv) |
3828
|
359 { |
4773
|
360 int dist = 1; |
|
361 for (int i = 0; i < rank; i++) |
|
362 dist *= dv(i); |
|
363 |
|
364 fftw_plan plan = fftw_planner.create_plan (FFTW_FORWARD, rank, dv, 1, 1, |
|
365 dist, in, out); |
|
366 |
|
367 fftw_execute_dft (plan, |
|
368 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
|
369 reinterpret_cast<fftw_complex *> (out)); |
|
370 |
|
371 return 0; |
|
372 } |
3828
|
373 |
4773
|
374 int |
|
375 octave_fftw::ifftNd (const Complex *in, Complex *out, const int rank, |
|
376 const dim_vector &dv) |
|
377 { |
|
378 int dist = 1; |
|
379 for (int i = 0; i < rank; i++) |
|
380 dist *= dv(i); |
|
381 |
|
382 fftw_plan plan = fftw_planner.create_plan (FFTW_BACKWARD, rank, dv, 1, 1, |
|
383 dist, in, out); |
|
384 |
|
385 fftw_execute_dft (plan, |
|
386 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
|
387 reinterpret_cast<fftw_complex *> (out)); |
|
388 |
|
389 const size_t npts = dv.numel (); |
3828
|
390 const Complex scale = npts; |
|
391 for (size_t i = 0; i < npts; i++) |
4773
|
392 out[i] /= scale; |
3828
|
393 |
|
394 return 0; |
|
395 } |
|
396 |
|
397 #endif |
|
398 |
|
399 /* |
|
400 ;;; Local Variables: *** |
|
401 ;;; mode: C++ *** |
|
402 ;;; End: *** |
|
403 */ |
|
404 |