LSSTApplications  10.0+286,10.0+36,10.0+46,10.0-2-g4f67435,10.1+152,10.1+37,11.0,11.0+1,11.0-1-g47edd16,11.0-1-g60db491,11.0-1-g7418c06,11.0-2-g04d2804,11.0-2-g68503cd,11.0-2-g818369d,11.0-2-gb8b8ce7
LSSTDataManagementBasePackage
FFTWTraits.h
Go to the documentation of this file.
1 // -*- lsst-c++ -*-
2 /*
3  * LSST Data Management System
4  * Copyright 2008, 2009, 2010, 2011 LSST Corporation.
5  *
6  * This product includes software developed by the
7  * LSST Project (http://www.lsst.org/).
8  *
9  * This program is free software: you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation, either version 3 of the License, or
12  * (at your option) any later version.
13  *
14  * This program is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the LSST License Statement and
20  * the GNU General Public License along with this program. If not,
21  * see <http://www.lsstcorp.org/LegalNotices/>.
22  */
23 #ifndef NDARRAY_FFT_FFTWTraits_h_INCLUDED
24 #define NDARRAY_FFT_FFTWTraits_h_INCLUDED
25 
32 #include <complex>
33 #include <fftw3.h>
35 
36 namespace ndarray {
38 namespace detail {
39 
44 template <typename T> struct FFTWTraits { BOOST_STATIC_ASSERT(sizeof(T) < 0); };
45 
47 
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,
58  unsigned flags) {
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,
63  flags);
64  }
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,
68  unsigned flags) {
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,
73  flags);
74  }
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) {
78  return OwnerX(
79  reinterpret_cast<ElementX*>(
80  fftwf_malloc(sizeof(ElementX)*n)
81  ),
82  fftwf_free
83  );
84  }
85  static inline OwnerK allocateK(int n) {
86  return OwnerK(
87  reinterpret_cast<ElementK*>(
88  fftwf_malloc(sizeof(ElementK)*n)
89  ),
90  fftwf_free
91  );
92  }
93  };
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,
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_FORWARD, flags);
110  }
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,
114  unsigned flags) {
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);
121  }
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) {
125  return OwnerX(
126  reinterpret_cast<ElementX*>(
127  fftwf_malloc(sizeof(ElementX)*n)
128  ),
129  fftwf_free
130  );
131  }
132  static inline OwnerK allocateK(int n) {
133  return OwnerK(
134  reinterpret_cast<ElementK*>(
135  fftwf_malloc(sizeof(ElementK)*n)
136  ),
137  fftwf_free
138  );
139  }
140  };
141 
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,
152  unsigned flags) {
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,
157  flags);
158  }
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,
162  unsigned flags) {
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,
167  flags);
168  }
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) {
172  return OwnerX(
173  reinterpret_cast<ElementX*>(
174  fftw_malloc(sizeof(ElementX)*n)
175  ),
176  fftw_free
177  );
178  }
179  static inline OwnerK allocateK(int n) {
180  return OwnerK(
181  reinterpret_cast<ElementK*>(
182  fftw_malloc(sizeof(ElementK)*n)
183  ),
184  fftw_free
185  );
186  }
187  };
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,
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_FORWARD, flags);
204  }
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,
208  unsigned flags) {
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);
215  }
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) {
219  return OwnerX(
220  reinterpret_cast<ElementX*>(
221  fftw_malloc(sizeof(ElementX)*n)
222  ),
223  fftw_free
224  );
225  }
226  static inline OwnerK allocateK(int n) {
227  return OwnerK(
228  reinterpret_cast<ElementK*>(
229  fftw_malloc(sizeof(ElementK)*n)
230  ),
231  fftw_free
232  );
233  }
234  };
235 #ifndef NDARRAY_FFT_NO_LONG_DOUBLE
236 
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,
247  unsigned flags) {
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,
252  flags);
253  }
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,
257  unsigned flags) {
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,
262  flags);
263  }
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) {
267  return OwnerX(
268  reinterpret_cast<ElementX*>(
269  fftwl_malloc(sizeof(ElementX)*n)
270  ),
271  fftwl_free
272  );
273  }
274  static inline OwnerK allocateK(int n) {
275  return OwnerK(
276  reinterpret_cast<ElementK*>(
277  fftwl_malloc(sizeof(ElementK)*n)
278  ),
279  fftwl_free
280  );
281  }
282  };
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,
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_FORWARD, flags);
299  }
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,
303  unsigned flags) {
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);
310  }
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) {
314  return OwnerX(
315  reinterpret_cast<ElementX*>(
316  fftwl_malloc(sizeof(ElementX)*n)
317  ),
318  fftwl_free
319  );
320  }
321  static inline OwnerK allocateK(int n) {
322  return OwnerK(
323  reinterpret_cast<ElementK*>(
324  fftwl_malloc(sizeof(ElementK)*n)
325  ),
326  fftwl_free
327  );
328  }
329  };
330 #endif
331 
333 } // namespace detail
335 } // namespace ndarray
336 
337 #endif // !NDARRAY_FFT_FFTWTraits_h_INCLUDED
Traits classes to handle real-data and complex-data FFTs in a template-friendly way.