More modular position handling.
[gromacs/qmmm-gamess-us.git] / include / gmx_parallel_3dfft.h
blob47b6e9afeb89f7de30acfb664e6c1bba948e87d0
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * Gromacs Copyright (c) 1991-2005
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
15 * And Hey:
16 * Gnomes, ROck Monsters And Chili Sauce
19 #ifndef _gmx_parallel_3dfft_h_
20 #define _gmx_parallel_3dfft_h_
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
26 #ifdef GMX_MPI
28 #include "types/simple.h"
29 #include "gmxcomplex.h"
30 #include "gmx_fft.h"
32 /* We NEED MPI here. */
33 #ifdef GMX_LIB_MPI
34 #include <mpi.h>
35 #endif
36 #ifdef GMX_THREADS
37 #include "tmpi.h"
38 #endif
40 typedef struct gmx_parallel_3dfft *
41 gmx_parallel_3dfft_t;
45 /*! \brief Initialize parallel MPI-based 3D-FFT.
47 * This routine performs real-to-complex and complex-to-real parallel 3D FFTs,
48 * but not complex-to-complex.
50 * The routine is optimized for small-to-medium size FFTs used for PME and
51 * PPPM algorithms, and do allocate extra workspace whenever it might improve
52 * performance.
54 * \param pfft_setup Pointer to parallel 3dfft setup structure, previously
55 * allocated or with automatic storage.
56 * \param ngridx Global number of grid cells in the x direction. Must be
57 * divisible by the number of nodes.
58 * \param ngridy Global number of grid cells in the y direction. Must be
59 * divisible by the number of nodes.
60 * \param ngridz Global number of grid cells in the z direction.
61 * \param node2slab Node id to slab index array, can be NULL.
62 * \param slab2grid_x Slab index to grid_x array (nnodes+1), can be NULL.
63 * \param comm MPI communicator, must have been initialized.
64 * \param bReproducible Try to avoid FFT timing optimizations and other stuff
65 * that could make results differ for two runs with
66 * identical input (reproducibility for debugging).
68 * \return 0 or a standard error code.
70 int
71 gmx_parallel_3dfft_init (gmx_parallel_3dfft_t * pfft_setup,
72 int ngridx,
73 int ngridy,
74 int ngridz,
75 int *node2slab,
76 int *slab2grid_x,
77 MPI_Comm comm,
78 bool bReproducible);
84 /*! \brief Get direct space grid index limits
86 * The z dimension is never distributed. In the direct space, the x dimension
87 * is distributed over nodes, and after the real-to-complex FFT we work with
88 * a transposed grid where the y dimension is partitioned over nodes.
90 * The node2slab array translates to node ids to slab indices,
91 * when NULL the slab ids are assumed to be identical to the node ids
92 * in the communicator comm.
94 int
95 gmx_parallel_3dfft_limits(gmx_parallel_3dfft_t pfft_setup,
96 int * local_x_start,
97 int * local_nx,
98 int * local_y_start,
99 int * local_ny);
103 gmx_parallel_transpose(t_complex * data,
104 t_complex * work,
105 int nx,
106 int ny,
107 int local_x_start,
108 int local_nx,
109 int local_y_start,
110 int local_ny,
111 int nelem,
112 int nnodes,
113 int *node2slab,
114 MPI_Comm comm);
117 /*! \brief Perform forward parallel MPI FFT.
119 * Direction is either GMX_FFT_REAL_TO_COMPLEX or GMX_FFT_COMPLEX_TO_REAL.
121 * If input and output arrays are separate there is no packing to consider.
122 * Input is simply nx*ny*nz in real, and output ny*nx*nzc in complex.
124 * In they are identical we need to make sure there is room for the complex
125 * (length nzc=nz/2+1) in the array, so the _real_ space dimensions is
126 * always padded to nzc*2.
127 * In this case, the real dimensions are nx*ny*(nzc*2) while the complex
128 * dimensions is ny*nx*nzc (of type complex).
130 * Note that the X and Y dimensions are transposed in the reciprocal space
131 * to avoid extra communication!
133 * The node2slab array translates to node ids to slab indices,
134 * when NULL the slab ids are assumed to be identical to the node ids
135 * in the communicator comm.
138 gmx_parallel_3dfft(gmx_parallel_3dfft_t pfft_setup,
139 enum gmx_fft_direction dir,
140 void * in_data,
141 void * out_data);
145 /*! \brief Release all data in parallel fft setup
147 * All temporary storage and FFT plans are released. The structure itself
148 * is not released, but the contents is invalid after this call.
150 * \param pfft_setup Parallel 3dfft setup.
152 * \return 0 or a standard error code.
155 gmx_parallel_3dfft_destroy(gmx_parallel_3dfft_t pfft_setup);
157 #endif /* GMX_MPI */
159 #endif /* _gmx_parallel_3dfft_h_ */