PSCF v1.3
|
Fourier transform wrapper. More...
#include <FFT.h>
Public Member Functions | |
FFT () | |
Default constructor. | |
virtual | ~FFT () |
Destructor. | |
void | setup (IntVec< D > const &meshDimensions) |
Setup grid dimensions, FFT plans and work space. | |
void | forwardTransform (RField< D > const &in, RFieldDft< D > &out) const |
Compute forward (real-to-complex) discrete Fourier transform (DFT). | |
void | inverseTransformUnsafe (RFieldDft< D > &in, RField< D > &out) const |
Compute inverse (complex-to-real) DFT, overwriting the input. | |
void | inverseTransformSafe (RFieldDft< D > const &in, RField< D > &out) const |
Compute inverse (complex-to-real) DFT without overwriting input. | |
void | forwardTransform (CField< D > const &in, CField< D > &out) const |
Compute forward (complex-to-complex) discrete Fourier transform. | |
void | inverseTransform (CField< D > const &in, CField< D > &out) const |
Compute complex-to-complex inverse Fourier transform. | |
IntVec< D > const & | meshDimensions () const |
Return the dimensions of the grid for which this was setup. | |
bool | isSetup () const |
Has this FFT object been setup? |
Static Public Member Functions | |
static void | computeKMesh (IntVec< D > const &rMeshDimensions, IntVec< D > &kMeshDimensions, int &kSize) |
Compute dimensions and size of k-space mesh for DFT of real data. | |
static bool | hasImplicitInverse (IntVec< D > const &wavevector, IntVec< D > const &meshDimensions) |
Does this wavevector have implicit inverse in DFT or real data? |
Fourier transform wrapper.
This class is a wrapper for plan creation and discrete Fourier transform (DFT) functions provided by the FFTW library, providing an interface to the field container classes RField<D>, RField<Dft>, and CField<D>.
Pscf::Prdc::Cuda::FFT< D >::FFT | ( | ) |
Default constructor.
Definition at line 57 of file cpu/FFT.tpp.
|
virtual |
Destructor.
Definition at line 72 of file cpu/FFT.tpp.
void Pscf::Prdc::Cuda::FFT< D >::setup | ( | IntVec< D > const & | meshDimensions | ) |
Setup grid dimensions, FFT plans and work space.
meshDimensions | Dimensions of real-space grid. |
Definition at line 92 of file cpu/FFT.tpp.
References Pscf::Prdc::Cpu::CField< D >::allocate(), Pscf::Prdc::Cpu::RField< D >::allocate(), Util::Array< Data >::capacity(), Pscf::Prdc::Cpu::CField< D >::meshDimensions(), meshDimensions(), Pscf::Prdc::Cpu::RField< D >::meshDimensions(), and UTIL_CHECK.
void Pscf::Prdc::Cuda::FFT< D >::forwardTransform | ( | RField< D > const & | in, |
RFieldDft< D > & | out ) const |
Compute forward (real-to-complex) discrete Fourier transform (DFT).
The resulting complex array is the output of FFTW's forward transform method scaled by a factor of 1/N (where N is the number of grid points). The scaling ensures that a round-trip Fourier transform (forward + inverse) reproduces the original input array.
This function does not overwrite or corrupt the input array.
in | array of real values on r-space grid |
out | array of complex values on k-space grid |
Definition at line 151 of file cpu/FFT.tpp.
References Util::Array< Data >::capacity(), Pscf::Prdc::Cpu::RField< D >::meshDimensions(), Pscf::Prdc::Cpu::RFieldDft< D >::meshDimensions(), and UTIL_CHECK.
void Pscf::Prdc::Cuda::FFT< D >::inverseTransformUnsafe | ( | RFieldDft< D > & | in, |
RField< D > & | out ) const |
Compute inverse (complex-to-real) DFT, overwriting the input.
The resulting real array is the unaltered output of FFTW's inverse transform function.
NOTE: The inverse transform generally overwrites and corrupts its input. This is the behavior of the complex-to-real (c2r) transform of the underlying FFTW library. See Sec. 2.3 of the FFTW 3 documentation at https://www.fftw.org/fftw3_doc, One-Dimensional DFTs of Real Data: "...the inverse transform (complex to real) has the side-effect of overwriting its input array, ..."
in | array of complex values on k-space grid (overwritten) |
out | array of real values on r-space grid |
Definition at line 179 of file cpu/FFT.tpp.
References Util::Array< Data >::capacity(), Pscf::Prdc::Cpu::RField< D >::meshDimensions(), Pscf::Prdc::Cpu::RFieldDft< D >::meshDimensions(), and UTIL_CHECK.
Referenced by inverseTransformSafe().
void Pscf::Prdc::Cuda::FFT< D >::inverseTransformSafe | ( | RFieldDft< D > const & | in, |
RField< D > & | out ) const |
Compute inverse (complex-to-real) DFT without overwriting input.
The resulting real array is the unaltered output of FFTW's inverse transform function.
This function makes a copy of the input data and passes the copy to inverseTransformUnsafe to avoid overwriting the input data.
in | array of complex values on k-space grid |
out | array of real values on r-space grid |
Definition at line 197 of file cpu/FFT.tpp.
References Util::Array< Data >::capacity(), inverseTransformUnsafe(), and UTIL_CHECK.
void Pscf::Prdc::Cuda::FFT< D >::forwardTransform | ( | CField< D > const & | in, |
CField< D > & | out ) const |
Compute forward (complex-to-complex) discrete Fourier transform.
The resulting complex array is the output of FFTW's forward transform method scaled by a factor of 1/N (where N is the number of grid points). The scaling ensures that a round-trip Fourier transform (forward + inverse) reproduces the original input array.
This function does not overwrite or corrupt the input array.
Note that, in this transform, the k-space grid is the same size as the r-space grid, which differs from the real-to-complex transform.
in | array of complex values on r-space grid |
out | array of complex values on k-space grid |
Definition at line 214 of file cpu/FFT.tpp.
References Util::Array< Data >::capacity(), Pscf::Prdc::Cpu::CField< D >::meshDimensions(), and UTIL_CHECK.
void Pscf::Prdc::Cuda::FFT< D >::inverseTransform | ( | CField< D > const & | in, |
CField< D > & | out ) const |
Compute complex-to-complex inverse Fourier transform.
The resulting complex array is the unaltered output of FFTW's inverse transform function.
Note that, in this transform, the k-space grid is the same size as the r-space grid, which differs from the complex-to-real transform.
in | array of complex values on k-space grid |
out | array of complex values on r-space grid |
Definition at line 241 of file cpu/FFT.tpp.
References Util::Array< Data >::capacity(), Pscf::Prdc::Cpu::CField< D >::meshDimensions(), and UTIL_CHECK.
|
inline |
Return the dimensions of the grid for which this was setup.
Definition at line 263 of file cpu/FFT.h.
Referenced by hasImplicitInverse(), and setup().
|
inline |
|
static |
Compute dimensions and size of k-space mesh for DFT of real data.
A corresponding function is not needed for complex-to-complex transforms because the real-space and Fourier-space grids for complex transforms have the same dimensions.
rMeshDimensions | dimensions of real space grid (real data) |
kMeshDimensions | dimensions of k-space grid (complex data) |
kSize | number of point in k-space grid |
Definition at line 262 of file cpu/FFT.tpp.
Referenced by Pscf::Prdc::Cpu::RFieldDft< D >::allocate(), Pscf::Prdc::Cpu::WaveList< D >::allocate(), Pscf::Rpc::Block< D >::allocate(), Pscf::Rpg::Block< D >::allocate(), Pscf::Rpc::IntraCorrelation< D >::computeIntraCorrelations(), Pscf::Rpg::IntraCorrelation< D >::computeIntraCorrelations(), Pscf::Rpc::BinaryStructureFactorGrid< D >::setup(), Pscf::Rpc::FourthOrderParameter< D >::setup(), Pscf::Rpc::MaxOrderParameter< D >::setup(), Pscf::Rpg::BinaryStructureFactorGrid< D >::setup(), Pscf::Rpg::FourthOrderParameter< D >::setup(), Pscf::Rpg::LrAmCompressor< D >::setup(), Pscf::Rpg::LrCompressor< D >::setup(), and Pscf::Rpg::MaxOrderParameter< D >::setup().
|
inlinestatic |
Does this wavevector have implicit inverse in DFT or real data?
The inverse of a wavevector within the DFT mesh used for a DFT of real data is "implicit" if its inverse is not equivalent to any wavevector in this mesh. The inverse is not implicit if its inverse is equivalent to a wavevector within this mesh.
To compute a Fourier sum using the DFT of real r-space data, for any summand that has inversion symmetry (i.e., the same value for every vector and its inverse), sum over all vectors in the k-space mesh and multiply each element for each wavevector that has an implicit inverse by 2.0, and all others by 1.0.
wavevector | integer indices of wavevector |
meshDimensions | dimensions of real-space mesh |
Definition at line 271 of file cpu/FFT.h.
References meshDimensions().
Referenced by Pscf::Prdc::Cpu::WaveList< D >::allocate().