Gromacs  2018.8
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Typedefs | Enumerations | Functions | Variables
fft.h File Reference
#include <stdio.h>
#include "gromacs/math/gmxcomplex.h"
#include "gromacs/utility/real.h"
+ Include dependency graph for fft.h:
+ This graph shows which files directly or indirectly include this file:

Description

Fast Fourier Transforms.

This file provides an abstract Gromacs interface to Fourier transforms, including multi-dimensional and real-to-complex transforms.

Internally it is implemented as wrappers to external libraries such as FFTW or the Intel Math Kernel Library, but we also have a built-in version of FFTPACK in case the faster alternatives are unavailable.

We also provide our own multi-dimensional transform setups even when the underlying library does not support it directly.

Typedefs

typedef struct gmx_fft * gmx_fft_t
 Datatype for FFT setup. More...
 
typedef int gmx_fft_flag
 Specifier for FFT flags. More...
 

Enumerations

enum  gmx_fft_direction { GMX_FFT_FORWARD, GMX_FFT_BACKWARD, GMX_FFT_REAL_TO_COMPLEX, GMX_FFT_COMPLEX_TO_REAL }
 Specifier for FFT direction. More...
 

Functions

int gmx_fft_init_1d (gmx_fft_t *fft, int nx, gmx_fft_flag flags)
 Setup a 1-dimensional complex-to-complex transform. More...
 
int gmx_fft_init_many_1d (gmx_fft_t *fft, int nx, int howmany, gmx_fft_flag flags)
 Setup multiple 1-dimensional complex-to-complex transform. More...
 
int gmx_fft_init_1d_real (gmx_fft_t *fft, int nx, gmx_fft_flag flags)
 Setup a 1-dimensional real-to-complex transform. More...
 
int gmx_fft_init_many_1d_real (gmx_fft_t *fft, int nx, int howmany, gmx_fft_flag flags)
 Setup multiple 1-dimensional real-to-complex transform. More...
 
int gmx_fft_init_2d_real (gmx_fft_t *fft, int nx, int ny, gmx_fft_flag flags)
 Setup a 2-dimensional real-to-complex transform. More...
 
int gmx_fft_1d (gmx_fft_t setup, enum gmx_fft_direction dir, void *in_data, void *out_data)
 Perform a 1-dimensional complex-to-complex transform. More...
 
int gmx_fft_many_1d (gmx_fft_t setup, enum gmx_fft_direction dir, void *in_data, void *out_data)
 Perform many 1-dimensional complex-to-complex transforms. More...
 
int gmx_fft_1d_real (gmx_fft_t setup, enum gmx_fft_direction dir, void *in_data, void *out_data)
 Perform a 1-dimensional real-to-complex transform. More...
 
int gmx_fft_many_1d_real (gmx_fft_t setup, enum gmx_fft_direction dir, void *in_data, void *out_data)
 Perform many 1-dimensional real-to-complex transforms. More...
 
int gmx_fft_2d_real (gmx_fft_t setup, enum gmx_fft_direction dir, void *in_data, void *out_data)
 Perform a 2-dimensional real-to-complex transform. More...
 
void gmx_fft_destroy (gmx_fft_t setup)
 Release an FFT setup structure. More...
 
void gmx_many_fft_destroy (gmx_fft_t setup)
 Release a many FFT setup structure. More...
 
int gmx_fft_transpose_2d (t_complex *in_data, t_complex *out_data, int nx, int ny)
 Transpose 2d complex matrix, in-place or out-of-place. More...
 
void gmx_fft_cleanup ()
 Cleanup global data of FFT. More...
 

Variables

static const int GMX_FFT_FLAG_NONE = 0
 Macro to indicate no special flags for FFT routines. More...
 
static const int GMX_FFT_FLAG_CONSERVATIVE = (1<<0)
 Flag to disable FFT optimizations based on timings, see gmx_fft_flag. More...
 

Typedef Documentation

typedef int gmx_fft_flag

Specifier for FFT flags.

Some FFT libraries (FFTW, in particular) can do timings and other tricks to try and optimize the FFT for the current architecture. However, this can also lead to results that differ between consecutive runs with identical input. To avoid this, the conservative flag will attempt to disable such optimization, but there are no guarantees since we cannot control what the FFT libraries do internally.

typedef struct gmx_fft* gmx_fft_t

Datatype for FFT setup.

The gmx_fft_t type contains all the setup information, e.g. twiddle factors, necessary to perform an FFT. Internally it is mapped to whatever FFT library we are using, or the built-in FFTPACK if no fast external library is available.

Since some of the libraries (e.g. MKL) store work array data in their handles this datatype should only be used for one thread at a time, i.e. they should allocate one instance each when executing in parallel.

Enumeration Type Documentation

Specifier for FFT direction.

The definition of the 1D forward transform from input x[] to output y[] is

