2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2005 David van der Spoel, Erik Lindahl, University of Groningen.
5 * Copyright (c) 2013, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
44 #include "gromacs/fft/parallel_3dfft.h"
45 #include "gromacs/fft/fft.h"
46 #include "gromacs/utility/gmxmpi.h"
49 #include "gmxcomplex.h"
50 #include "gmx_fatal.h"
54 struct gmx_parallel_3dfft {
59 gmx_parallel_3dfft_init (gmx_parallel_3dfft_t * pfft_setup,
62 t_complex ** complex_data,
64 gmx_bool bReproducible,
67 int rN = ndata[2], M = ndata[1], K = ndata[0];
68 int flags = FFT5D_REALCOMPLEX | FFT5D_ORDER_YZ; /* FFT5D_DEBUG */
69 MPI_Comm rcomm[] = {comm[1], comm[0]};
70 int Nb, Mb, Kb; /* dimension for backtransform (in starting order) */
71 t_complex *buf1, *buf2; /*intermediate buffers - used internally.*/
76 flags |= FFT5D_NOMEASURE;
79 if (!(flags&FFT5D_ORDER_YZ))
81 Nb = M; Mb = K; Kb = rN;
85 Nb = K; Mb = rN; Kb = M; /* currently always true because ORDER_YZ always set */
88 (*pfft_setup)->p1 = fft5d_plan_3d(rN, M, K, rcomm, flags, (t_complex**)real_data, complex_data, &buf1, &buf2, nthreads);
90 (*pfft_setup)->p2 = fft5d_plan_3d(Nb, Mb, Kb, rcomm,
91 (flags|FFT5D_BACKWARD|FFT5D_NOMALLOC)^FFT5D_ORDER_YZ, complex_data, (t_complex**)real_data, &buf1, &buf2, nthreads);
93 return (*pfft_setup)->p1 != 0 && (*pfft_setup)->p2 != 0;
98 fft5d_limits(fft5d_plan p,
104 local_offset[1] = p->oM[0]; /*=p->coor[0]*p->MG/p->P[0]; */
105 local_offset[0] = p->oK[0]; /*=p->coor[1]*p->KG/p->P[1]; */
107 local_ndata[2] = p->rC[0];
108 local_ndata[1] = p->pM[0];
109 local_ndata[0] = p->pK[0];
111 if ((!(p->flags&FFT5D_BACKWARD)) && (p->flags&FFT5D_REALCOMPLEX))
113 //C is length in multiples of complex local_size in multiples of real
114 local_size[2] = p->C[0]*2;
118 local_size[2] = p->C[0];
120 local_size[1] = p->pM[0];
121 local_size[0] = p->pK[0];
126 gmx_parallel_3dfft_real_limits(gmx_parallel_3dfft_t pfft_setup,
131 return fft5d_limits(pfft_setup->p1, local_ndata, local_offset, local_size);
134 static void reorder_ivec_yzx(ivec v)
145 gmx_parallel_3dfft_complex_limits(gmx_parallel_3dfft_t pfft_setup,
153 /* For now everything is in-order, but prepare to save communication by avoiding transposes */
154 complex_order[0] = 0;
155 complex_order[1] = 1;
156 complex_order[2] = 2;
158 ret = fft5d_limits(pfft_setup->p2, local_ndata, local_offset, local_size);
160 reorder_ivec_yzx(local_ndata);
161 reorder_ivec_yzx(local_offset);
162 reorder_ivec_yzx(local_size);
169 gmx_parallel_3dfft_execute(gmx_parallel_3dfft_t pfft_setup,
170 enum gmx_fft_direction dir,
172 gmx_wallcycle_t wcycle)
174 if ((!(pfft_setup->p1->flags&FFT5D_REALCOMPLEX)) ^ (dir == GMX_FFT_FORWARD || dir == GMX_FFT_BACKWARD))
176 gmx_fatal(FARGS, "Invalid transform. Plan and execution don't match regarding reel/complex");
178 if (dir == GMX_FFT_FORWARD || dir == GMX_FFT_REAL_TO_COMPLEX)
180 fft5d_execute(pfft_setup->p1, thread, wcycle);
184 fft5d_execute(pfft_setup->p2, thread, wcycle);
190 gmx_parallel_3dfft_destroy(gmx_parallel_3dfft_t pfft_setup)
194 fft5d_destroy(pfft_setup->p2);
195 fft5d_destroy(pfft_setup->p1);