mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
fftT.hpp
Go to the documentation of this file.
1/** \file fftT.hpp
2 * \brief The Fast Fourier Transform interface
3 * \ingroup ft_files
4 * \author Jared R. Males (jaredmales@gmail.com)
5 *
6 */
7
8//***********************************************************************//
9// Copyright 2015-2025 Jared R. Males (jaredmales@gmail.com)
10//
11// This file is part of mxlib.
12//
13// mxlib is free software: you can redistribute it and/or modify
14// it under the terms of the GNU General Public License as published by
15// the Free Software Foundation, either version 3 of the License, or
16// (at your option) any later version.
17//
18// mxlib is distributed in the hope that it will be useful,
19// but WITHOUT ANY WARRANTY; without even the implied warranty of
20// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21// GNU General Public License for more details.
22//
23// You should have received a copy of the GNU General Public License
24// along with mxlib. If not, see <http://www.gnu.org/licenses/>.
25//***********************************************************************//
26
27#ifndef fftT_hpp
28#define fftT_hpp
29
30#pragma GCC system_header
31#include <Eigen/Dense>
32
33#include <fftw3.h>
34
35#include "fftwTemplates.hpp"
36#include "../../meta/trueFalseT.hpp"
37
38#include "ftTypes.hpp"
39
40namespace mx
41{
42namespace math
43{
44namespace ft
45{
46
47template <typename _inputT, typename _outputT, size_t _rank, int _cudaGPU = 0>
48class fftT;
49
50template <int rank>
51std::vector<int> fftwDimVec( int szX, int szY, int szZ );
52
53/// Fast Fourier Transforms using the \ref fftw "FFTW library"
54/** The \ref fftw_templates "FFTW Templates" type resolution system is used to allow the compiler
55 * to access the right plan and types for the transforms based on inputT and outputT.
56 *
57 * Planning is done with "measure" to achieve optimum performance at the cost of extra time during
58 * the planning step. This can be mitigated using "wisdom", managed by \ref fftwEnvironment.
59 *
60 * Calls to the FFTW plan functions are protected by `\#pragma omp critical` directives
61 * unless MX_FFTW_NOOMP is defined prior to the first inclusion of this file. This means that
62 * you do not need to use the `fftw_make_planner_thread_safe` facility in FFTW.
63 *
64 * \note I recommend not using fftw_make_planner_thread_safe or specifying number of threads in
65 * \ref fftwEnvironment.
66 *
67 * \todo add execute interface with fftw like signature
68 * \todo add plan interface where user passes in pointers (to avoid allocations)
69 *
70 * \tparam _inputT is the input type of the transform, can be either real or complex
71 * \tparam _outputT is the output type of the transform, can be either real or complex
72 * \tparam _rank is the rank of the transform. Limited to 1, 2, or 3.
73 *
74 * \ingroup fft
75 */
76template <typename _inputT, typename _outputT, size_t _rank>
77class fftT<_inputT, _outputT, _rank, 0>
78{
79
80 public:
81 typedef _inputT inputT;
82 typedef _outputT outputT;
83
84 static const size_t rank = _rank;
85
86 typedef typename fftwTypeSpec<inputT, outputT>::realT realT;
87
88 typedef typename fftwTypeSpec<inputT, outputT>::complexT complexT;
89
90 typedef Eigen::Array<inputT, -1, -1> eigenArrayInT;
91
92 typedef Eigen::Array<outputT, -1, -1> eigenArrayOutT;
93
94 typedef typename fftwPlanSpec<realT>::planT planT;
95
96 protected:
97 dir m_direction{ dir::forward }; ///< Direction of this FFT, either dir::forward (default) or dir::backward
98
99 int m_szX{ 0 }; ///< Size of the x dimension
100 int m_szY{ 0 }; ///< Size of the y dimension
101 int m_szZ{ 0 }; ///< size of the z dimension
102
103 planT m_plan{ nullptr }; ///< The FFTW plan object. This is a pointer, allocated by FFTW library calls.
104
105 public:
106 /// Default c'tor
108
109 /// Constructor for rank 1 FFT.
110 template <int crank = _rank>
111 fftT( int nx, ///< [in] the desired size of the FFT
112 dir ndir = dir::forward, ///< [in] [optional] direction of this FFT, either dir::forward (default) or
113 ///< dir::backward
114 bool inPlace = false, /**< [in] [optional] whether or not this is an in-place transform.
115 Default is false, out-of-place.*/
116 typename std::enable_if<crank == 1>::type * = 0 );
117
118 /// Constructor for rank 2 FFT.
119 template <int crank = _rank>
120 fftT( int nx, ///< [in] the desired x size of the FFT
121 int ny, ///< [in] the desired y size of the FFT
122 dir ndir = dir::forward, /**< [in] [optional] direction of this FFT, either dir::forward
123 (default) or dir::backward */
124 bool inPlace = false, /**< [in] [optional] whether or not this is an in-place transform.
125 Default is false, out-of-place. */
126 typename std::enable_if<crank == 2>::type * = 0 );
127
128 /// Constructor for rank 3 FFT.
129 template <int crank = _rank>
130 fftT( int nx, ///< [in] the desired x size of the FFT
131 int ny, ///< [in] the desired y size of the FFT
132 int nz, ///< [in] the desired z size of the FFT
133 dir ndir = dir::forward, /**< [in] [optional] direction of this FFT, either dir::forward (default)
134 or dir::backward*/
135 bool inPlace = false, /**< [in] [optional] whether or not this is an in-place
136 transform. Default is false, out-of-place. */
137 typename std::enable_if<crank == 3>::type * = 0 );
138
139 /// Destructor
141
142 /// Destroy (de-allocate) the plan
144
145 /// Get the direction of this FFT
146 /** The direction is either dir::forward or dir::backward.
147 *
148 * \returns the current value of m_direction.
149 */
151
152protected:
153 /// Call the FFTW planning routine
154 void doPlan( bool inPlace /**< [in] whether or not the tranform is in-place */ );
155
156public:
157
158 /// Planning routine for rank 1 transforms.
159 template <int crank = _rank>
160 void plan( int nx, ///< [in] the desired size of the FFT
161 ft::dir ndir = dir::forward, /**< [in] [optional] direction of this FFT, either dir::forward (default)
162 or dir::backward */
163 bool inPlace = false, /**< [in] [optional] whether or not this is an in-place transform.
164 Default is false, out-of-place. */
165 typename std::enable_if<crank == 1>::type * = 0 );
166
167 /// Planning routine for rank 2 transforms.
168 template <int crank = _rank>
169 void plan( int nx, ///< [in] the desired x size of the FFT
170 int ny, ///< [in] the desired y size of the FFT
171 ft::dir ndir = dir::forward, /**< [in] [optional] direction of this FFT, either dir::forward (default)
172 or dir::backward */
173 bool inPlace = false, /**< [in] [optional] whether or not this is an in-place transform.
174 Default is false, out-of-place. */
175 typename std::enable_if<crank == 2>::type * = 0 );
176
177 /// Planning routine for rank 3 transforms.
178 template <int crank = _rank>
179 void plan( int nx, ///< [in] the desired x size of the FFT
180 int ny, ///< [in] the desired y size of the FFT
181 int nz, ///< [in] the desired z size of the FFT
182 ft::dir ndir = dir::forward, /**< [in] [optional] direction of this FFT, either dir::forward
183 (default) or dir::backward */
184 bool inPlace = false, /**< [in] [optional] whether or not this is an in-place transform.
185 Default is false, out-of-place. */
186 typename std::enable_if<crank == 3>::type * = 0 );
187
188 /// Conduct the FFT
189 void operator()( outputT *out, ///< [out] the output of the FFT, must be pre-allocated
190 inputT *in ///< [in] the input to the FFT
191 ) const;
192
193 /// Conduct the FFT
194 void operator()( eigenArrayOutT &out, /**< [out] the output of the DFT */
195 eigenArrayInT &in /**< [in] the input to the DFT */ ) const;
196};
197
198template <typename inputT, typename outputT, size_t rank>
199fftT<inputT, outputT, rank, 0>::fftT()
200{
201}
202
203template <typename inputT, typename outputT, size_t rank>
204template <int crank>
205fftT<inputT, outputT, rank, 0>::fftT( int nx, ft::dir ndir, bool inPlace, typename std::enable_if<crank == 1>::type * )
206{
207 m_direction = ndir;
208
209 plan( nx, ndir, inPlace );
210}
211
212template <typename inputT, typename outputT, size_t rank>
213template <int crank>
214fftT<inputT, outputT, rank, 0>::fftT(
215 int nx, int ny, ft::dir ndir, bool inPlace, typename std::enable_if<crank == 2>::type * )
216{
217 m_direction = ndir;
218
219 plan( nx, ny, ndir, inPlace );
220}
221
222template <typename inputT, typename outputT, size_t rank>
223template <int crank>
224fftT<inputT, outputT, rank, 0>::fftT(
225 int nx, int ny, int nz, ft::dir ndir, bool inPlace, typename std::enable_if<crank == 3>::type * )
226{
227 m_direction = ndir;
228
229 plan( nx, ny, nz, ndir, inPlace );
230}
231
232template <typename inputT, typename outputT, size_t rank>
233fftT<inputT, outputT, rank, 0>::~fftT()
234{
235 destroyPlan();
236}
237
238template <typename inputT, typename outputT, size_t rank>
239void fftT<inputT, outputT, rank, 0>::destroyPlan()
240{
241 if( m_plan )
242 {
243 fftw_destroy_plan<realT>( m_plan );
244 }
245
246 m_plan = 0;
247
248 m_szX = 0;
249 m_szY = 0;
250 m_szZ = 0;
251}
252
253template <typename inputT, typename outputT, size_t rank>
254ft::dir fftT<inputT, outputT, rank, 0>::direction()
255{
256 return m_direction;
257}
258
259template <typename inputT, typename outputT, size_t rank>
260void fftT<inputT, outputT, rank, 0>::doPlan( bool inPlace )
261{
262
263 inputT *forplan1 = nullptr;
264 outputT *forplan2 = nullptr;
265
266 int sz;
267
268 if( rank == 1 )
269 {
270 sz = m_szX;
271 }
272 if( rank == 2 )
273 {
274 sz = m_szX * m_szY;
275 }
276 if( rank == 3 )
277 {
278 sz = m_szX * m_szY * m_szZ;
279 }
280
282
283 if( inPlace )
284 {
285 forplan2 = reinterpret_cast<outputT *>(forplan1);
286 }
287 else
288 {
290 }
291
292 int pdir = FFTW_FORWARD;
293 if( m_direction == dir::backward )
294 {
296 }
297
298 // clang-format off
299 #ifndef MX_FFTW_NOOMP
300 #pragma omp critical
301 {
302 #endif // clang-format on
303
304 m_plan = fftw_plan_dft<inputT, outputT>( fftwDimVec<rank>( m_szX, m_szY, m_szZ ),
305 forplan1,
306 forplan2,
307 pdir,
308 FFTW_MEASURE );
309 // clang-format off
310 #ifndef MX_FFTW_NOOMP
311 }
312 #endif // clang-format on
313
314 //Only free forplan2 if it's distinct. Have to do comparison as same type.
315 if( forplan2 != nullptr && reinterpret_cast<char *>(forplan2) != reinterpret_cast<char *>(forplan1) )
316 {
318 }
319
320 if( forplan1 != nullptr )
321 {
323 }
324}
325
326template <typename inputT, typename outputT, size_t rank>
327template <int crank>
328void fftT<inputT, outputT, rank, 0>::plan( int nx,
330 bool inPlace,
331 typename std::enable_if<crank == 1>::type * )
332{
333 if( m_szX == nx && m_direction == ndir && m_plan )
334 {
335 return;
336 }
337
338 destroyPlan();
339
340 m_direction = ndir;
341
342 m_szX = nx;
343 m_szY = 0;
344 m_szZ = 0;
345
346 doPlan( inPlace );
347}
348
349template <typename inputT, typename outputT, size_t rank>
350template <int crank>
351void fftT<inputT, outputT, rank, 0>::plan(
352 int nx, int ny, ft::dir ndir, bool inPlace, typename std::enable_if<crank == 2>::type * )
353{
354 if( m_szX == nx && m_szY == ny && m_direction == ndir && m_plan )
355 {
356 return;
357 }
358
359 destroyPlan();
360
361 m_direction = ndir;
362
363 m_szX = nx;
364 m_szY = ny;
365 m_szZ = 0;
366
367 doPlan( inPlace );
368}
369
370template <typename inputT, typename outputT, size_t rank>
371template <int crank>
372void fftT<inputT, outputT, rank, 0>::plan(
373 int nx, int ny, int nz, ft::dir ndir, bool inPlace, typename std::enable_if<crank == 3>::type * )
374{
375 if( m_szX == nx && m_szY == ny && m_szZ == nz && m_direction == ndir && m_plan )
376 {
377 return;
378 }
379
380 destroyPlan();
381
382 m_direction = ndir;
383
384 m_szX = nx;
385 m_szY = ny;
386 m_szZ = nz;
387
388 doPlan( inPlace );
389}
390
391template <typename inputT, typename outputT, size_t rank>
392void fftT<inputT, outputT, rank, 0>::operator()( outputT *out, inputT *in ) const
393{
395}
396
397template <typename inputT, typename outputT, size_t rank>
398void fftT<inputT, outputT, rank, 0>::operator()( eigenArrayOutT &out, eigenArrayInT &in ) const
399{
400 fftw_execute_dft<inputT, outputT>( m_plan, in.data(), out.data() );
401}
402
403template <>
404std::vector<int> fftwDimVec<1>( int szX, int szY, int szZ );
405
406template <>
407std::vector<int> fftwDimVec<2>( int szX, int szY, int szZ );
408
409template <>
410std::vector<int> fftwDimVec<3>( int szX, int szY, int szZ );
411
412
413} // namespace ft
414} // namespace math
415} // namespace mx
416
417#ifdef MXLIB_CUDA
418
419 #include "fftTcuda.hpp"
420
421#endif
422
423#endif // fftT_hpp
void plan(int nx, int ny, ft::dir ndir=dir::forward, bool inPlace=false, typename std::enable_if< crank==2 >::type *=0)
Planning routine for rank 2 transforms.
void operator()(outputT *out, inputT *in) const
Conduct the FFT.
void operator()(eigenArrayOutT &out, eigenArrayInT &in) const
Conduct the FFT.
void plan(int nx, ft::dir ndir=dir::forward, bool inPlace=false, typename std::enable_if< crank==1 >::type *=0)
Planning routine for rank 1 transforms.
fftT(int nx, int ny, int nz, dir ndir=dir::forward, bool inPlace=false, typename std::enable_if< crank==3 >::type *=0)
Constructor for rank 3 FFT.
void destroyPlan()
Destroy (de-allocate) the plan.
void doPlan(bool inPlace)
Call the FFTW planning routine.
fftT(int nx, int ny, dir ndir=dir::forward, bool inPlace=false, typename std::enable_if< crank==2 >::type *=0)
Constructor for rank 2 FFT.
ft::dir direction()
Get the direction of this FFT.
void plan(int nx, int ny, int nz, ft::dir ndir=dir::forward, bool inPlace=false, typename std::enable_if< crank==3 >::type *=0)
Planning routine for rank 3 transforms.
fftT(int nx, dir ndir=dir::forward, bool inPlace=false, typename std::enable_if< crank==1 >::type *=0)
Constructor for rank 1 FFT.
Declares and defines templatized wrappers for the fftw library.
Fourier Transform Types.
dir
Directions of the Fourier Transform.
Definition ftTypes.hpp:41
@ backward
Specifies the backward transform.
@ forward
Specifies the forward transform.
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
The mxlib c++ namespace.
Definition mxlib.hpp:37
fftwX_plan planT
Specializations typedef planT as fftwf_plan, fftw_plan, fftwl_plan, or fftwq_plan.
std::complex< realT > complexT
The complex data type.
_realT realT
The real data type (_realT is actually defined in specializations).