\[ y_{k} = \sum_{j=0}^{N-1} x_{j} \exp{-i 2 \pi j k /N} \]

while the corresponding backward transform is

\[ y_{k} = \sum_{j=0}^{N-1} x_{j} \exp{i 2 \pi j k /N} \]

A forward-backward transform pair will this result in data scaled by N.

For complex-to-complex transforms you can only use one of GMX_FFT_FORWARD or GMX_FFT_BACKWARD, and for real-complex transforms you can only use GMX_FFT_REAL_TO_COMPLEX or GMX_FFT_COMPLEX_TO_REAL.

Enumerator
GMX_FFT_FORWARD 

Forward complex-to-complex transform.

GMX_FFT_BACKWARD 

Backward complex-to-complex transform.

GMX_FFT_REAL_TO_COMPLEX 

Real-to-complex valued FFT.

GMX_FFT_COMPLEX_TO_REAL 

Complex-to-real valued FFT.

Function Documentation

int gmx_fft_1d ( gmx_fft_t  setup,
enum gmx_fft_direction  dir,
void *  in_data,
void *  out_data 
)

Perform a 1-dimensional complex-to-complex transform.

Performs an instance of a transform previously initiated.

Parameters
setupSetup returned from gmx_fft_init_1d()
dirForward or Backward
in_dataInput grid data. This should be allocated with gmx_new() to make it 16-byte aligned for better performance.
out_dataOutput grid data. This should be allocated with gmx_new() to make it 16-byte aligned for better performance. You can provide the same pointer for in_data and out_data to perform an in-place transform.
Returns
0 on success, or an error code.
Note
Data pointers are declared as void, to avoid casting pointers depending on your grid type.
int gmx_fft_1d_real ( gmx_fft_t  setup,
enum gmx_fft_direction  dir,
void *  in_data,
void *  out_data 
)

Perform a 1-dimensional real-to-complex transform.

Performs an instance of a transform previously initiated.

Parameters
setupSetup returned from gmx_fft_init_1d_real()
dirReal-to-complex or complex-to-real
in_dataInput grid data. This should be allocated with gmx_new() to make it 16-byte aligned for better performance.
out_dataOutput grid data. This should be allocated with gmx_new() to make it 16-byte aligned for better performance. You can provide the same pointer for in_data and out_data to perform an in-place transform.

