PSCF v1.3
Pscf::Prdc::Cuda::FFT< D > Class Template Reference

Fourier transform wrapper for real or complex data. More...

#include <FFT.h>

Public Member Functions

 FFT ()
 Default constructor.
virtual ~FFT ()
 Destructor.
void setup (IntVec< D > const &meshDimensions)
 Setup grid dimensions, plans and work space.
void forwardTransform (RField< D > const &rField, RFieldDft< D > &kField) const
 Compute forward (real-to-complex) discrete Fourier transform.
void inverseTransformUnsafe (RFieldDft< D > &kField, RField< D > &rField) const
 Compute inverse (complex-to-real) DFT, overwriting the input.
void inverseTransformSafe (RFieldDft< D > const &kField, RField< D > &rField) const
 Compute inverse (complex-to-real) DFT without overwriting input.
void forwardTransform (CField< D > const &rField, CField< D > &kField) const
 Compute forward (complex-to-complex) discrete Fourier transform.
void inverseTransform (CField< D > const &kField, CField< D > &rField) const
 Compute inverse (complex-to-complex) discrete Fourier transform.
const IntVec< D > & meshDimensions () const
 Return the dimensions of the grid for which this was allocated.
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 a wavevector have an implicit inverse in DFT for real data?

Detailed Description

template<int D>
class Pscf::Prdc::Cuda::FFT< D >

Fourier transform wrapper for real or complex data.

This class is a wrapper for plan creation and discrete Fourier transform (DFT) functions provided by the NVIDIA cufft library, providing an interface to the field container classes RField<D>, RField<Dft>, and CField<D> in namespace Pscf::Prdc::Cuda.

Definition at line 38 of file cuda/FFT.h.

Constructor & Destructor Documentation

◆ FFT()

template<int D>
Pscf::Prdc::Cuda::FFT< D >::FFT ( )

Default constructor.

◆ ~FFT()

template<int D>
virtual Pscf::Prdc::Cuda::FFT< D >::~FFT ( )
virtual

Destructor.

Member Function Documentation

◆ setup()

template<int D>
void Pscf::Prdc::Cuda::FFT< D >::setup ( IntVec< D > const & meshDimensions)

Setup grid dimensions, plans and work space.

Parameters
meshDimensionsDimensions of real-space grid.

References meshDimensions().

◆ forwardTransform() [1/2]

template<int D>
void Pscf::Prdc::Cuda::FFT< D >::forwardTransform ( RField< D > const & rField,
RFieldDft< D > & kField ) const

Compute forward (real-to-complex) discrete Fourier transform.

The resulting complex array is the output of cuFFT'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 transform is inherently "safe", meaning that the input array will not be modified during the DFT.

Parameters
rFieldreal values on r-space grid (input, gpu mem)
kFieldcomplex values on k-space grid (output, gpu mem)

◆ inverseTransformUnsafe()

template<int D>
void Pscf::Prdc::Cuda::FFT< D >::inverseTransformUnsafe ( RFieldDft< D > & kField,
RField< D > & rField ) const

Compute inverse (complex-to-real) DFT, overwriting the input.

The resulting real array is the unaltered output of cuFFT's inverse transform function.

This transform is inherently "unsafe", meaning that the input array will be overwritten during the DFT. Accordingly, the input array kField is not const like it is in all the other transform methods.

Parameters
kFieldcomplex values on k-space grid (input, gpu mem)
rFieldreal values on r-space grid (output, gpu mem)

◆ inverseTransformSafe()

template<int D>
void Pscf::Prdc::Cuda::FFT< D >::inverseTransformSafe ( RFieldDft< D > const & kField,
RField< D > & rField ) const

Compute inverse (complex-to-real) DFT without overwriting input.

This method makes a copy of the input kField array before performing the DFT, thus preserving the input.

Parameters
kFieldcomplex values on k-space grid (input, gpu mem)
rFieldreal values on r-space grid (output, gpu mem)

◆ forwardTransform() [2/2]

template<int D>
void Pscf::Prdc::Cuda::FFT< D >::forwardTransform ( CField< D > const & rField,
CField< D > & kField ) const

Compute forward (complex-to-complex) discrete Fourier transform.

The resulting complex array is the output of cuFFT'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 transform is inherently "safe", meaning that the input array will not be modified during the DFT.

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.

Parameters
rFieldcomplex values on r-space grid (input, gpu mem)
kFieldcomplex values on k-space grid (output, gpu mem)

◆ inverseTransform()

template<int D>
void Pscf::Prdc::Cuda::FFT< D >::inverseTransform ( CField< D > const & kField,
CField< D > & rField ) const

Compute inverse (complex-to-complex) discrete Fourier transform.

The resulting complex array is the unaltered output of cuFFT's inverse transform function.

This transform is inherently "safe", meaning that the input array will not be modified during the DFT.

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.

Parameters
kFieldcomplex values on k-space grid (input, gpu mem)
rFieldcomplex values on r-space grid (output, gpu mem)

◆ meshDimensions()

template<int D>
const IntVec< D > & Pscf::Prdc::Cuda::FFT< D >::meshDimensions ( ) const

Return the dimensions of the grid for which this was allocated.

Referenced by hasImplicitInverse(), and setup().

◆ isSetup()

template<int D>
bool Pscf::Prdc::Cuda::FFT< D >::isSetup ( ) const

Has this FFT object been setup?

◆ computeKMesh()

template<int D>
void Pscf::Prdc::Cuda::FFT< D >::computeKMesh ( IntVec< D > const & rMeshDimensions,
IntVec< D > & kMeshDimensions,
int & kSize )
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 have the same dimensions in this case.

Parameters
rMeshDimensionsdimensions of real space grid (real data)
kMeshDimensionsdimensions of k-space grid (complex data)
kSizenumber of point in k-space grid

Referenced by Pscf::Prdc::Cuda::RFieldDft< D >::associate().

◆ hasImplicitInverse()

template<int D>
bool Pscf::Prdc::Cuda::FFT< D >::hasImplicitInverse ( IntVec< D > const & wavevector,
IntVec< D > const & meshDimensions )
static

Does a wavevector have an implicit inverse in DFT for 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 k-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.

Parameters
wavevectorinteger indices of wavevector (in)
meshDimensionsdimensions of real-space mesh (in)
Returns
true iff inverse is implicit (i.e., not in DFT k-mesh)

References meshDimensions().


The documentation for this class was generated from the following file: