Fix confusing gmx_header_config.h include
[gromacs.git] / include / gmx_fft.h
blob858597f0bb9b6ef8aaa94aa2dc0647654688002f
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Gromacs 4.0 Copyright (c) 1991-2003
5 * David van der Spoel, Erik Lindahl, University of Groningen.
6 * Copyright (c) 2012, by the GROMACS development team, led by
7 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
8 * others, as listed in the AUTHORS file in the top-level source
9 * directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifndef _GMX_FFT_H_
39 #define _GMX_FFT_H_
41 /*! \file gmx_fft.h
42 * \brief Fast Fourier Transforms.
44 * This file provides an abstract Gromacs interface to Fourier transforms,
45 * including multi-dimensional and real-to-complex transforms.
47 * Internally it is implemented as wrappers to external libraries such
48 * as FFTW or the Intel Math Kernel Library, but we also have a built-in
49 * version of FFTPACK in case the faster alternatives are unavailable.
51 * We also provide our own multi-dimensional transform setups even when
52 * the underlying library does not support it directly.
56 #include <stdio.h>
57 #include "visibility.h"
58 #include "types/simple.h"
59 #include "gmxcomplex.h"
62 #ifdef __cplusplus
63 extern "C" {
64 #endif
65 #if 0
66 } /* fixes auto-indentation problems */
67 #endif
71 /*! \brief Datatype for FFT setup
73 * The gmx_fft_t type contains all the setup information, e.g. twiddle
74 * factors, necessary to perform an FFT. Internally it is mapped to
75 * whatever FFT library we are using, or the built-in FFTPACK if no fast
76 * external library is available.
78 * Since some of the libraries (e.g. MKL) store work array data in their
79 * handles this datatype should only be used for one thread at a time, i.e.
80 * they should allocate one instance each when executing in parallel.
82 typedef struct gmx_fft *
83 gmx_fft_t;
88 /*! \brief Specifier for FFT direction.
90 * The definition of the 1D forward transform from input x[] to output y[] is
91 * \f[
92 * y_{k} = \sum_{j=0}^{N-1} x_{j} \exp{-i 2 \pi j k /N}
93 * \f]
95 * while the corresponding backward transform is
97 * \f[
98 * y_{k} = \sum_{j=0}^{N-1} x_{j} \exp{i 2 \pi j k /N}
99 * \f]
101 * A forward-backward transform pair will this result in data scaled by N.
103 * For complex-to-complex transforms you can only use one of
104 * GMX_FFT_FORWARD or GMX_FFT_BACKWARD, and for real-complex transforms you
105 * can only use GMX_FFT_REAL_TO_COMPLEX or GMX_FFT_COMPLEX_TO_REAL.
107 typedef enum gmx_fft_direction
109 GMX_FFT_FORWARD, /*!< Forward complex-to-complex transform */
110 GMX_FFT_BACKWARD, /*!< Backward complex-to-complex transform */
111 GMX_FFT_REAL_TO_COMPLEX, /*!< Real-to-complex valued fft */
112 GMX_FFT_COMPLEX_TO_REAL /*!< Complex-to-real valued fft */
113 } gmx_fft_direction;
115 /*! \brief Specifier for FFT flags.
117 * Some FFT libraries (FFTW, in particular) can do timings and other
118 * tricks to try and optimize the FFT for the current architecture. However,
119 * this can also lead to results that differ between consecutive runs with
120 * identical input.
121 * To avoid this, the conservative flag will attempt to disable such
122 * optimization, but there are no guarantees since we cannot control what
123 * the FFT libraries do internally.
126 typedef int gmx_fft_flag;
127 static const int GMX_FFT_FLAG_NONE = 0;
128 static const int GMX_FFT_FLAG_CONSERVATIVE = (1<<0);
130 /*! \brief Setup a 1-dimensional complex-to-complex transform
132 * \param fft Pointer to opaque Gromacs FFT datatype
133 * \param nx Length of transform
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 (gmx_fft_t * fft,
144 int nx,
145 gmx_fft_flag flags);
148 /*! \brief Setup multiple 1-dimensional complex-to-complex transform
150 * \param fft Pointer to opaque Gromacs FFT datatype
151 * \param nx Length of transform
152 * \param howmany Howmany 1D FFT
153 * \param flags FFT options
155 * \return status - 0 or a standard error message.
157 * \note Since some of the libraries (e.g. MKL) store work array data in their
158 * handles this datatype should only be used for one thread at a time,
159 * i.e. you should create one copy per thread when executing in parallel.
162 gmx_fft_init_many_1d (gmx_fft_t * fft,
163 int nx,
164 int howmany,
165 gmx_fft_flag flags);
168 /*! \brief Setup a 1-dimensional real-to-complex transform
170 * \param fft Pointer to opaque Gromacs FFT datatype
171 * \param nx Length of transform in real space
172 * \param flags FFT options
174 * \return status - 0 or a standard error message.
176 * \note Since some of the libraries (e.g. MKL) store work array data in their
177 * handles this datatype should only be used for one thread at a time,
178 * i.e. you should create one copy per thread when executing in parallel.
180 GMX_LIBMD_EXPORT
182 gmx_fft_init_1d_real (gmx_fft_t * fft,
183 int nx,
184 gmx_fft_flag flags);
187 /*! \brief Setup multiple 1-dimensional real-to-complex transform
189 * \param fft Pointer to opaque Gromacs FFT datatype
190 * \param nx Length of transform in real space
191 * \param howmany Homany 1D FFTs
192 * \param flags FFT options
194 * \return status - 0 or a standard error message.
196 * \note Since some of the libraries (e.g. MKL) store work array data in their
197 * handles this datatype should only be used for one thread at a time,
198 * i.e. you should create one copy per thread when executing in parallel.
201 gmx_fft_init_many_1d_real (gmx_fft_t * fft,
202 int nx,
203 int howmany,
204 gmx_fft_flag flags);
208 /*! \brief Setup a 2-dimensional complex-to-complex transform
210 * \param fft Pointer to opaque Gromacs FFT datatype
211 * \param nx Length of transform in first dimension
212 * \param ny Length of transform in second dimension
213 * \param flags FFT options
215 * \return status - 0 or a standard error message.
217 * \note Since some of the libraries (e.g. MKL) store work array data in their
218 * handles this datatype should only be used for one thread at a time,
219 * i.e. you should create one copy per thread when executing in parallel.
222 gmx_fft_init_2d (gmx_fft_t * fft,
223 int nx,
224 int ny,
225 gmx_fft_flag flags);
228 /*! \brief Setup a 2-dimensional real-to-complex transform
230 * \param fft Pointer to opaque Gromacs FFT datatype
231 * \param nx Length of transform in first dimension
232 * \param ny Length of transform in second dimension
233 * \param flags FFT options
235 * The normal space is assumed to be real, while the values in
236 * frequency space are complex.
238 * \return status - 0 or a standard error message.
240 * \note Since some of the libraries (e.g. MKL) store work array data in their
241 * handles this datatype should only be used for one thread at a time,
242 * i.e. you should create one copy per thread when executing in parallel.
244 GMX_LIBMD_EXPORT
246 gmx_fft_init_2d_real (gmx_fft_t * fft,
247 int nx,
248 int ny,
249 gmx_fft_flag flags);
252 /*! \brief Setup a 3-dimensional complex-to-complex transform
254 * \param fft Pointer to opaque Gromacs FFT datatype
255 * \param nx Length of transform in first dimension
256 * \param ny Length of transform in second dimension
257 * \param nz Length of transform in third dimension
258 * \param flags FFT options
260 * \return status - 0 or a standard error message.
262 * \note Since some of the libraries (e.g. MKL) store work array data in their
263 * handles this datatype should only be used for one thread at a time,
264 * i.e. you should create one copy per thread when executing in parallel.
267 gmx_fft_init_3d (gmx_fft_t * fft,
268 int nx,
269 int ny,
270 int nz,
271 gmx_fft_flag flags);
274 /*! \brief Setup a 3-dimensional real-to-complex transform
276 * \param fft Pointer to opaque Gromacs FFT datatype
277 * \param nx Length of transform in first dimension
278 * \param ny Length of transform in second dimension
279 * \param nz Length of transform in third dimension
280 * \param flags FFT options
282 * The normal space is assumed to be real, while the values in
283 * frequency space are complex.
285 * \return status - 0 or a standard error message.
287 * \note Since some of the libraries (e.g. MKL) store work array data in their
288 * handles this datatype should only be used for one thread at a time,
289 * i.e. you should create one copy per thread when executing in parallel.
292 gmx_fft_init_3d_real (gmx_fft_t * fft,
293 int nx,
294 int ny,
295 int nz,
296 gmx_fft_flag flags);
300 /*! \brief Perform a 1-dimensional complex-to-complex transform
302 * Performs an instance of a transform previously initiated.
304 * \param setup Setup returned from gmx_fft_init_1d()
305 * \param dir Forward or Backward
306 * \param in_data Input grid data. This should be allocated with gmx_new()
307 * to make it 16-byte aligned for better performance.
308 * \param out_data Output grid data. This should be allocated with gmx_new()
309 * to make it 16-byte aligned for better performance.
310 * You can provide the same pointer for in_data and out_data
311 * to perform an in-place transform.
313 * \return 0 on success, or an error code.
315 * \note Data pointers are declared as void, to avoid casting pointers
316 * depending on your grid type.
318 int
319 gmx_fft_1d (gmx_fft_t setup,
320 enum gmx_fft_direction dir,
321 void * in_data,
322 void * out_data);
325 /*! \brief Perform many 1-dimensional complex-to-complex transforms
327 * Performs many instances of a transform previously initiated.
329 * \param setup Setup returned from gmx_fft_init_1d()
330 * \param dir Forward or Backward
331 * \param in_data Input grid data. This should be allocated with gmx_new()
332 * to make it 16-byte aligned for better performance.
333 * \param out_data Output grid data. This should be allocated with gmx_new()
334 * to make it 16-byte aligned for better performance.
335 * You can provide the same pointer for in_data and out_data
336 * to perform an in-place transform.
338 * \return 0 on success, or an error code.
340 * \note Data pointers are declared as void, to avoid casting pointers
341 * depending on your grid type.
343 int
344 gmx_fft_many_1d (gmx_fft_t setup,
345 enum gmx_fft_direction dir,
346 void * in_data,
347 void * out_data);
350 /*! \brief Perform a 1-dimensional real-to-complex transform
352 * Performs an instance of a transform previously initiated.
354 * \param setup Setup returned from gmx_fft_init_1d_real()
355 * \param dir Real-to-complex or complex-to-real
356 * \param in_data Input grid data. This should be allocated with gmx_new()
357 * to make it 16-byte aligned for better performance.
358 * \param out_data Output grid data. This should be allocated with gmx_new()
359 * to make it 16-byte aligned for better performance.
360 * You can provide the same pointer for in_data and out_data
361 * to perform an in-place transform.
363 * If you are doing an in-place transform, the array must be padded up to
364 * an even integer length so n/2 complex numbers can fit. Out-of-place arrays
365 * should not be padded (although it doesn't matter in 1d).
367 * \return 0 on success, or an error code.
369 * \note Data pointers are declared as void, to avoid casting pointers
370 * depending on transform direction.
372 GMX_LIBMD_EXPORT
373 int
374 gmx_fft_1d_real (gmx_fft_t setup,
375 enum gmx_fft_direction dir,
376 void * in_data,
377 void * out_data);
379 /*! \brief Perform many 1-dimensional real-to-complex transforms
381 * Performs many instances of a transform previously initiated.
383 * \param setup Setup returned from gmx_fft_init_1d_real()
384 * \param dir Real-to-complex or complex-to-real
385 * \param in_data Input grid data. This should be allocated with gmx_new()
386 * to make it 16-byte aligned for better performance.
387 * \param out_data Output grid data. This should be allocated with gmx_new()
388 * to make it 16-byte aligned for better performance.
389 * You can provide the same pointer for in_data and out_data
390 * to perform an in-place transform.
392 * If you are doing an in-place transform, the array must be padded up to
393 * an even integer length so n/2 complex numbers can fit. Out-of-place arrays
394 * should not be padded (although it doesn't matter in 1d).
396 * \return 0 on success, or an error code.
398 * \note Data pointers are declared as void, to avoid casting pointers
399 * depending on transform direction.
401 int
402 gmx_fft_many_1d_real (gmx_fft_t setup,
403 enum gmx_fft_direction dir,
404 void * in_data,
405 void * out_data);
408 /*! \brief Perform a 2-dimensional complex-to-complex transform
410 * Performs an instance of a transform previously initiated.
412 * \param setup Setup returned from gmx_fft_init_1d()
413 * \param dir Forward or Backward
414 * \param in_data Input grid data. This should be allocated with gmx_new()
415 * to make it 16-byte aligned for better performance.
416 * \param out_data Output grid data. This should be allocated with gmx_new()
417 * to make it 16-byte aligned for better performance.
418 * You can provide the same pointer for in_data and out_data
419 * to perform an in-place transform.
421 * \return 0 on success, or an error code.
423 * \note Data pointers are declared as void, to avoid casting pointers
424 * depending on your grid type.
426 int
427 gmx_fft_2d (gmx_fft_t setup,
428 enum gmx_fft_direction dir,
429 void * in_data,
430 void * out_data);
433 /*! \brief Perform a 2-dimensional real-to-complex transform
435 * Performs an instance of a transform previously initiated.
437 * \param setup Setup returned from gmx_fft_init_1d_real()
438 * \param dir Real-to-complex or complex-to-real
439 * \param in_data Input grid data. This should be allocated with gmx_new()
440 * to make it 16-byte aligned for better performance.
441 * \param out_data Output grid data. This should be allocated with gmx_new()
442 * to make it 16-byte aligned for better performance.
443 * You can provide the same pointer for in_data and out_data
444 * to perform an in-place transform.
446 * \return 0 on success, or an error code.
448 * \note If you are doing an in-place transform, the last dimension of the
449 * array MUST be padded up to an even integer length so n/2 complex numbers can
450 * fit. Thus, if the real grid e.g. has dimension 5*3, you must allocate it as
451 * a 5*4 array, where the last element in the second dimension is padding.
452 * The complex data will be written to the same array, but since that dimension
453 * is 5*2 it will now fill the entire array. Reverse complex-to-real in-place
454 * transformation will produce the same sort of padded array.
456 * The padding does NOT apply to out-of-place transformation. In that case the
457 * input array will simply be 5*3 of real, while the output is 5*2 of complex.
459 * \note Data pointers are declared as void, to avoid casting pointers
460 * depending on transform direction.
462 GMX_LIBMD_EXPORT
464 gmx_fft_2d_real (gmx_fft_t setup,
465 enum gmx_fft_direction dir,
466 void * in_data,
467 void * out_data);
470 /*! \brief Perform a 3-dimensional complex-to-complex transform
472 * Performs an instance of a transform previously initiated.
474 * \param setup Setup returned from gmx_fft_init_1d()
475 * \param dir Forward or Backward
476 * \param in_data Input grid data. This should be allocated with gmx_new()
477 * to make it 16-byte aligned for better performance.
478 * \param out_data Output grid data. This should be allocated with gmx_new()
479 * to make it 16-byte aligned for better performance.
480 * You can provide the same pointer for in_data and out_data
481 * to perform an in-place transform.
483 * \return 0 on success, or an error code.
485 * \note Data pointers are declared as void, to avoid casting pointers
486 * depending on your grid type.
488 int
489 gmx_fft_3d (gmx_fft_t setup,
490 enum gmx_fft_direction dir,
491 void * in_data,
492 void * out_data);
495 /*! \brief Perform a 3-dimensional real-to-complex transform
497 * Performs an instance of a transform previously initiated.
499 * \param setup Setup returned from gmx_fft_init_1d_real()
500 * \param dir Real-to-complex or complex-to-real
501 * \param in_data Input grid data. This should be allocated with gmx_new()
502 * to make it 16-byte aligned for better performance.
503 * \param out_data Output grid data. This should be allocated with gmx_new()
504 * to make it 16-byte aligned for better performance.
505 * You can provide the same pointer for in_data and out_data
506 * to perform an in-place transform.
508 * \return 0 on success, or an error code.
510 * \note If you are doing an in-place transform, the last dimension of the
511 * array MUST be padded up to an even integer length so n/2 complex numbers can
512 * fit. Thus, if the real grid e.g. has dimension 7*5*3, you must allocate it as
513 * a 7*5*4 array, where the last element in the second dimension is padding.
514 * The complex data will be written to the same array, but since that dimension
515 * is 7*5*2 it will now fill the entire array. Reverse complex-to-real in-place
516 * transformation will produce the same sort of padded array.
518 * The padding does NOT apply to out-of-place transformation. In that case the
519 * input will simply be 7*5*3 of real, while the output is 7*5*2 of complex.
521 * \note Data pointers are declared as void, to avoid casting pointers
522 * depending on transform direction.
524 int
525 gmx_fft_3d_real (gmx_fft_t setup,
526 enum gmx_fft_direction dir,
527 void * in_data,
528 void * out_data);
531 /*! \brief Release an FFT setup structure
533 * Destroy setup and release all allocated memory.
535 * \param setup Setup returned from gmx_fft_init_1d(), or one
536 * of the other initializers.
539 GMX_LIBMD_EXPORT
540 void
541 gmx_fft_destroy (gmx_fft_t setup);
543 /*! \brief Release a many FFT setup structure
545 * Destroy setup and release all allocated memory.
547 * \param setup Setup returned from gmx_fft_init_1d(), or one
548 * of the other initializers.
551 void
552 gmx_many_fft_destroy (gmx_fft_t setup);
555 /*! \brief Transpose 2d complex matrix, in-place or out-of-place.
557 * This routines works when the matrix is non-square, i.e. nx!=ny too,
558 * without allocating an entire matrix of work memory, which is important
559 * for huge FFT grids.
561 * \param in_data Input data, to be transposed
562 * \param out_data Output, transposed data. If this is identical to
563 * in_data, an in-place transpose is performed.
564 * \param nx Number of rows before transpose
565 * \param ny Number of columns before transpose
567 * \return GMX_SUCCESS, or an error code from gmx_errno.h
570 gmx_fft_transpose_2d (t_complex * in_data,
571 t_complex * out_data,
572 int nx,
573 int ny);
576 /*! \brief Transpose 2d multi-element matrix
578 * This routine is very similar to gmx_fft_transpose_2d(), but it
579 * supports matrices with more than one data value for each position.
580 * It is extremely useful when transposing the x/y dimensions of a 3d
581 * matrix - in that case you just set nelem to nz, and the routine will do
582 * and x/y transpose where it moves entire columns of z data
584 * This routines works when the matrix is non-square, i.e. nx!=ny too,
585 * without allocating an entire matrix of work memory, which is important
586 * for huge FFT grid.
588 * For performance reasons you need to provide a \a small workarray
589 * with length at least 2*nelem (note that the type is char, not t_complex).
591 * \param in_data Input data, to be transposed
592 * \param out_data Output, transposed data. If this is identical to
593 * in_data, an in-place transpose is performed.
594 * \param nx Number of rows before transpose
595 * \param ny Number of columns before transpose
596 * \param nelem Number of t_complex values in each position. If this
597 * is 1 it is faster to use gmx_fft_transpose_2d() directly.
598 * \param work Work array of length 2*nelem, type t_complex.
600 * \return GMX_SUCCESS, or an error code from gmx_errno.h
603 gmx_fft_transpose_2d_nelem(t_complex * in_data,
604 t_complex * out_data,
605 int nx,
606 int ny,
607 int nelem,
608 t_complex * work);
613 #ifdef __cplusplus
615 #endif
617 #endif /* _GMX_FFT_H_ */