Introduction
Wrappers Reference
1. Wrappers for Using Plans
2. Basic Interface Wrappers
3. Advanced Interface Wrappers
4. Guru Interface Wrappers
5. Wisdom Wrappers
6. Memory Allocation
Parallel Mode
Calling Wrappers from Fortran
Installation
Creating a Wrapper Library
Application Assembling
Running Examples
Technical Support
This document describes a collection of C routines – wrappers that allow the FFTW interface to call the Intel® Math Kernel Library (Intel® MKL) discrete Fourier transform interface (DFTI). These wrappers correspond to the FFTW version 3.x and the Intel MKL versions 7.0 and later.
The purpose of this set of wrappers is to enable developers whose programs
currently use FFTW to achieve the performance of the Intel MKL Fourier
transforms without changing the program source code. The only change that is
required is to modify the header file fftw3.h
(see
Creating a Wrapper Library).
Because of differences between FFTW and Intel MKL DFTI functionalities, there
are a lot of restrictions on using wrappers instead of FTTW functions. However,
many typical DFTs can be computed using these wrappers.
Please refer to the Intel MKL DFTI documentation for better understanding the effects from the use of the wrappers.
Additional wrappers may be added in the future to extend FFTW functionality available with Intel MKL.
Each FFTW function has its own wrapper. Some of them are empty and do nothing, but they are still needed to avoid compilation errors and satisfy the function calls.
Note that Intel MKL DFTI operates on float and double-precision data types and does not support the long-double data type used by the FFTW functions.
void fftw_execute(const fftw_plan plan);
void fftw_destroy_plan(const fftw_plan plan);
void fftwf_execute(const fftw_plan plan);
void fftwf_destroy_plan(const fftw_plan plan);
Wrappers for execution and plan destruction functions are listed in Wrappers for Using Plans.
fftw_plan fftw_plan_dft_1d(int n, fftw_complex
*in, fftw_complex *out, int sign, unsigned flags);
fftw_plan fftw_plan_dft_2d(int nx, int ny,
fftw_complex *in, fftw_complex *out, int sign, unsigned flags);
fftw_plan fftw_plan_dft_3d(int nx, int ny, int
nz, fftw_complex *in, fftw_complex *out, int sign, unsigned flags);
fftw_plan fftw_plan_dft(int rank, const int
*n, fftw_complex *in, fftw_complex *out, int sign, unsigned flags);
fftwf_plan fftwf_plan_dft_1d(int n,
fftwf_complex *in, fftwf_complex *out, int sign, unsigned flags);
fftwf_plan fftwf_plan_dft_2d(int nx, int ny,
fftwf_complex *in, fftwf_complex *out, int sign, unsigned flags);
fftwf_plan fftwf_plan_dft_3d(int nx, int ny,
int nz, fftwf_complex *in, fftwf_complex *out, int sign, unsigned flags);
fftwf_plan fftwf_plan_dft(int rank, const int
*n, fftwf_complex *in, fftwf_complex *out, int sign, unsigned flags);
Argument restrictions. The
same algorithm corresponds to all values of the flags
parameter.
fftw_plan fftw_plan_dft_r2c(int rank, const
int *n, double *in, fftw_complex *out, unsigned flags);
fftw_plan fftw_plan_dft_r2c_1d(int n, double
*in, fftw_complex *out, unsigned flags);
fftw_plan fftw_plan_dft_r2c_2d(int nx, int ny,
double *in, fftw_complex *out, unsigned flags);
fftw_plan fftw_plan_dft_r2c_3d(int nx, int ny,
int nz, double *in, fftw_complex *out, unsigned flags);
fftw_plan fftw_plan_dft_c2r(int rank, const
int *n, fftw_complex *in, double *out, unsigned flags);
fftw_plan fftw_plan_dft_c2r_1d(int n,
fftw_complex *in, double *out, unsigned flags);
fftw_plan fftw_plan_dft_c2r_2d(int nx, int ny,
fftw_complex *in, double *out, unsigned flags);
fftw_plan fftw_plan_dft_c2r_3d(int nx, int ny,
int nz, fftw_complex *in, double *out, unsigned flags);
fftwf_plan fftwf_plan_dft_r2c(int rank, const
int *n, float *in, fftwf_complex *out, unsigned flags);
fftwf_plan fftwf_plan_dft_r2c_1d(int n, float
*in, fftwf_complex *out, unsigned flags);
fftwf_plan fftwf_plan_dft_r2c_2d(int nx, int
ny, float *in, fftwf_complex *out, unsigned flags);
fftwf_plan fftwf_plan_dft_r2c_3d(int nx, int
ny, int nz, float *in, fftwf_complex *out, unsigned flags);
fftwf_plan fftwf_plan_dft_c2r(int rank, const
int *n, fftwf_complex *in, float *out, unsigned flags);
fftwf_plan fftwf_plan_dft_c2r_1d(int n,
fftwf_complex *in, float *out, unsigned flags);
fftwf_plan fftwf_plan_dft_c2r_2d(int nx, int
ny, fftwf_complex *in, float *out, unsigned flags);
fftwf_plan fftwf_plan_dft_c2r_3d(int nx, int
ny, int nz, fftwf_complex *in, float *out, unsigned flags);
Argument restrictions. The
same algorithm corresponds to all values of the flags
parameter.
In-place
transforms for 2D and 3D DFT are not supported.
All wrappers are empty and do nothing, as Intel MKL DFTI does not currently support this functionality.
Wrappers for execution and plan destruction functions are listed in Wrappers for Using Plans.
fftw_plan fftw_plan_many_dft(int rank, const
int *n, int howmany, fftw_complex *in, const int *inembed, int istride, int
idist, fftw_complex *out, const int *onembed, int ostride, int odist, int sign,
unsigned flags);
fftwf_plan fftwf_plan_many_dft(int rank, const
int *n, int howmany, fftwf_complex *in, const int *inembed, int istride, int
idist, fftwf_complex *out, const int *onembed, int ostride, int odist, int sign,
unsigned flags);
Argument restrictions. The
same algorithm corresponds to all values of the flags
parameter.
fftw_plan fftw_plan_many_dft_r2c(int rank,
const int *n, int howmany, double* in, const int *inembed, int istride, int
idist, fftw_complex *out, const int *onembed, int ostride, int odist, unsigned
flags);
fftwf_plan fftwf_plan_many_dft_r2c(int rank,
const int *n, int howmany, float* in, const int *inembed, int istride, int
idist, fftwf_complex *out, const int *onembed, int ostride, int odist, unsigned
flags);
fftw_plan fftw_plan_many_dft_c2r(int rank,
const int *n, int howmany, fftw_complex * in, const int *inembed, int istride,
int idist, double *out, const int *onembed, int ostride, int odist, unsigned
flags);
fftwf_plan fftwf_plan_many_dft_c2r(int rank,
const int *n, int howmany, fftwf_complex* in, const int *inembed, int istride,
int idist, float *out, const int *onembed, int ostride, int odist, unsigned
flags);
All wrappers are empty and do nothing, as the Intel MKL DFTI does not currently support this functionality.
fftw_plan fftw_plan_guru_dft(int rank, const
fftw_iodim *dims, int howmany_rank, const fftw_iodim *howmany_dims, fftw_complex
*in, fftw_complex *out, int sign, unsigned flags);
fftwf_plan fftwf_plan_guru_dft(int rank, const
fftwf_iodim *dims, int howmany_rank, const fftwf_iodim *howmany_dims,
fftwf_complex *in, fftwf_complex *out, int sign, unsigned flags);
flags
parameter.
The rest of the wrappers are empty and do nothing, as the Intel MKL DFTI currently does not support split arrays.
All wrappers are empty and do nothing.
Real-data wrappers (without support of split arrays) will be added in later versions of the Intel MKL.
All wrappers are empty and do nothing, as the Intel MKL DFTI currently does not support these functionalities.
void fftw_execute_dft(const fftw_plan p,
fftw_complex *in, fftw_complex *out);
void fftw_execute_dft_r2c(const fftw_plan p,
double *in, fftw_complex *out);
void fftw_execute_dft_c2r(const fftw_plan p,
fftw_complex *in, double *out);
void fftwf_execute_dft(const fftwf_plan p,
fftwf_complex *in, fftwf_complex *out);
void fftwf_execute_dft_r2c(const fftwf_plan p,
float *in, fftwf_complex *out);
void fftwf_execute_dft_c2r(const fftwf_plan p,
fftwf_complex *in, float *out);
The rest of the wrappers are empty and do nothing.
Real-data wrappers (without support of split arrays) will be added in later versions of the Intel MKL.
Wrappers for more execution and plan destruction functions are listed in Wrappers for Using Plans.
All wrappers are empty and do nothing, as the Intel MKL DFTI currently does not support these functionalities.
All wrappers in this
section are empty and do nothing. The Intel MKL allocation function cannot
align the allocatable array. It is recommended to allocate memory using
malloc
or similar functions and
align your array yourself. To do that, it is necessary to allocate extra memory
and shift the array address for the DFT data. See also the
Performance section in the Intel MKL documentation file
mkluse.htm
All wrappers in this section are empty and do nothing, as the Intel MKL DFTI implements a different mechanism of parallelization. If you want to use Intel MKL DFTI routines in parallel mode or call wrappers from a multi-threaded application, please refer to the Intel MKL documentation to learn how to manage the number of threads.
There are no wrappers in this section yet.
Wrappers are delivered as the source code, which must be compiled by a user to build the wrapper library. Then the FFTW library can be substituted by the wrapper and Intel MKL libraries. The source code for the wrappers and makefiles with the wrapper list files are located in the \interfaces\fftw3xc sub-directory in the Intel MKL directory for C wrappers. Both __release_lnx and __release_win directories have the "interfaces" sub-directory.
Two header files are used
to compile the wrapper library:
fftw3_mkl.h
and
fftw3.h
.
The
fftw3_mkl.h
file is located in the \interfaces\fftw3xc\wrappers sub-directory in the Intel
MKL directory. The original FFTW (www.fftw.org) header file
fftw3.h
is slightly
modified (all rows containing calls to the fftw3.lib
are
commented) and placed in the \include\fftw sub-directory in the Intel MKL directory.
Makefiles contain the following parameters: platform (required), compiler, and function. Description of these parameters can be found in the makefile comment heading.
Examples
The commandmake lib64
make lib32 F=intel7
In the wrapper library names, the suffix corresponds to the used compiler and the underscore is preceded with letter
"с" (meaning C wrappers).
For example,
fftw3xc_intel.lib
(libfftw2xc_intel.a
for Linux*)
fftw3xc_ms.lib
(libfftw2xc_gnu.a
for Linux).
The adapted
fftw3.h
header file (see above)
should be used when you build applications.
There are some examples that demonstrate how to use the wrapper library. The source code for the examples, makefiles used to run them, and the example list files are located in the \examples\fftw3xc sub-directory in the Intel MKL directory.
Example makefile parameters are the same as wrapper library makefile parameters. Example makefiles run examples. However, if the appropriate wrapper library is not yet created, the makefile will first build it in the same way as the wrapper library makefile does and then proceed to examples.
If the parameter
function=<example_name>
is defined, then only the specified example will run.
Otherwise, all examples from the \examples\fftw3xc\source sub-directory will run. The
sub-directory \_results will be created, and the results will be stored there in
the example_name.res
files.
Intel, the Intel logo, Intel SpeedStep, Intel
NetBurst, Intel NetStructure, MMX, Intel386, Intel486, Celeron, Intel Centrino,
Intel Xeon, Intel XScale, Itanium, Pentium, Pentium II Xeon, Pentium III Xeon, Pentium M, and VTune are trademarks or registered
trademarks of Intel Corporation or its subsidiaries in the United States and
other countries.
* Other names and brands may be claimed as the property
of others.
Copyright © 2005 - 2006, Intel Corporation.