30#pragma GCC system_header
36#include "../../meta/trueFalseT.hpp"
47template <
typename _inputT,
typename _outputT,
size_t _rank,
int _cudaGPU = 0>
51std::vector<int> fftwDimVec(
int szX,
int szY,
int szZ );
76template <
typename _inputT,
typename _outputT,
size_t _rank>
84 static const size_t rank =
_rank;
90 typedef Eigen::Array<inputT, -1, -1> eigenArrayInT;
92 typedef Eigen::Array<outputT, -1, -1> eigenArrayOutT;
103 planT m_plan{
nullptr };
110 template <
int crank = _rank>
116 typename std::enable_if<crank == 1>::type * = 0 );
119 template <
int crank = _rank>
126 typename std::enable_if<crank == 2>::type * = 0 );
129 template <
int crank = _rank>
137 typename std::enable_if<crank == 3>::type * = 0 );
159 template <
int crank = _rank>
165 typename std::enable_if<crank == 1>::type * = 0 );
168 template <
int crank = _rank>
175 typename std::enable_if<crank == 2>::type * = 0 );
178 template <
int crank = _rank>
186 typename std::enable_if<crank == 3>::type * = 0 );
195 eigenArrayInT &
in )
const;
198template <
typename inputT,
typename outputT,
size_t rank>
199fftT<inputT, outputT, rank, 0>::fftT()
203template <
typename inputT,
typename outputT,
size_t rank>
205fftT<inputT, outputT, rank, 0>::fftT(
int nx,
ft::dir ndir,
bool inPlace,
typename std::enable_if<crank == 1>::type * )
212template <
typename inputT,
typename outputT,
size_t rank>
214fftT<inputT, outputT, rank, 0>::fftT(
215 int nx,
int ny,
ft::dir ndir,
bool inPlace,
typename std::enable_if<crank == 2>::type * )
222template <
typename inputT,
typename outputT,
size_t rank>
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 * )
232template <
typename inputT,
typename outputT,
size_t rank>
233fftT<inputT, outputT, rank, 0>::~fftT()
238template <
typename inputT,
typename outputT,
size_t rank>
239void fftT<inputT, outputT, rank, 0>::destroyPlan()
250template <
typename inputT,
typename outputT,
size_t rank>
251ft::dir fftT<inputT, outputT, rank, 0>::dir()
256template <
typename inputT,
typename outputT,
size_t rank>
257void fftT<inputT, outputT, rank, 0>::doPlan(
const meta::trueFalseT<false> &
inPlace )
271 sz = m_szX * m_szY * m_szZ;
292template <
typename inputT,
typename outputT,
size_t rank>
293void fftT<inputT, outputT, rank, 0>::doPlan(
const meta::trueFalseT<true> &
inPlace )
306 sz = m_szX * m_szY * m_szZ;
319 reinterpret_cast<inputT *
>(
forplan ),
320 reinterpret_cast<outputT *
>(
forplan ),
328template <
typename inputT,
typename outputT,
size_t rank>
330void fftT<inputT, outputT, rank, 0>::plan(
int nx,
ft::dir ndir,
bool inPlace,
typename std::enable_if<crank == 1>::type * )
332 if( m_szX == nx && m_dir ==
ndir && m_plan )
347 doPlan( meta::trueFalseT<false>() );
351 doPlan( meta::trueFalseT<true>() );
355template <
typename inputT,
typename outputT,
size_t rank>
357void fftT<inputT, outputT, rank, 0>::plan(
358 int nx,
int ny,
ft::dir ndir,
bool inPlace,
typename std::enable_if<crank == 2>::type * )
360 if( m_szX == nx && m_szY == ny && m_dir ==
ndir && m_plan )
375 doPlan( meta::trueFalseT<false>() );
379 doPlan( meta::trueFalseT<true>() );
383template <
typename inputT,
typename outputT,
size_t rank>
385void fftT<inputT, outputT, rank, 0>::plan(
386 int nx,
int ny,
int nz,
ft::dir ndir,
bool inPlace,
typename std::enable_if<crank == 3>::type * )
388 if( m_szX == nx && m_szY == ny && m_szZ ==
nz && m_dir ==
ndir && m_plan )
403 doPlan( meta::trueFalseT<false>() );
407 doPlan( meta::trueFalseT<true>() );
411template <
typename inputT,
typename outputT,
size_t rank>
412void fftT<inputT, outputT, rank, 0>::operator()( outputT *
out, inputT *
in )
const
417template <
typename inputT,
typename outputT,
size_t rank>
418void fftT<inputT, outputT, rank, 0>::operator()( eigenArrayOutT &
out, eigenArrayInT &
in )
const
424std::vector<int> fftwDimVec<1>(
int szX,
int szY,
int szZ );
427std::vector<int> fftwDimVec<2>(
int szX,
int szY,
int szZ );
430std::vector<int> fftwDimVec<3>(
int szX,
int szY,
int szZ );
435class fft<std::complex<float>, std::complex<float>, 2, 1>
450 fft(
int nx,
int ny )
458 void plan(
int nx,
int ny )
460 if( m_szX == nx && m_szY == ny )
468 cufftPlan2d( &m_plan, m_szX, m_szY, CUFFT_C2C );
471 void fwdFFT( std::complex<float> *in, std::complex<float> *out )
473 cufftExecC2C( m_plan, (cufftComplex *)in, (cufftComplex *)out, CUFFT_FORWARD );
Declares and defines templatized wrappers for the fftw library.
dir
Directions of the Fourier Transform.
@ backward
Specifies the backward transform.
@ forward
Specifies the forward transform.
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
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).