1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * Gromacs 4.0 Copyright (c) 1991-2003
5 * David van der Spoel, Erik Lindahl, University of Groningen.
7 * This program is free software; you can redistribute it and/or
8 * modify it under the terms of the GNU General Public License
9 * as published by the Free Software Foundation; either version 2
10 * of the License, or (at your option) any later version.
12 * To help us fund GROMACS development, we humbly ask that you cite
13 * the research papers on the package. Check out http://www.gromacs.org
16 * Gnomes, ROck Monsters And Chili Sauce
23 * \brief Fast Fourier Transforms.
25 * This file provides an abstract Gromacs interface to Fourier transforms,
26 * including multi-dimensional and real-to-complex transforms.
28 * Internally it is implemented as wrappers to external libraries such
29 * as FFTW or the Intel Math Kernel Library, but we also have a built-in
30 * version of FFTPACK in case the faster alternatives are unavailable.
32 * We also provide our own multi-dimensional transform setups even when
33 * the underlying library does not support it directly.
39 #include "types/simple.h"
40 #include "gmxcomplex.h"
47 } /* fixes auto-indentation problems */
52 /*! \brief Datatype for FFT setup
54 * The gmx_fft_t type contains all the setup information, e.g. twiddle
55 * factors, necessary to perform an FFT. Internally it is mapped to
56 * whatever FFT library we are using, or the built-in FFTPACK if no fast
57 * external library is available.
59 * Since some of the libraries (e.g. MKL) store work array data in their
60 * handles this datatype should only be used for one thread at a time, i.e.
61 * they should allocate one instance each when executing in parallel.
63 typedef struct gmx_fft
*
69 /*! \brief Specifier for FFT direction.
71 * The definition of the 1D forward transform from input x[] to output y[] is
73 * y_{k} = \sum_{j=0}^{N-1} x_{j} \exp{-i 2 \pi j k /N}
76 * while the corresponding backward transform is
79 * y_{k} = \sum_{j=0}^{N-1} x_{j} \exp{i 2 \pi j k /N}
82 * A forward-backward transform pair will this result in data scaled by N.
84 * For complex-to-complex transforms you can only use one of
85 * GMX_FFT_FORWARD or GMX_FFT_BACKWARD, and for real-complex transforms you
86 * can only use GMX_FFT_REAL_TO_COMPLEX or GMX_FFT_COMPLEX_TO_REAL.
88 typedef enum gmx_fft_direction
90 GMX_FFT_FORWARD
, /*!< Forward complex-to-complex transform */
91 GMX_FFT_BACKWARD
, /*!< Backward complex-to-complex transform */
92 GMX_FFT_REAL_TO_COMPLEX
, /*!< Real-to-complex valued fft */
93 GMX_FFT_COMPLEX_TO_REAL
/*!< Complex-to-real valued fft */
96 /*! \brief Specifier for FFT flags.
98 * Some FFT libraries (FFTW, in particular) can do timings and other
99 * tricks to try and optimize the FFT for the current architecture. However,
100 * this can also lead to results that differ between consecutive runs with
102 * To avoid this, the conservative flag will attempt to disable such
103 * optimization, but there are no guarantees since we cannot control what
104 * the FFT libraries do internally.
107 typedef int gmx_fft_flag
;
108 static const int GMX_FFT_FLAG_NONE
= 0;
109 static const int GMX_FFT_FLAG_CONSERVATIVE
= (1<<0);
111 /*! \brief Setup a 1-dimensional complex-to-complex transform
113 * \param fft Pointer to opaque Gromacs FFT datatype
114 * \param nx Length of transform
115 * \param flags FFT options
117 * \return status - 0 or a standard error message.
119 * \note Since some of the libraries (e.g. MKL) store work array data in their
120 * handles this datatype should only be used for one thread at a time,
121 * i.e. you should create one copy per thread when executing in parallel.
124 gmx_fft_init_1d (gmx_fft_t
* fft
,
130 /*! \brief Setup a 1-dimensional real-to-complex transform
132 * \param fft Pointer to opaque Gromacs FFT datatype
133 * \param nx Length of transform in real space
134 * \param flags FFT options
136 * \return status - 0 or a standard error message.
138 * \note Since some of the libraries (e.g. MKL) store work array data in their
139 * handles this datatype should only be used for one thread at a time,
140 * i.e. you should create one copy per thread when executing in parallel.
143 gmx_fft_init_1d_real (gmx_fft_t
* fft
,
149 /*! \brief Setup a 2-dimensional complex-to-complex transform
151 * \param fft Pointer to opaque Gromacs FFT datatype
152 * \param nx Length of transform in first dimension
153 * \param ny Length of transform in second dimension
154 * \param flags FFT options
156 * \return status - 0 or a standard error message.
158 * \note Since some of the libraries (e.g. MKL) store work array data in their
159 * handles this datatype should only be used for one thread at a time,
160 * i.e. you should create one copy per thread when executing in parallel.
163 gmx_fft_init_2d (gmx_fft_t
* fft
,
169 /*! \brief Setup a 2-dimensional real-to-complex transform
171 * \param fft Pointer to opaque Gromacs FFT datatype
172 * \param nx Length of transform in first dimension
173 * \param ny Length of transform in second dimension
174 * \param flags FFT options
176 * The normal space is assumed to be real, while the values in
177 * frequency space are complex.
179 * \return status - 0 or a standard error message.
181 * \note Since some of the libraries (e.g. MKL) store work array data in their
182 * handles this datatype should only be used for one thread at a time,
183 * i.e. you should create one copy per thread when executing in parallel.
186 gmx_fft_init_2d_real (gmx_fft_t
* fft
,
192 /*! \brief Setup a 3-dimensional complex-to-complex transform
194 * \param fft Pointer to opaque Gromacs FFT datatype
195 * \param nx Length of transform in first dimension
196 * \param ny Length of transform in second dimension
197 * \param nz Length of transform in third dimension
198 * \param flags FFT options
200 * \return status - 0 or a standard error message.
202 * \note Since some of the libraries (e.g. MKL) store work array data in their
203 * handles this datatype should only be used for one thread at a time,
204 * i.e. you should create one copy per thread when executing in parallel.
207 gmx_fft_init_3d (gmx_fft_t
* fft
,
214 /*! \brief Setup a 3-dimensional real-to-complex transform
216 * \param fft Pointer to opaque Gromacs FFT datatype
217 * \param nx Length of transform in first dimension
218 * \param ny Length of transform in second dimension
219 * \param nz Length of transform in third dimension
220 * \param flags FFT options
222 * The normal space is assumed to be real, while the values in
223 * frequency space are complex.
225 * \return status - 0 or a standard error message.
227 * \note Since some of the libraries (e.g. MKL) store work array data in their
228 * handles this datatype should only be used for one thread at a time,
229 * i.e. you should create one copy per thread when executing in parallel.
232 gmx_fft_init_3d_real (gmx_fft_t
* fft
,
240 /*! \brief Perform a 1-dimensional complex-to-complex transform
242 * Performs an instance of a transform previously initiated.
244 * \param setup Setup returned from gmx_fft_init_1d()
245 * \param dir Forward or Backward
246 * \param in_data Input grid data. This should be allocated with gmx_new()
247 * to make it 16-byte aligned for better performance.
248 * \param out_data Output grid data. This should be allocated with gmx_new()
249 * to make it 16-byte aligned for better performance.
250 * You can provide the same pointer for in_data and out_data
251 * to perform an in-place transform.
253 * \return 0 on success, or an error code.
255 * \note Data pointers are declared as void, to avoid casting pointers
256 * depending on your grid type.
259 gmx_fft_1d (gmx_fft_t setup
,
260 enum gmx_fft_direction dir
,
265 /*! \brief Perform a 1-dimensional real-to-complex transform
267 * Performs an instance of a transform previously initiated.
269 * \param setup Setup returned from gmx_fft_init_1d_real()
270 * \param dir Real-to-complex or complex-to-real
271 * \param in_data Input grid data. This should be allocated with gmx_new()
272 * to make it 16-byte aligned for better performance.
273 * \param out_data Output grid data. This should be allocated with gmx_new()
274 * to make it 16-byte aligned for better performance.
275 * You can provide the same pointer for in_data and out_data
276 * to perform an in-place transform.
278 * If you are doing an in-place transform, the array must be padded up to
279 * an even integer length so n/2 complex numbers can fit. Out-of-place arrays
280 * should not be padded (although it doesn't matter in 1d).
282 * \return 0 on success, or an error code.
284 * \note Data pointers are declared as void, to avoid casting pointers
285 * depending on transform direction.
288 gmx_fft_1d_real (gmx_fft_t setup
,
289 enum gmx_fft_direction dir
,
294 /*! \brief Perform a 2-dimensional complex-to-complex transform
296 * Performs an instance of a transform previously initiated.
298 * \param setup Setup returned from gmx_fft_init_1d()
299 * \param dir Forward or Backward
300 * \param in_data Input grid data. This should be allocated with gmx_new()
301 * to make it 16-byte aligned for better performance.
302 * \param out_data Output grid data. This should be allocated with gmx_new()
303 * to make it 16-byte aligned for better performance.
304 * You can provide the same pointer for in_data and out_data
305 * to perform an in-place transform.
307 * \return 0 on success, or an error code.
309 * \note Data pointers are declared as void, to avoid casting pointers
310 * depending on your grid type.
313 gmx_fft_2d (gmx_fft_t setup
,
314 enum gmx_fft_direction dir
,
319 /*! \brief Perform a 2-dimensional real-to-complex transform
321 * Performs an instance of a transform previously initiated.
323 * \param setup Setup returned from gmx_fft_init_1d_real()
324 * \param dir Real-to-complex or complex-to-real
325 * \param in_data Input grid data. This should be allocated with gmx_new()
326 * to make it 16-byte aligned for better performance.
327 * \param out_data Output grid data. This should be allocated with gmx_new()
328 * to make it 16-byte aligned for better performance.
329 * You can provide the same pointer for in_data and out_data
330 * to perform an in-place transform.
332 * \return 0 on success, or an error code.
334 * \note If you are doing an in-place transform, the last dimension of the
335 * array MUST be padded up to an even integer length so n/2 complex numbers can
336 * fit. Thus, if the real grid e.g. has dimension 5*3, you must allocate it as
337 * a 5*4 array, where the last element in the second dimension is padding.
338 * The complex data will be written to the same array, but since that dimension
339 * is 5*2 it will now fill the entire array. Reverse complex-to-real in-place
340 * transformation will produce the same sort of padded array.
342 * The padding does NOT apply to out-of-place transformation. In that case the
343 * input array will simply be 5*3 of real, while the output is 5*2 of complex.
345 * \note Data pointers are declared as void, to avoid casting pointers
346 * depending on transform direction.
349 gmx_fft_2d_real (gmx_fft_t setup
,
350 enum gmx_fft_direction dir
,
355 /*! \brief Perform a 3-dimensional complex-to-complex transform
357 * Performs an instance of a transform previously initiated.
359 * \param setup Setup returned from gmx_fft_init_1d()
360 * \param dir Forward or Backward
361 * \param in_data Input grid data. This should be allocated with gmx_new()
362 * to make it 16-byte aligned for better performance.
363 * \param out_data Output grid data. This should be allocated with gmx_new()
364 * to make it 16-byte aligned for better performance.
365 * You can provide the same pointer for in_data and out_data
366 * to perform an in-place transform.
368 * \return 0 on success, or an error code.
370 * \note Data pointers are declared as void, to avoid casting pointers
371 * depending on your grid type.
374 gmx_fft_3d (gmx_fft_t setup
,
375 enum gmx_fft_direction dir
,
380 /*! \brief Perform a 3-dimensional real-to-complex transform
382 * Performs an instance of a transform previously initiated.
384 * \param setup Setup returned from gmx_fft_init_1d_real()
385 * \param dir Real-to-complex or complex-to-real
386 * \param in_data Input grid data. This should be allocated with gmx_new()
387 * to make it 16-byte aligned for better performance.
388 * \param out_data Output grid data. This should be allocated with gmx_new()
389 * to make it 16-byte aligned for better performance.
390 * You can provide the same pointer for in_data and out_data
391 * to perform an in-place transform.
393 * \return 0 on success, or an error code.
395 * \note If you are doing an in-place transform, the last dimension of the
396 * array MUST be padded up to an even integer length so n/2 complex numbers can
397 * fit. Thus, if the real grid e.g. has dimension 7*5*3, you must allocate it as
398 * a 7*5*4 array, where the last element in the second dimension is padding.
399 * The complex data will be written to the same array, but since that dimension
400 * is 7*5*2 it will now fill the entire array. Reverse complex-to-real in-place
401 * transformation will produce the same sort of padded array.
403 * The padding does NOT apply to out-of-place transformation. In that case the
404 * input will simply be 7*5*3 of real, while the output is 7*5*2 of complex.
406 * \note Data pointers are declared as void, to avoid casting pointers
407 * depending on transform direction.
410 gmx_fft_3d_real (gmx_fft_t setup
,
411 enum gmx_fft_direction dir
,
416 /*! \brief Release an FFT setup structure
418 * Destroy setup and release all allocated memory.
420 * \param setup Setup returned from gmx_fft_init_1d(), or one
421 * of the other initializers.
425 gmx_fft_destroy (gmx_fft_t setup
);
428 /*! \brief Transpose 2d complex matrix, in-place or out-of-place.
430 * This routines works when the matrix is non-square, i.e. nx!=ny too,
431 * without allocating an entire matrix of work memory, which is important
432 * for huge FFT grids.
434 * \param in_data Input data, to be transposed
435 * \param out_data Output, transposed data. If this is identical to
436 * in_data, an in-place transpose is performed.
437 * \param nx Number of rows before transpose
438 * \param ny Number of columns before transpose
440 * \return GMX_SUCCESS, or an error code from gmx_errno.h
443 gmx_fft_transpose_2d (t_complex
* in_data
,
444 t_complex
* out_data
,
449 /*! \brief Transpose 2d multi-element matrix
451 * This routine is very similar to gmx_fft_transpose_2d(), but it
452 * supports matrices with more than one data value for each position.
453 * It is extremely useful when transposing the x/y dimensions of a 3d
454 * matrix - in that case you just set nelem to nz, and the routine will do
455 * and x/y transpose where it moves entire columns of z data
457 * This routines works when the matrix is non-square, i.e. nx!=ny too,
458 * without allocating an entire matrix of work memory, which is important
461 * For performance reasons you need to provide a \a small workarray
462 * with length at least 2*nelem (note that the type is char, not t_complex).
464 * \param in_data Input data, to be transposed
465 * \param out_data Output, transposed data. If this is identical to
466 * in_data, an in-place transpose is performed.
467 * \param nx Number of rows before transpose
468 * \param ny Number of columns before transpose
469 * \param nelem Number of t_complex values in each position. If this
470 * is 1 it is faster to use gmx_fft_transpose_2d() directly.
471 * \param work Work array of length 2*nelem, type t_complex.
473 * \return GMX_SUCCESS, or an error code from gmx_errno.h
476 gmx_fft_transpose_2d_nelem(t_complex
* in_data
,
477 t_complex
* out_data
,
490 #endif /* _GMX_FFT_H_ */