If you are doing an in-place transform, the array must be padded up to an even integer length so n/2 complex numbers can fit. Out-of-place arrays should not be padded (although it doesn't matter in 1d).

Returns
0 on success, or an error code.
Note
Data pointers are declared as void, to avoid casting pointers depending on transform direction.
int gmx_fft_2d_real ( gmx_fft_t  setup,
enum gmx_fft_direction  dir,
void *  in_data,
void *  out_data 
)

Perform a 2-dimensional real-to-complex transform.

Performs an instance of a transform previously initiated.

Parameters
setupSetup returned from gmx_fft_init_1d_real()
dirReal-to-complex or complex-to-real
in_dataInput grid data. This should be allocated with gmx_new() to make it 16-byte aligned for better performance.
out_dataOutput grid data. This should be allocated with gmx_new() to make it 16-byte aligned for better performance. You can provide the same pointer for in_data and out_data to perform an in-place transform.
Returns
0 on success, or an error code.
Note
If you are doing an in-place transform, the last dimension of the array MUST be padded up to an even integer length so n/2 complex numbers can fit. Thus, if the real grid e.g. has dimension 5*3, you must allocate it as a 5*4 array, where the last element in the second dimension is padding. The complex data will be written to the same array, but since that dimension is 5*2 it will now fill the entire array. Reverse complex-to-real in-place transformation will produce the same sort of padded array.

The padding does NOT apply to out-of-place transformation. In that case the input array will simply be 5*3 of real, while the output is 5*2 of complex.

Note
Data pointers are declared as void, to avoid casting pointers depending on transform direction.
void gmx_fft_cleanup ( )

Cleanup global data of FFT.

Any plans are invalid after this function. Should be called after all plans have been destroyed.

void gmx_fft_destroy ( gmx_fft_t  setup)

Release an FFT setup structure.

Destroy setup and release all allocated memory.

Parameters
setupSetup returned from gmx_fft_init_1d(), or one of the other initializers.
int gmx_fft_init_1d ( gmx_fft_t fft,
int  nx,
gmx_fft_flag  flags 
)

Setup a 1-dimensional complex-to-complex transform.

Parameters
fftPointer to opaque Gromacs FFT datatype
nxLength of transform
flagsFFT options
Returns
status - 0 or a standard error message.
Note
Since some of the libraries (e.g. MKL) store work array data in their handles this datatype should only be used for one thread at a time, i.e. you should create one copy per thread when executing in parallel.
int gmx_fft_init_1d_real ( gmx_fft_t fft,
int  nx,
gmx_fft_flag  flags 
)

Setup a 1-dimensional real-to-complex transform.

Parameters
fftPointer to opaque Gromacs FFT datatype
nxLength of transform in real space
flagsFFT options
Returns
status - 0 or a standard error message.
Note
Since some of the libraries (e.g. MKL) store work array data in their handles this datatype should only be used for one thread at a time, i.e. you should create one copy per thread when executing in parallel.
int gmx_fft_init_2d_real ( gmx_fft_t fft,
int  nx,
int  ny,
gmx_fft_flag  flags 
)

Setup a 2-dimensional real-to-complex transform.

Parameters
fftPointer to opaque Gromacs FFT datatype
nxLength of transform in first dimension
nyLength of transform in second dimension
flagsFFT options

The normal space is assumed to be real, while the values in frequency space are complex.

Returns
status - 0 or a standard error message.
Note
Since some of the libraries (e.g. MKL) store work array data in their handles this datatype should only be used for one thread at a time, i.e. you should create one copy per thread when executing in parallel.
int gmx_fft_init_many_1d ( gmx_fft_t fft,
int  nx,
int  howmany,
gmx_fft_flag  flags 
)

Setup multiple 1-dimensional complex-to-complex transform.

Parameters
fftPointer to opaque Gromacs FFT datatype
nxLength of transform
howmanyHowmany 1D FFT
flagsFFT options
Returns
status - 0 or a standard error message.
Note
Since some of the libraries (e.g. MKL) store work array data in their handles this datatype should only be used for one thread at a time, i.e. you should create one copy per thread when executing in parallel.
int gmx_fft_init_many_1d_real ( gmx_fft_t fft,
int  nx,
int  howmany,
gmx_fft_flag  flags 
)

Setup multiple 1-dimensional real-to-complex transform.

Parameters
fftPointer to opaque Gromacs FFT datatype
nxLength of transform in real space
howmanyHomany 1D FFTs
flagsFFT options
Returns
status - 0 or a standard error message.
Note
Since some of the libraries (e.g. MKL) store work array data in their handles this datatype should only be used for one thread at a time, i.e. you should create one copy per thread when executing in parallel.
int gmx_fft_many_1d ( gmx_fft_t  setup,
enum gmx_fft_direction  dir,
void *  in_data,
void *  out_data 
)

Perform many 1-dimensional complex-to-complex transforms.

Performs many instances of a transform previously initiated.

Parameters
setupSetup returned from gmx_fft_init_1d()
dirForward or Backward
in_dataInput grid data. This should be allocated with gmx_new() to make it 16-byte aligned for better performance.
out_dataOutput grid data. This should be allocated with gmx_new() to make it 16-byte aligned for better performance. You can provide the same pointer for in_data and out_data to perform an in-place transform.
Returns
0 on success, or an error code.
Note
Data pointers are declared as void, to avoid casting pointers depending on your grid type.
int gmx_fft_many_1d_real ( gmx_fft_t  setup,
enum gmx_fft_direction  dir,
void *  in_data,
void *  out_data 
)

Perform many 1-dimensional real-to-complex transforms.

Performs many instances of a transform previously initiated.

Parameters
setupSetup returned from gmx_fft_init_1d_real()
dirReal-to-complex or complex-to-real
in_dataInput grid data. This should be allocated with gmx_new() to make it 16-byte aligned for better performance.
out_dataOutput grid data. This should be allocated with gmx_new() to make it 16-byte aligned for better performance. You can provide the same pointer for in_data and out_data to perform an in-place transform.

If you are doing an in-place transform, the array must be padded up to an even integer length so n/2 complex numbers can fit. Out-of-place arrays should not be padded (although it doesn't matter in 1d).

Returns
0 on success, or an error code.
Note
Data pointers are declared as void, to avoid casting pointers depending on transform direction.
int gmx_fft_transpose_2d ( t_complex *  in_data,
t_complex *  out_data,
int  nx,
int  ny 
)

Transpose 2d complex matrix, in-place or out-of-place.

This routines works when the matrix is non-square, i.e. nx!=ny too, without allocating an entire matrix of work memory, which is important for huge FFT grids.

Parameters
in_dataInput data, to be transposed
out_dataOutput, transposed data. If this is identical to in_data, an in-place transpose is performed.
nxNumber of rows before transpose
nyNumber of columns before transpose
Returns
GMX_SUCCESS, or an error code from gmx_errno.h
void gmx_many_fft_destroy ( gmx_fft_t  setup)

Release a many FFT setup structure.

Destroy setup and release all allocated memory.

Parameters
setupSetup returned from gmx_fft_init_1d(), or one of the other initializers.

Variable Documentation

const int GMX_FFT_FLAG_CONSERVATIVE = (1<<0)
static

Flag to disable FFT optimizations based on timings, see gmx_fft_flag.

const int GMX_FFT_FLAG_NONE = 0
static

Macro to indicate no special flags for FFT routines.