Merge branch 'release-4-6' into master
[alexxy/gromacs.git] / src / gromacs / fft / parallel_3dfft.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  * Gromacs                               Copyright (c) 1991-2005
5  * David van der Spoel, Erik Lindahl, University of Groningen.
6  *
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.
11  *
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
14  *
15  * And Hey:
16  * Gnomes, ROck Monsters And Chili Sauce
17  */
18 #ifdef HAVE_CONFIG_H
19 #include <config.h>
20 #endif
21
22 #include <stdlib.h>
23 #include <string.h>
24 #include <errno.h>
25
26 #include "gromacs/fft/parallel_3dfft.h"
27 #include "gromacs/fft/fft.h"
28 #include "gromacs/utility/gmxmpi.h"
29
30 #include "smalloc.h"
31 #include "gmxcomplex.h"
32 #include "gmx_fatal.h"
33
34 #include "fft5d.h"
35
36 struct gmx_parallel_3dfft  {
37     fft5d_plan p1, p2;
38 };
39
40 int
41 gmx_parallel_3dfft_init   (gmx_parallel_3dfft_t     *    pfft_setup,
42                            ivec                          ndata,
43                            real     **                   real_data,
44                            t_complex     **              complex_data,
45                            MPI_Comm                      comm[2],
46                            int     *                     slab2index_major,
47                            int     *                     slab2index_minor,
48                            gmx_bool                      bReproducible,
49                            int                           nthreads)
50 {
51     int        rN      = ndata[2], M = ndata[1], K = ndata[0];
52     int        flags   = FFT5D_REALCOMPLEX | FFT5D_ORDER_YZ; /* FFT5D_DEBUG */
53     MPI_Comm   rcomm[] = {comm[1], comm[0]};
54     int        Nb, Mb, Kb;                                   /* dimension for backtransform (in starting order) */
55     t_complex *buf1, *buf2;                                  /*intermediate buffers - used internally.*/
56
57     snew(*pfft_setup, 1);
58     if (bReproducible)
59     {
60         flags |= FFT5D_NOMEASURE;
61     }
62
63     if (!(flags&FFT5D_ORDER_YZ))
64     {
65         Nb = M; Mb = K; Kb = rN;
66     }
67     else
68     {
69         Nb = K; Mb = rN; Kb = M;  /* currently always true because ORDER_YZ always set */
70     }
71
72     (*pfft_setup)->p1 = fft5d_plan_3d(rN, M, K, rcomm, flags, (t_complex**)real_data, complex_data, &buf1, &buf2, nthreads);
73
74     (*pfft_setup)->p2 = fft5d_plan_3d(Nb, Mb, Kb, rcomm,
75                                       (flags|FFT5D_BACKWARD|FFT5D_NOMALLOC)^FFT5D_ORDER_YZ, complex_data, (t_complex**)real_data, &buf1, &buf2, nthreads);
76
77     return (*pfft_setup)->p1 != 0 && (*pfft_setup)->p2 != 0;
78 }
79
80
81 static int
82 fft5d_limits(fft5d_plan                p,
83              ivec                      local_ndata,
84              ivec                      local_offset,
85              ivec                      local_size)
86 {
87     local_offset[2] = 0;
88     local_offset[1] = p->oM[0];  /*=p->coor[0]*p->MG/p->P[0]; */
89     local_offset[0] = p->oK[0];  /*=p->coor[1]*p->KG/p->P[1]; */
90
91     local_ndata[2] = p->rC[0];
92     local_ndata[1] = p->pM[0];
93     local_ndata[0] = p->pK[0];
94
95     if ((!(p->flags&FFT5D_BACKWARD)) && (p->flags&FFT5D_REALCOMPLEX))
96     {
97         //C is length in multiples of complex local_size in multiples of real
98         local_size[2] = p->C[0]*2;
99     }
100     else
101     {
102         local_size[2] = p->C[0];
103     }
104     local_size[1] = p->pM[0];
105     local_size[0] = p->pK[0];
106     return 0;
107 }
108
109 int
110 gmx_parallel_3dfft_real_limits(gmx_parallel_3dfft_t      pfft_setup,
111                                ivec                      local_ndata,
112                                ivec                      local_offset,
113                                ivec                      local_size)
114 {
115     return fft5d_limits(pfft_setup->p1, local_ndata, local_offset, local_size);
116 }
117
118 static void reorder_ivec_yzx(ivec v)
119 {
120     real tmp;
121
122     tmp   = v[0];
123     v[XX] = v[2];
124     v[ZZ] = v[1];
125     v[YY] = tmp;
126 }
127
128 int
129 gmx_parallel_3dfft_complex_limits(gmx_parallel_3dfft_t      pfft_setup,
130                                   ivec                      complex_order,
131                                   ivec                      local_ndata,
132                                   ivec                      local_offset,
133                                   ivec                      local_size)
134 {
135     int ret;
136
137     /* For now everything is in-order, but prepare to save communication by avoiding transposes */
138     complex_order[0] = 0;
139     complex_order[1] = 1;
140     complex_order[2] = 2;
141
142     ret = fft5d_limits(pfft_setup->p2, local_ndata, local_offset, local_size);
143
144     reorder_ivec_yzx(local_ndata);
145     reorder_ivec_yzx(local_offset);
146     reorder_ivec_yzx(local_size);
147
148     return ret;
149 }
150
151
152 int
153 gmx_parallel_3dfft_execute(gmx_parallel_3dfft_t    pfft_setup,
154                            enum gmx_fft_direction  dir,
155                            void *                  in_data,
156                            void *                  out_data,
157                            int                     thread,
158                            gmx_wallcycle_t         wcycle)
159 {
160     if ((!(pfft_setup->p1->flags&FFT5D_REALCOMPLEX)) ^ (dir == GMX_FFT_FORWARD || dir == GMX_FFT_BACKWARD))
161     {
162         gmx_fatal(FARGS, "Invalid transform. Plan and execution don't match regarding reel/complex");
163     }
164     if (dir == GMX_FFT_FORWARD || dir == GMX_FFT_REAL_TO_COMPLEX)
165     {
166         fft5d_execute(pfft_setup->p1, thread, wcycle);
167     }
168     else
169     {
170         fft5d_execute(pfft_setup->p2, thread, wcycle);
171     }
172     return 0;
173 }
174
175 int
176 gmx_parallel_3dfft_destroy(gmx_parallel_3dfft_t    pfft_setup)
177 {
178     if (pfft_setup)
179     {
180         fft5d_destroy(pfft_setup->p2);
181         fft5d_destroy(pfft_setup->p1);
182         sfree(pfft_setup);
183     }
184     return 0;
185 }