ndarray
NumPy-friendly multidimensional arrays in C++
Loading...
Searching...
No Matches
FFTWTraits.h
Go to the documentation of this file.
1// -*- c++ -*-
2/*
3 * Copyright (c) 2010-2012, Jim Bosch
4 * All rights reserved.
5 *
6 * ndarray is distributed under a simple BSD-like license;
7 * see the LICENSE file that should be present in the root
8 * of the source distribution, or alternately available at:
9 * https://github.com/ndarray/ndarray
10 */
11// THIS FILE IS MACHINE GENERATED BY SCONS. DO NOT EDIT MANUALLY.
12#ifndef NDARRAY_FFT_FFTWTraits_h_INCLUDED
13#define NDARRAY_FFT_FFTWTraits_h_INCLUDED
14
21#include <complex>
22#include <fftw3.h>
24
25namespace ndarray {
27namespace detail {
28
33template <typename T> struct FFTWTraits { BOOST_STATIC_ASSERT(sizeof(T) < 0); };
34
36
37 template <> struct FFTWTraits<float> {
38 BOOST_STATIC_ASSERT((!boost::is_const<float>::value));
39 typedef fftwf_plan Plan;
40 typedef FourierTraits<float>::ElementX ElementX;
41 typedef FourierTraits<float>::ElementK ElementK;
42 typedef boost::shared_ptr<ElementX> OwnerX;
43 typedef boost::shared_ptr<ElementK> OwnerK;
44 static inline Plan forward(int rank, const int *n, int howmany,
45 ElementX *in, const int *inembed, int istride, int idist,
46 ElementK *out, const int *onembed, int ostride, int odist,
47 unsigned flags) {
48 return fftwf_plan_many_dft_r2c(rank, n, howmany,
49 in, inembed, istride, idist,
50 reinterpret_cast<fftwf_complex*>(out),
51 onembed, ostride, odist,
52 flags);
53 }
54 static inline Plan inverse(int rank, const int *n, int howmany,
55 ElementK *in, const int *inembed, int istride, int idist,
56 ElementX *out, const int *onembed, int ostride, int odist,
57 unsigned flags) {
58 return fftwf_plan_many_dft_c2r(rank, n, howmany,
59 reinterpret_cast<fftwf_complex*>(in),
60 inembed, istride, idist,
61 out, onembed, ostride, odist,
62 flags);
63 }
64 static inline void destroy(Plan p) { fftwf_destroy_plan(p); }
65 static inline void execute(Plan p) { fftwf_execute(p); }
66 static inline OwnerX allocateX(int n) {
67 return OwnerX(
68 reinterpret_cast<ElementX*>(
69 fftwf_malloc(sizeof(ElementX)*n)
70 ),
71 fftwf_free
72 );
73 }
74 static inline OwnerK allocateK(int n) {
75 return OwnerK(
76 reinterpret_cast<ElementK*>(
77 fftwf_malloc(sizeof(ElementK)*n)
78 ),
79 fftwf_free
80 );
81 }
82 };
83 template <> struct FFTWTraits< std::complex<float> > {
84 typedef fftwf_plan Plan;
85 typedef FourierTraits< std::complex<float> >::ElementX ElementX;
86 typedef FourierTraits< std::complex<float> >::ElementK ElementK;
87 typedef boost::shared_ptr<ElementX> OwnerX;
88 typedef boost::shared_ptr<ElementK> OwnerK;
89 static inline Plan forward(int rank, const int *n, int howmany,
90 ElementX *in, const int *inembed, int istride, int idist,
91 ElementK *out, const int *onembed, int ostride, int odist,
92 unsigned flags) {
93 return fftwf_plan_many_dft(rank, n, howmany,
94 reinterpret_cast<fftwf_complex*>(in),
95 inembed, istride, idist,
96 reinterpret_cast<fftwf_complex*>(out),
97 onembed, ostride, odist,
98 FFTW_FORWARD, flags);
99 }
100 static inline Plan inverse(int rank, const int *n, int howmany,
101 ElementK *in, const int *inembed, int istride, int idist,
102 ElementX *out, const int *onembed, int ostride, int odist,
103 unsigned flags) {
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_BACKWARD,flags);
110 }
111 static inline void destroy(Plan p) { fftwf_destroy_plan(p); }
112 static inline void execute(Plan p) { fftwf_execute(p); }
113 static inline OwnerX allocateX(int n) {
114 return OwnerX(
115 reinterpret_cast<ElementX*>(
116 fftwf_malloc(sizeof(ElementX)*n)
117 ),
118 fftwf_free
119 );
120 }
121 static inline OwnerK allocateK(int n) {
122 return OwnerK(
123 reinterpret_cast<ElementK*>(
124 fftwf_malloc(sizeof(ElementK)*n)
125 ),
126 fftwf_free
127 );
128 }
129 };
130
131 template <> struct FFTWTraits<double> {
132 BOOST_STATIC_ASSERT((!boost::is_const<double>::value));
133 typedef fftw_plan Plan;
134 typedef FourierTraits<double>::ElementX ElementX;
135 typedef FourierTraits<double>::ElementK ElementK;
136 typedef boost::shared_ptr<ElementX> OwnerX;
137 typedef boost::shared_ptr<ElementK> OwnerK;
138 static inline Plan forward(int rank, const int *n, int howmany,
139 ElementX *in, const int *inembed, int istride, int idist,
140 ElementK *out, const int *onembed, int ostride, int odist,
141 unsigned flags) {
142 return fftw_plan_many_dft_r2c(rank, n, howmany,
143 in, inembed, istride, idist,
144 reinterpret_cast<fftw_complex*>(out),
145 onembed, ostride, odist,
146 flags);
147 }
148 static inline Plan inverse(int rank, const int *n, int howmany,
149 ElementK *in, const int *inembed, int istride, int idist,
150 ElementX *out, const int *onembed, int ostride, int odist,
151 unsigned flags) {
152 return fftw_plan_many_dft_c2r(rank, n, howmany,
153 reinterpret_cast<fftw_complex*>(in),
154 inembed, istride, idist,
155 out, onembed, ostride, odist,
156 flags);
157 }
158 static inline void destroy(Plan p) { fftw_destroy_plan(p); }
159 static inline void execute(Plan p) { fftw_execute(p); }
160 static inline OwnerX allocateX(int n) {
161 return OwnerX(
162 reinterpret_cast<ElementX*>(
163 fftw_malloc(sizeof(ElementX)*n)
164 ),
165 fftw_free
166 );
167 }
168 static inline OwnerK allocateK(int n) {
169 return OwnerK(
170 reinterpret_cast<ElementK*>(
171 fftw_malloc(sizeof(ElementK)*n)
172 ),
173 fftw_free
174 );
175 }
176 };
177 template <> struct FFTWTraits< std::complex<double> > {
178 typedef fftw_plan Plan;
179 typedef FourierTraits< std::complex<double> >::ElementX ElementX;
180 typedef FourierTraits< std::complex<double> >::ElementK ElementK;
181 typedef boost::shared_ptr<ElementX> OwnerX;
182 typedef boost::shared_ptr<ElementK> OwnerK;
183 static inline Plan forward(int rank, const int *n, int howmany,
184 ElementX *in, const int *inembed, int istride, int idist,
185 ElementK *out, const int *onembed, int ostride, int odist,
186 unsigned flags) {
187 return fftw_plan_many_dft(rank, n, howmany,
188 reinterpret_cast<fftw_complex*>(in),
189 inembed, istride, idist,
190 reinterpret_cast<fftw_complex*>(out),
191 onembed, ostride, odist,
192 FFTW_FORWARD, flags);
193 }
194 static inline Plan inverse(int rank, const int *n, int howmany,
195 ElementK *in, const int *inembed, int istride, int idist,
196 ElementX *out, const int *onembed, int ostride, int odist,
197 unsigned flags) {
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_BACKWARD,flags);
204 }
205 static inline void destroy(Plan p) { fftw_destroy_plan(p); }
206 static inline void execute(Plan p) { fftw_execute(p); }
207 static inline OwnerX allocateX(int n) {
208 return OwnerX(
209 reinterpret_cast<ElementX*>(
210 fftw_malloc(sizeof(ElementX)*n)
211 ),
212 fftw_free
213 );
214 }
215 static inline OwnerK allocateK(int n) {
216 return OwnerK(
217 reinterpret_cast<ElementK*>(
218 fftw_malloc(sizeof(ElementK)*n)
219 ),
220 fftw_free
221 );
222 }
223 };
224#ifndef NDARRAY_FFT_NO_LONG_DOUBLE
225
226 template <> struct FFTWTraits<long double> {
227 BOOST_STATIC_ASSERT((!boost::is_const<long double>::value));
228 typedef fftwl_plan Plan;
229 typedef FourierTraits<long double>::ElementX ElementX;
230 typedef FourierTraits<long double>::ElementK ElementK;
231 typedef boost::shared_ptr<ElementX> OwnerX;
232 typedef boost::shared_ptr<ElementK> OwnerK;
233 static inline Plan forward(int rank, const int *n, int howmany,
234 ElementX *in, const int *inembed, int istride, int idist,
235 ElementK *out, const int *onembed, int ostride, int odist,
236 unsigned flags) {
237 return fftwl_plan_many_dft_r2c(rank, n, howmany,
238 in, inembed, istride, idist,
239 reinterpret_cast<fftwl_complex*>(out),
240 onembed, ostride, odist,
241 flags);
242 }
243 static inline Plan inverse(int rank, const int *n, int howmany,
244 ElementK *in, const int *inembed, int istride, int idist,
245 ElementX *out, const int *onembed, int ostride, int odist,
246 unsigned flags) {
247 return fftwl_plan_many_dft_c2r(rank, n, howmany,
248 reinterpret_cast<fftwl_complex*>(in),
249 inembed, istride, idist,
250 out, onembed, ostride, odist,
251 flags);
252 }
253 static inline void destroy(Plan p) { fftwl_destroy_plan(p); }
254 static inline void execute(Plan p) { fftwl_execute(p); }
255 static inline OwnerX allocateX(int n) {
256 return OwnerX(
257 reinterpret_cast<ElementX*>(
258 fftwl_malloc(sizeof(ElementX)*n)
259 ),
260 fftwl_free
261 );
262 }
263 static inline OwnerK allocateK(int n) {
264 return OwnerK(
265 reinterpret_cast<ElementK*>(
266 fftwl_malloc(sizeof(ElementK)*n)
267 ),
268 fftwl_free
269 );
270 }
271 };
272 template <> struct FFTWTraits< std::complex<long double> > {
273 typedef fftwl_plan Plan;
274 typedef FourierTraits< std::complex<long double> >::ElementX ElementX;
275 typedef FourierTraits< std::complex<long double> >::ElementK ElementK;
276 typedef boost::shared_ptr<ElementX> OwnerX;
277 typedef boost::shared_ptr<ElementK> OwnerK;
278 static inline Plan forward(int rank, const int *n, int howmany,
279 ElementX *in, const int *inembed, int istride, int idist,
280 ElementK *out, const int *onembed, int ostride, int odist,
281 unsigned flags) {
282 return fftwl_plan_many_dft(rank, n, howmany,
283 reinterpret_cast<fftwl_complex*>(in),
284 inembed, istride, idist,
285 reinterpret_cast<fftwl_complex*>(out),
286 onembed, ostride, odist,
287 FFTW_FORWARD, flags);
288 }
289 static inline Plan inverse(int rank, const int *n, int howmany,
290 ElementK *in, const int *inembed, int istride, int idist,
291 ElementX *out, const int *onembed, int ostride, int odist,
292 unsigned flags) {
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_BACKWARD,flags);
299 }
300 static inline void destroy(Plan p) { fftwl_destroy_plan(p); }
301 static inline void execute(Plan p) { fftwl_execute(p); }
302 static inline OwnerX allocateX(int n) {
303 return OwnerX(
304 reinterpret_cast<ElementX*>(
305 fftwl_malloc(sizeof(ElementX)*n)
306 ),
307 fftwl_free
308 );
309 }
310 static inline OwnerK allocateK(int n) {
311 return OwnerK(
312 reinterpret_cast<ElementK*>(
313 fftwl_malloc(sizeof(ElementK)*n)
314 ),
315 fftwl_free
316 );
317 }
318 };
319#endif
321
322} // namespace detail
324} // namespace ndarray
325
326#endif // !NDARRAY_FFT_FFTWTraits_h_INCLUDED
Traits classes to handle real-data and complex-data FFTs in a template-friendly way.