3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gromacs Runs On Most of All Computer Systems
45 #include "gmxcomplex.h"
50 #include "gmx_parallel_3dfft.h"
55 #define INDEX(i,j,k) ((i)*la12+(j)*la2+(k))
58 int local_nx
,local_x_start
,local_ny_after_transpose
;
59 int local_y_start_after_transpose
;
66 int nx
,ny
,nz
,la2r
,la2c
,la12r
,la12c
;
70 gmx_parallel_3dfft_t mpi_fft_setup
;
75 extern t_fftgrid
*mk_fftgrid(int nx
,
83 /* Create an FFT grid (1 Dimensional), to be indexed by the INDEX macro
84 * Setup FFT plans and extract local sizes for the grid.
85 * If the file pointer is given, information is printed to it.
86 * If cr is non-NULL and cr->nnodes>1, a parallel grid and FFT will be created.
87 * The node2slab array translates to node ids to slab indices,
88 * when NULL the slab ids are assumed to be identical to the node ids.
89 * The slab2grid_x array determines which grid x-indices beling to which slab
90 * (array size nnodes+1), when NULL this is determined automatically.
91 * Set bReproducible to avoid FFTW timing and other optimizations that
92 * could affect reproducibility of simulations.
95 extern void pr_fftgrid(FILE *fp
,char *title
,t_fftgrid
*grid
);
96 /* Dump a grid to a file */
98 extern void done_fftgrid(t_fftgrid
*grid
);
99 /* And throw it away again */
101 extern void gmxfft3D(t_fftgrid
*grid
,enum gmx_fft_direction dir
,t_commrec
*cr
);
102 /* Do the FFT, direction may be either
103 * FFTW_FORWARD (sign -1) for real -> complex transform
104 * FFTW_BACKWARD (sign 1) for complex -> real transform
107 extern void clear_fftgrid(t_fftgrid
*grid
);
110 extern void unpack_fftgrid(t_fftgrid
*grid
,int *nx
,int *ny
,int *nz
,
111 int *nx2
,int *ny2
,int *nz2
,
112 int *la2
, int *la12
,bool bReal
, real
**ptr
);
114 /* Get the values for the constants into local copies */
119 /************************************************************************
121 * For backward compatibility (for testing the ewald code vs. PPPM etc)
122 * some old grid routines are retained here.
124 ************************************************************************/
126 extern real
***mk_rgrid(int nx
,int ny
,int nz
);
128 extern void free_rgrid(real
***grid
,int nx
,int ny
);
130 extern real
print_rgrid(FILE *fp
,char *title
,int nx
,int ny
,int nz
,
133 extern void print_rgrid_pdb(char *fn
,int nx
,int ny
,int nz
,real
***grid
);
135 extern t_complex
***mk_cgrid(int nx
,int ny
,int nz
);
137 extern void free_cgrid(t_complex
***grid
,int nx
,int ny
);
139 extern t_complex
print_cgrid(FILE *fp
,char *title
,int nx
,int ny
,int nz
,
142 extern void clear_cgrid(int nx
,int ny
,int nz
,t_complex
***grid
);
144 extern void clear_rgrid(int nx
,int ny
,int nz
,real
***grid
);