23 #ifndef NDARRAY_FFT_FFTWTraits_h_INCLUDED
24 #define NDARRAY_FFT_FFTWTraits_h_INCLUDED
44 template <
typename T>
struct FFTWTraits { BOOST_STATIC_ASSERT(
sizeof(T) < 0); };
48 template <>
struct FFTWTraits<float> {
49 BOOST_STATIC_ASSERT((!boost::is_const<float>::value));
50 typedef fftwf_plan Plan;
51 typedef FourierTraits<float>::ElementX ElementX;
52 typedef FourierTraits<float>::ElementK ElementK;
53 typedef boost::shared_ptr<ElementX> OwnerX;
54 typedef boost::shared_ptr<ElementK> OwnerK;
55 static inline Plan forward(
int rank,
const int *n,
int howmany,
56 ElementX *in,
const int *inembed,
int istride,
int idist,
57 ElementK *out,
const int *onembed,
int ostride,
int odist,
59 return fftwf_plan_many_dft_r2c(rank, n, howmany,
60 in, inembed, istride, idist,
61 reinterpret_cast<fftwf_complex*>(out),
62 onembed, ostride, odist,
65 static inline Plan inverse(
int rank,
const int *n,
int howmany,
66 ElementK *in,
const int *inembed,
int istride,
int idist,
67 ElementX *out,
const int *onembed,
int ostride,
int odist,
69 return fftwf_plan_many_dft_c2r(rank, n, howmany,
70 reinterpret_cast<fftwf_complex*>(in),
71 inembed, istride, idist,
72 out, onembed, ostride, odist,
75 static inline void destroy(Plan p) { fftwf_destroy_plan(p); }
76 static inline void execute(Plan p) { fftwf_execute(p); }
77 static inline OwnerX allocateX(
int n) {
79 reinterpret_cast<ElementX*>(
80 fftwf_malloc(
sizeof(ElementX)*n)
85 static inline OwnerK allocateK(
int n) {
87 reinterpret_cast<ElementK*>(
88 fftwf_malloc(
sizeof(ElementK)*n)
94 template <>
struct FFTWTraits< std::complex<float> > {
95 typedef fftwf_plan Plan;
96 typedef FourierTraits< std::complex<float> >::ElementX ElementX;
97 typedef FourierTraits< std::complex<float> >::ElementK ElementK;
98 typedef boost::shared_ptr<ElementX> OwnerX;
99 typedef boost::shared_ptr<ElementK> OwnerK;
100 static inline Plan forward(
int rank,
const int *n,
int howmany,
101 ElementX *in,
const int *inembed,
int istride,
int idist,
102 ElementK *out,
const int *onembed,
int ostride,
int odist,
104 return fftwf_plan_many_dft(rank, n, howmany,
105 reinterpret_cast<fftwf_complex*>(in),
106 inembed, istride, idist,
107 reinterpret_cast<fftwf_complex*>(out),
108 onembed, ostride, odist,
109 FFTW_FORWARD, flags);
111 static inline Plan inverse(
int rank,
const int *n,
int howmany,
112 ElementK *in,
const int *inembed,
int istride,
int idist,
113 ElementX *out,
const int *onembed,
int ostride,
int odist,
115 return fftwf_plan_many_dft(rank, n, howmany,
116 reinterpret_cast<fftwf_complex*>(in),
117 inembed, istride, idist,
118 reinterpret_cast<fftwf_complex*>(out),
119 onembed, ostride, odist,
120 FFTW_BACKWARD,flags);
122 static inline void destroy(Plan p) { fftwf_destroy_plan(p); }
123 static inline void execute(Plan p) { fftwf_execute(p); }
124 static inline OwnerX allocateX(
int n) {
126 reinterpret_cast<ElementX*>(
127 fftwf_malloc(
sizeof(ElementX)*n)
132 static inline OwnerK allocateK(
int n) {
134 reinterpret_cast<ElementK*>(
135 fftwf_malloc(
sizeof(ElementK)*n)
142 template <>
struct FFTWTraits<double> {
143 BOOST_STATIC_ASSERT((!boost::is_const<double>::value));
144 typedef fftw_plan Plan;
145 typedef FourierTraits<double>::ElementX ElementX;
146 typedef FourierTraits<double>::ElementK ElementK;
147 typedef boost::shared_ptr<ElementX> OwnerX;
148 typedef boost::shared_ptr<ElementK> OwnerK;
149 static inline Plan forward(
int rank,
const int *n,
int howmany,
150 ElementX *in,
const int *inembed,
int istride,
int idist,
151 ElementK *out,
const int *onembed,
int ostride,
int odist,
153 return fftw_plan_many_dft_r2c(rank, n, howmany,
154 in, inembed, istride, idist,
155 reinterpret_cast<fftw_complex*>(out),
156 onembed, ostride, odist,
159 static inline Plan inverse(
int rank,
const int *n,
int howmany,
160 ElementK *in,
const int *inembed,
int istride,
int idist,
161 ElementX *out,
const int *onembed,
int ostride,
int odist,
163 return fftw_plan_many_dft_c2r(rank, n, howmany,
164 reinterpret_cast<fftw_complex*>(in),
165 inembed, istride, idist,
166 out, onembed, ostride, odist,
169 static inline void destroy(Plan p) { fftw_destroy_plan(p); }
170 static inline void execute(Plan p) { fftw_execute(p); }
171 static inline OwnerX allocateX(
int n) {
173 reinterpret_cast<ElementX*>(
174 fftw_malloc(
sizeof(ElementX)*n)
179 static inline OwnerK allocateK(
int n) {
181 reinterpret_cast<ElementK*>(
182 fftw_malloc(
sizeof(ElementK)*n)
188 template <>
struct FFTWTraits< std::complex<double> > {
189 typedef fftw_plan Plan;
190 typedef FourierTraits< std::complex<double> >::ElementX ElementX;
191 typedef FourierTraits< std::complex<double> >::ElementK ElementK;
192 typedef boost::shared_ptr<ElementX> OwnerX;
193 typedef boost::shared_ptr<ElementK> OwnerK;
194 static inline Plan forward(
int rank,
const int *n,
int howmany,
195 ElementX *in,
const int *inembed,
int istride,
int idist,
196 ElementK *out,
const int *onembed,
int ostride,
int odist,
198 return fftw_plan_many_dft(rank, n, howmany,
199 reinterpret_cast<fftw_complex*>(in),
200 inembed, istride, idist,
201 reinterpret_cast<fftw_complex*>(out),
202 onembed, ostride, odist,
203 FFTW_FORWARD, flags);
205 static inline Plan inverse(
int rank,
const int *n,
int howmany,
206 ElementK *in,
const int *inembed,
int istride,
int idist,
207 ElementX *out,
const int *onembed,
int ostride,
int odist,
209 return fftw_plan_many_dft(rank, n, howmany,
210 reinterpret_cast<fftw_complex*>(in),
211 inembed, istride, idist,
212 reinterpret_cast<fftw_complex*>(out),
213 onembed, ostride, odist,
214 FFTW_BACKWARD,flags);
216 static inline void destroy(Plan p) { fftw_destroy_plan(p); }
217 static inline void execute(Plan p) { fftw_execute(p); }
218 static inline OwnerX allocateX(
int n) {
220 reinterpret_cast<ElementX*>(
221 fftw_malloc(
sizeof(ElementX)*n)
226 static inline OwnerK allocateK(
int n) {
228 reinterpret_cast<ElementK*>(
229 fftw_malloc(
sizeof(ElementK)*n)
235 #ifndef NDARRAY_FFT_NO_LONG_DOUBLE
237 template <>
struct FFTWTraits<long double> {
238 BOOST_STATIC_ASSERT((!boost::is_const<long double>::value));
239 typedef fftwl_plan Plan;
240 typedef FourierTraits<long double>::ElementX ElementX;
241 typedef FourierTraits<long double>::ElementK ElementK;
242 typedef boost::shared_ptr<ElementX> OwnerX;
243 typedef boost::shared_ptr<ElementK> OwnerK;
244 static inline Plan forward(
int rank,
const int *n,
int howmany,
245 ElementX *in,
const int *inembed,
int istride,
int idist,
246 ElementK *out,
const int *onembed,
int ostride,
int odist,
248 return fftwl_plan_many_dft_r2c(rank, n, howmany,
249 in, inembed, istride, idist,
250 reinterpret_cast<fftwl_complex*>(out),
251 onembed, ostride, odist,
254 static inline Plan inverse(
int rank,
const int *n,
int howmany,
255 ElementK *in,
const int *inembed,
int istride,
int idist,
256 ElementX *out,
const int *onembed,
int ostride,
int odist,
258 return fftwl_plan_many_dft_c2r(rank, n, howmany,
259 reinterpret_cast<fftwl_complex*>(in),
260 inembed, istride, idist,
261 out, onembed, ostride, odist,
264 static inline void destroy(Plan p) { fftwl_destroy_plan(p); }
265 static inline void execute(Plan p) { fftwl_execute(p); }
266 static inline OwnerX allocateX(
int n) {
268 reinterpret_cast<ElementX*>(
269 fftwl_malloc(
sizeof(ElementX)*n)
274 static inline OwnerK allocateK(
int n) {
276 reinterpret_cast<ElementK*>(
277 fftwl_malloc(
sizeof(ElementK)*n)
283 template <>
struct FFTWTraits< std::complex<long double> > {
284 typedef fftwl_plan Plan;
285 typedef FourierTraits< std::complex<long double> >::ElementX ElementX;
286 typedef FourierTraits< std::complex<long double> >::ElementK ElementK;
287 typedef boost::shared_ptr<ElementX> OwnerX;
288 typedef boost::shared_ptr<ElementK> OwnerK;
289 static inline Plan forward(
int rank,
const int *n,
int howmany,
290 ElementX *in,
const int *inembed,
int istride,
int idist,
291 ElementK *out,
const int *onembed,
int ostride,
int odist,
293 return fftwl_plan_many_dft(rank, n, howmany,
294 reinterpret_cast<fftwl_complex*>(in),
295 inembed, istride, idist,
296 reinterpret_cast<fftwl_complex*>(out),
297 onembed, ostride, odist,
298 FFTW_FORWARD, flags);
300 static inline Plan inverse(
int rank,
const int *n,
int howmany,
301 ElementK *in,
const int *inembed,
int istride,
int idist,
302 ElementX *out,
const int *onembed,
int ostride,
int odist,
304 return fftwl_plan_many_dft(rank, n, howmany,
305 reinterpret_cast<fftwl_complex*>(in),
306 inembed, istride, idist,
307 reinterpret_cast<fftwl_complex*>(out),
308 onembed, ostride, odist,
309 FFTW_BACKWARD,flags);
311 static inline void destroy(Plan p) { fftwl_destroy_plan(p); }
312 static inline void execute(Plan p) { fftwl_execute(p); }
313 static inline OwnerX allocateX(
int n) {
315 reinterpret_cast<ElementX*>(
316 fftwl_malloc(
sizeof(ElementX)*n)
321 static inline OwnerK allocateK(
int n) {
323 reinterpret_cast<ElementK*>(
324 fftwl_malloc(
sizeof(ElementK)*n)
337 #endif // !NDARRAY_FFT_FFTWTraits_h_INCLUDED
Traits classes to handle real-data and complex-data FFTs in a template-friendly way.