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 |
4786
|
30 #include "lo-error.h" |
3828
|
31 #include "oct-fftw.h" |
4786
|
32 #include "quit.h" |
3828
|
33 |
4809
|
34 // Helper class to create and cache fftw plans for both 1d and |
|
35 // 2d. This implementation uses FFTW_ESTIMATE to create the plans, |
|
36 // which in theory is suboptimal, but provides quite reasonable |
|
37 // performance. |
4773
|
38 |
|
39 // Also note that if FFTW_ESTIMATE is not used the planner in FFTW3 |
4809
|
40 // destroys the input and output arrays. So with the form of the |
|
41 // current code we definitely want FFTW_ESTIMATE!! However, we use any |
|
42 // wsidom that is available, either in a FFTW3 system wide file or as |
|
43 // supplied by the user. |
4773
|
44 |
4809
|
45 // XXX FIXME XXX -- if we can ensure 16 byte alignment in Array<T> |
|
46 // (<T> *data) the FFTW3 can use SIMD instructions for further |
|
47 // acceleration. |
4773
|
48 |
4809
|
49 // Note that it is profitable to store the FFTW3 plans, for small |
|
50 // ffts. |
3828
|
51 |
|
52 class |
|
53 octave_fftw_planner |
|
54 { |
|
55 public: |
4808
|
56 |
|
57 octave_fftw_planner (void); |
3828
|
58 |
4773
|
59 fftw_plan create_plan (int dir, const int rank, const dim_vector dims, |
|
60 int howmany, int stride, int dist, |
|
61 const Complex *in, Complex *out); |
4808
|
62 |
4773
|
63 fftw_plan create_plan (const int rank, const dim_vector dims, |
|
64 int howmany, int stride, int dist, |
|
65 const double *in, Complex *out); |
3828
|
66 |
|
67 private: |
4808
|
68 |
3828
|
69 int plan_flags; |
|
70 |
4809
|
71 // XXX FIXME XXX -- perhaps this should be split into two classes? |
|
72 |
4773
|
73 // Plan for fft and ifft of complex values |
3828
|
74 fftw_plan plan[2]; |
4809
|
75 |
|
76 // dist |
|
77 int d[2]; |
|
78 |
|
79 // stride |
|
80 int s[2]; |
|
81 |
|
82 // rank |
|
83 int r[2]; |
|
84 |
|
85 // howmany |
|
86 int h[2]; |
|
87 |
|
88 // dims |
|
89 dim_vector n[2]; |
|
90 |
4808
|
91 bool simd_align[2]; |
3828
|
92 |
4809
|
93 |
4773
|
94 // Plan for fft of real values |
|
95 fftw_plan rplan; |
4809
|
96 |
|
97 // dist |
|
98 int rd; |
|
99 |
|
100 // stride |
|
101 int rs; |
|
102 |
|
103 // rank |
|
104 int rr; |
|
105 |
|
106 // howmany |
|
107 int rh; |
|
108 |
|
109 // dims |
|
110 dim_vector rn; |
|
111 |
4808
|
112 bool rsimd_align; |
3828
|
113 }; |
|
114 |
4808
|
115 octave_fftw_planner::octave_fftw_planner (void) |
3828
|
116 { |
|
117 plan_flags = FFTW_ESTIMATE; |
|
118 |
|
119 plan[0] = plan[1] = 0; |
4773
|
120 d[0] = d[1] = s[0] = s[1] = r[0] = r[1] = h[0] = h[1] = 0; |
4808
|
121 simd_align[0] = simd_align[1] = false; |
|
122 n[0] = n[1] = dim_vector (); |
4773
|
123 |
|
124 rplan = 0; |
|
125 rd = rs = rr = rh = 0; |
4808
|
126 rsimd_align = false; |
4773
|
127 rn = dim_vector (); |
|
128 |
4809
|
129 // If we have a system wide wisdom file, import it. |
4808
|
130 fftw_import_system_wisdom (); |
3828
|
131 } |
|
132 |
4808
|
133 #define CHECK_SIMD_ALIGNMENT(x) \ |
|
134 ((reinterpret_cast<ptrdiff_t> (x)) & 0xF == 0) |
|
135 |
3828
|
136 fftw_plan |
4773
|
137 octave_fftw_planner::create_plan (int dir, const int rank, |
|
138 const dim_vector dims, int howmany, |
|
139 int stride, int dist, |
|
140 const Complex *in, Complex *out) |
3828
|
141 { |
4773
|
142 int which = (dir == FFTW_FORWARD) ? 0 : 1; |
3828
|
143 fftw_plan *cur_plan_p = &plan[which]; |
|
144 bool create_new_plan = false; |
4808
|
145 bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out); |
3828
|
146 |
4809
|
147 // Don't create a new plan if we have a non SIMD plan already but |
|
148 // can do SIMD. This prevents endlessly recreating plans if we |
|
149 // change the alignment. |
|
150 |
4783
|
151 if (plan[which] == 0 || d[which] != dist || s[which] != stride |
|
152 || r[which] != rank || h[which] != howmany |
4808
|
153 || ((ioalign != simd_align[which]) ? !ioalign : false)) |
4773
|
154 create_new_plan = true; |
|
155 else |
4809
|
156 { |
|
157 // We still might not have the same shape of array. |
|
158 |
|
159 for (int i = 0; i < rank; i++) |
|
160 if (dims(i) != n[which](i)) |
|
161 { |
|
162 create_new_plan = true; |
|
163 break; |
|
164 } |
|
165 } |
3828
|
166 |
|
167 if (create_new_plan) |
|
168 { |
4773
|
169 d[which] = dist; |
|
170 s[which] = stride; |
|
171 r[which] = rank; |
|
172 h[which] = howmany; |
4808
|
173 simd_align[which] = ioalign; |
4773
|
174 n[which] = dims; |
|
175 |
4808
|
176 if (ioalign) |
|
177 plan_flags &= ~FFTW_UNALIGNED; |
|
178 else |
|
179 plan_flags |= FFTW_UNALIGNED; |
|
180 |
3828
|
181 if (*cur_plan_p) |
|
182 fftw_destroy_plan (*cur_plan_p); |
|
183 |
4809
|
184 // Note reversal of dimensions for column major storage in FFTW. |
|
185 |
4773
|
186 OCTAVE_LOCAL_BUFFER (int, tmp, rank); |
4809
|
187 |
4773
|
188 for (int i = 0, j = rank-1; i < rank; i++, j--) |
|
189 tmp[i] = dims(j); |
|
190 |
|
191 *cur_plan_p = |
|
192 fftw_plan_many_dft (rank, tmp, howmany, |
|
193 reinterpret_cast<fftw_complex *> (const_cast<Complex *> (in)), |
4774
|
194 0, stride, dist, reinterpret_cast<fftw_complex *> (out), |
|
195 0, stride, dist, dir, plan_flags); |
3828
|
196 |
|
197 if (*cur_plan_p == 0) |
|
198 (*current_liboctave_error_handler) ("Error creating fftw plan"); |
|
199 } |
|
200 |
|
201 return *cur_plan_p; |
|
202 } |
|
203 |
4773
|
204 fftw_plan |
|
205 octave_fftw_planner::create_plan (const int rank, const dim_vector dims, |
|
206 int howmany, int stride, int dist, |
|
207 const double *in, Complex *out) |
3828
|
208 { |
4773
|
209 fftw_plan *cur_plan_p = &rplan; |
3828
|
210 bool create_new_plan = false; |
4808
|
211 bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out); |
3828
|
212 |
4809
|
213 // Don't create a new plan if we have a non SIMD plan already but |
|
214 // can do SIMD. This prevents endlessly recreating plans if we |
|
215 // change the alignment. |
|
216 |
4783
|
217 if (rplan == 0 || rd != dist || rs != stride || rr != rank |
4808
|
218 || rh != howmany || ((ioalign != rsimd_align) ? !ioalign : false)) |
4773
|
219 create_new_plan = true; |
|
220 else |
4809
|
221 { |
|
222 // We still might not have the same shape of array. |
|
223 |
|
224 for (int i = 0; i < rank; i++) |
|
225 if (dims(i) != rn(i)) |
|
226 { |
|
227 create_new_plan = true; |
|
228 break; |
|
229 } |
|
230 } |
3828
|
231 |
|
232 if (create_new_plan) |
|
233 { |
4773
|
234 rd = dist; |
|
235 rs = stride; |
|
236 rr = rank; |
|
237 rh = howmany; |
4808
|
238 rsimd_align = ioalign; |
4773
|
239 rn = dims; |
|
240 |
4808
|
241 if (ioalign) |
|
242 plan_flags &= ~FFTW_UNALIGNED; |
|
243 else |
|
244 plan_flags |= FFTW_UNALIGNED; |
|
245 |
3828
|
246 if (*cur_plan_p) |
4773
|
247 fftw_destroy_plan (*cur_plan_p); |
3828
|
248 |
4809
|
249 // Note reversal of dimensions for column major storage in FFTW. |
|
250 |
4773
|
251 OCTAVE_LOCAL_BUFFER (int, tmp, rank); |
4809
|
252 |
4773
|
253 for (int i = 0, j = rank-1; i < rank; i++, j--) |
|
254 tmp[i] = dims(j); |
|
255 |
|
256 *cur_plan_p = |
|
257 fftw_plan_many_dft_r2c (rank, tmp, howmany, |
|
258 (const_cast<double *> (in)), |
4774
|
259 0, stride, dist, reinterpret_cast<fftw_complex *> (out), |
|
260 0, stride, dist, plan_flags); |
3828
|
261 |
|
262 if (*cur_plan_p == 0) |
4773
|
263 (*current_liboctave_error_handler) ("Error creating fftw plan"); |
3828
|
264 } |
|
265 |
|
266 return *cur_plan_p; |
|
267 } |
|
268 |
|
269 static octave_fftw_planner fftw_planner; |
|
270 |
4775
|
271 static inline void |
|
272 convert_packcomplex_1d (Complex *out, size_t nr, size_t nc, |
|
273 int stride, int dist) |
4773
|
274 { |
4785
|
275 OCTAVE_QUIT; |
|
276 |
|
277 // Fill in the missing data. |
|
278 |
4773
|
279 for (size_t i = 0; i < nr; i++) |
|
280 for (size_t j = nc/2+1; j < nc; j++) |
|
281 out[j*stride + i*dist] = conj(out[(nc - j)*stride + i*dist]); |
4785
|
282 |
|
283 OCTAVE_QUIT; |
4773
|
284 } |
|
285 |
4775
|
286 static inline void |
|
287 convert_packcomplex_Nd (Complex *out, const dim_vector &dv) |
3828
|
288 { |
4773
|
289 size_t nc = dv(0); |
|
290 size_t nr = dv(1); |
4808
|
291 size_t np = (dv.length () > 2 ? dv.numel () / nc / nr : 1); |
4773
|
292 size_t nrp = nr * np; |
|
293 Complex *ptr1, *ptr2; |
|
294 |
4785
|
295 OCTAVE_QUIT; |
|
296 |
|
297 // Create space for the missing elements. |
|
298 |
4773
|
299 for (size_t i = 0; i < nrp; i++) |
|
300 { |
|
301 ptr1 = out + i * (nc/2 + 1) + nrp*((nc-1)/2); |
|
302 ptr2 = out + i * nc; |
|
303 for (size_t j = 0; j < nc/2+1; j++) |
|
304 *ptr2++ = *ptr1++; |
|
305 } |
|
306 |
4785
|
307 OCTAVE_QUIT; |
|
308 |
|
309 // Fill in the missing data for the rank = 2 case directly for speed. |
|
310 |
4773
|
311 for (size_t i = 0; i < np; i++) |
|
312 { |
|
313 for (size_t j = 1; j < nr; j++) |
|
314 for (size_t k = nc/2+1; k < nc; k++) |
|
315 out[k + (j + i*nr)*nc] = conj(out[nc - k + ((i+1)*nr - j)*nc]); |
|
316 |
|
317 for (size_t j = nc/2+1; j < nc; j++) |
|
318 out[j + i*nr*nc] = conj(out[(i*nr+1)*nc - j]); |
|
319 } |
|
320 |
4785
|
321 OCTAVE_QUIT; |
|
322 |
|
323 // Now do the permutations needed for rank > 2 cases. |
|
324 |
4773
|
325 size_t jstart = dv(0) * dv(1); |
|
326 size_t kstep = dv(0); |
|
327 size_t nel = dv.numel (); |
4785
|
328 |
4808
|
329 for (int inner = 2; inner < dv.length (); inner++) |
4773
|
330 { |
|
331 size_t jmax = jstart * dv(inner); |
|
332 for (size_t i = 0; i < nel; i+=jmax) |
|
333 for (size_t j = jstart, jj = jmax-jstart; j < jj; |
|
334 j+=jstart, jj-=jstart) |
|
335 for (size_t k = 0; k < jstart; k+= kstep) |
|
336 for (size_t l = nc/2+1; l < nc; l++) |
|
337 { |
|
338 Complex tmp = out[i+ j + k + l]; |
|
339 out[i + j + k + l] = out[i + jj + k + l]; |
|
340 out[i + jj + k + l] = tmp; |
|
341 } |
|
342 jstart = jmax; |
|
343 } |
4785
|
344 |
|
345 OCTAVE_QUIT; |
4773
|
346 } |
|
347 |
|
348 int |
|
349 octave_fftw::fft (const double *in, Complex *out, size_t npts, |
|
350 size_t nsamples, int stride, int dist) |
|
351 { |
|
352 dist = (dist < 0 ? npts : dist); |
|
353 |
|
354 dim_vector dv (npts); |
|
355 fftw_plan plan = fftw_planner.create_plan (1, dv, nsamples, stride, dist, |
|
356 in, out); |
|
357 |
|
358 fftw_execute_dft_r2c (plan, (const_cast<double *>(in)), |
|
359 reinterpret_cast<fftw_complex *> (out)); |
|
360 |
4809
|
361 // Need to create other half of the transform. |
|
362 |
4773
|
363 convert_packcomplex_1d (out, nsamples, npts, stride, dist); |
3828
|
364 |
|
365 return 0; |
|
366 } |
|
367 |
|
368 int |
4773
|
369 octave_fftw::fft (const Complex *in, Complex *out, size_t npts, |
|
370 size_t nsamples, int stride, int dist) |
3828
|
371 { |
4773
|
372 dist = (dist < 0 ? npts : dist); |
|
373 |
|
374 dim_vector dv (npts); |
|
375 fftw_plan plan = fftw_planner.create_plan (FFTW_FORWARD, 1, dv, nsamples, |
|
376 stride, dist, in, out); |
|
377 |
|
378 fftw_execute_dft (plan, |
|
379 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
|
380 reinterpret_cast<fftw_complex *> (out)); |
|
381 |
|
382 return 0; |
|
383 } |
|
384 |
|
385 int |
|
386 octave_fftw::ifft (const Complex *in, Complex *out, size_t npts, |
|
387 size_t nsamples, int stride, int dist) |
|
388 { |
|
389 dist = (dist < 0 ? npts : dist); |
|
390 |
|
391 dim_vector dv (npts); |
|
392 fftw_plan plan = fftw_planner.create_plan (FFTW_BACKWARD, 1, dv, nsamples, |
|
393 stride, dist, in, out); |
|
394 |
|
395 fftw_execute_dft (plan, |
|
396 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
|
397 reinterpret_cast<fftw_complex *> (out)); |
3828
|
398 |
|
399 const Complex scale = npts; |
4773
|
400 for (size_t j = 0; j < nsamples; j++) |
|
401 for (size_t i = 0; i < npts; i++) |
|
402 out[i*stride + j*dist] /= scale; |
3828
|
403 |
|
404 return 0; |
|
405 } |
|
406 |
|
407 int |
4773
|
408 octave_fftw::fftNd (const double *in, Complex *out, const int rank, |
|
409 const dim_vector &dv) |
3828
|
410 { |
4773
|
411 int dist = 1; |
|
412 for (int i = 0; i < rank; i++) |
|
413 dist *= dv(i); |
|
414 |
|
415 // Fool with the position of the start of the output matrix, so that |
4809
|
416 // creating other half of the matrix won't cause cache problems. |
|
417 |
4773
|
418 int offset = (dv.numel () / dv(0)) * ((dv(0) - 1) / 2); |
|
419 |
|
420 fftw_plan plan = fftw_planner.create_plan (rank, dv, 1, 1, dist, |
|
421 in, out + offset); |
|
422 |
|
423 fftw_execute_dft_r2c (plan, (const_cast<double *>(in)), |
|
424 reinterpret_cast<fftw_complex *> (out+ offset)); |
|
425 |
4809
|
426 // Need to create other half of the transform. |
|
427 |
4773
|
428 convert_packcomplex_Nd (out, dv); |
3828
|
429 |
|
430 return 0; |
|
431 } |
|
432 |
|
433 int |
4773
|
434 octave_fftw::fftNd (const Complex *in, Complex *out, const int rank, |
|
435 const dim_vector &dv) |
3828
|
436 { |
4773
|
437 int dist = 1; |
|
438 for (int i = 0; i < rank; i++) |
|
439 dist *= dv(i); |
|
440 |
|
441 fftw_plan plan = fftw_planner.create_plan (FFTW_FORWARD, rank, dv, 1, 1, |
|
442 dist, in, out); |
|
443 |
|
444 fftw_execute_dft (plan, |
|
445 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
|
446 reinterpret_cast<fftw_complex *> (out)); |
|
447 |
|
448 return 0; |
|
449 } |
3828
|
450 |
4773
|
451 int |
|
452 octave_fftw::ifftNd (const Complex *in, Complex *out, const int rank, |
4784
|
453 const dim_vector &dv) |
4773
|
454 { |
|
455 int dist = 1; |
|
456 for (int i = 0; i < rank; i++) |
|
457 dist *= dv(i); |
|
458 |
|
459 fftw_plan plan = fftw_planner.create_plan (FFTW_BACKWARD, rank, dv, 1, 1, |
|
460 dist, in, out); |
|
461 |
|
462 fftw_execute_dft (plan, |
|
463 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
|
464 reinterpret_cast<fftw_complex *> (out)); |
|
465 |
|
466 const size_t npts = dv.numel (); |
3828
|
467 const Complex scale = npts; |
|
468 for (size_t i = 0; i < npts; i++) |
4773
|
469 out[i] /= scale; |
3828
|
470 |
|
471 return 0; |
|
472 } |
|
473 |
|
474 #endif |
|
475 |
|
476 /* |
|
477 ;;; Local Variables: *** |
|
478 ;;; mode: C++ *** |
|
479 ;;; End: *** |
|
480 */ |
|
481 |