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 |
|
25 #ifdef HAVE_FFTW |
|
26 |
|
27 #include "oct-fftw.h" |
|
28 #include "lo-error.h" |
|
29 |
|
30 |
|
31 // Helper class to create and cache fftw plans for both 1d and 2d. This |
|
32 // implementation uses FFTW_ESTIMATE to create the plans, which in theory |
|
33 // is suboptimal, but provides quite reasonable performance. Future |
|
34 // enhancement will be to add a dynamically loadable interface ("fftw") |
|
35 // to manipulate fftw wisdom so that users may choose the appropriate |
|
36 // planner. |
|
37 |
|
38 class |
|
39 octave_fftw_planner |
|
40 { |
|
41 public: |
|
42 octave_fftw_planner (); |
|
43 |
|
44 fftw_plan create_plan (fftw_direction, size_t); |
|
45 fftwnd_plan create_plan2d (fftw_direction, size_t, size_t); |
|
46 |
|
47 private: |
|
48 int plan_flags; |
|
49 |
|
50 fftw_plan plan[2]; |
|
51 fftwnd_plan plan2d[2]; |
|
52 |
|
53 size_t n[2]; |
|
54 size_t n2d[2][2]; |
|
55 }; |
|
56 |
|
57 octave_fftw_planner::octave_fftw_planner () |
|
58 { |
|
59 plan_flags = FFTW_ESTIMATE; |
|
60 |
|
61 plan[0] = plan[1] = 0; |
|
62 plan2d[0] = plan2d[1] = 0; |
|
63 |
|
64 n[0] = n[1] = 0; |
|
65 n2d[0][0] = n2d[0][1] = n2d[1][0] = n2d[1][1] = 0; |
|
66 } |
|
67 |
|
68 fftw_plan |
|
69 octave_fftw_planner::create_plan (fftw_direction dir, size_t npts) |
|
70 { |
|
71 size_t which = (dir == FFTW_FORWARD) ? 0 : 1; |
|
72 fftw_plan *cur_plan_p = &plan[which]; |
|
73 bool create_new_plan = false; |
|
74 |
|
75 if (plan[which] == 0 || n[which] != npts) |
|
76 { |
|
77 create_new_plan = true; |
|
78 n[which] = npts; |
|
79 } |
|
80 |
|
81 if (create_new_plan) |
|
82 { |
|
83 if (*cur_plan_p) |
|
84 fftw_destroy_plan (*cur_plan_p); |
|
85 |
|
86 *cur_plan_p = fftw_create_plan (npts, dir, plan_flags); |
|
87 |
|
88 if (*cur_plan_p == 0) |
|
89 (*current_liboctave_error_handler) ("Error creating fftw plan"); |
|
90 } |
|
91 |
|
92 return *cur_plan_p; |
|
93 } |
|
94 |
|
95 fftwnd_plan |
|
96 octave_fftw_planner::create_plan2d (fftw_direction dir, |
|
97 size_t nrows, size_t ncols) |
|
98 { |
|
99 size_t which = (dir == FFTW_FORWARD) ? 0 : 1; |
|
100 fftwnd_plan *cur_plan_p = &plan2d[which]; |
|
101 bool create_new_plan = false; |
|
102 |
|
103 if (plan2d[which] == 0 || n2d[which][0] != nrows || n2d[which][1] != ncols) |
|
104 { |
|
105 create_new_plan = true; |
|
106 |
|
107 n2d[which][0] = nrows; |
|
108 n2d[which][1] = ncols; |
|
109 } |
|
110 |
|
111 if (create_new_plan) |
|
112 { |
|
113 if (*cur_plan_p) |
|
114 fftwnd_destroy_plan (*cur_plan_p); |
|
115 |
|
116 *cur_plan_p = fftw2d_create_plan (nrows, ncols, dir, |
|
117 plan_flags | FFTW_IN_PLACE); |
|
118 |
|
119 if (*cur_plan_p == 0) |
|
120 (*current_liboctave_error_handler) ("Error creating 2d fftw plan"); |
|
121 } |
|
122 |
|
123 return *cur_plan_p; |
|
124 } |
|
125 |
|
126 static octave_fftw_planner fftw_planner; |
|
127 |
|
128 int |
|
129 octave_fftw::fft (const Complex *in, Complex *out, size_t npts) |
|
130 { |
|
131 fftw_one (fftw_planner.create_plan (FFTW_FORWARD, npts), |
|
132 reinterpret_cast<fftw_complex *> (const_cast<Complex *> (in)), |
|
133 reinterpret_cast<fftw_complex *> (out)); |
|
134 |
|
135 return 0; |
|
136 } |
|
137 |
|
138 int |
|
139 octave_fftw::ifft (const Complex *in, Complex *out, size_t npts) |
|
140 { |
|
141 fftw_one (fftw_planner.create_plan (FFTW_BACKWARD, npts), |
|
142 reinterpret_cast<fftw_complex *> (const_cast<Complex *> (in)), |
|
143 reinterpret_cast<fftw_complex *> (out)); |
|
144 |
|
145 const Complex scale = npts; |
|
146 for (size_t i = 0; i < npts; i++) |
|
147 out[i] /= scale; |
|
148 |
|
149 return 0; |
|
150 } |
|
151 |
|
152 int |
|
153 octave_fftw::fft2d (Complex *inout, size_t nr, size_t nc) |
|
154 { |
|
155 fftwnd_one (fftw_planner.create_plan2d (FFTW_FORWARD, nr, nc), |
|
156 reinterpret_cast<fftw_complex *> (inout), |
3874
|
157 0); |
3828
|
158 |
|
159 return 0; |
|
160 } |
|
161 |
|
162 int |
|
163 octave_fftw::ifft2d (Complex *inout, size_t nr, size_t nc) |
|
164 { |
|
165 fftwnd_one (fftw_planner.create_plan2d (FFTW_BACKWARD, nr, nc), |
|
166 reinterpret_cast<fftw_complex *> (inout), |
3874
|
167 0); |
3828
|
168 |
|
169 const size_t npts = nr * nc; |
|
170 const Complex scale = npts; |
|
171 for (size_t i = 0; i < npts; i++) |
|
172 inout[i] /= scale; |
|
173 |
|
174 return 0; |
|
175 } |
|
176 |
|
177 #endif |
|
178 |
|
179 /* |
|
180 ;;; Local Variables: *** |
|
181 ;;; mode: C++ *** |
|
182 ;;; End: *** |
|
183 */ |
|
184 |