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                            gmx_bool                      bReproducible,
47                            int                           nthreads)
48 {
49     int        rN      = ndata[2], M = ndata[1], K = ndata[0];
50     int        flags   = FFT5D_REALCOMPLEX | FFT5D_ORDER_YZ; /* FFT5D_DEBUG */
51     MPI_Comm   rcomm[] = {comm[1], comm[0]};
52     int        Nb, Mb, Kb;                                   /* dimension for backtransform (in starting order) */
53     t_complex *buf1, *buf2;                                  /*intermediate buffers - used internally.*/
54
55     snew(*pfft_setup, 1);
56     if (bReproducible)
57     {
58         flags |= FFT5D_NOMEASURE;
59     }
60
61     if (!(flags&FFT5D_ORDER_YZ))
62     {
63         Nb = M; Mb = K; Kb = rN;
64     }
65     else
66     {
67         Nb = K; Mb = rN; Kb = M;  /* currently always true because ORDER_YZ always set */
68     }
69
70     (*pfft_setup)->p1 = fft5d_plan_3d(rN, M, K, rcomm, flags, (t_complex**)real_data, complex_data, &buf1, &buf2, nthreads);
71
72     (*pfft_setup)->p2 = fft5d_plan_3d(Nb, Mb, Kb, rcomm,
73                                       (flags|FFT5D_BACKWARD|FFT5D_NOMALLOC)^FFT5D_ORDER_YZ, complex_data, (t_complex**)real_data, &buf1, &buf2, nthreads);
74
75     return (*pfft_setup)->p1 != 0 && (*pfft_setup)->p2 != 0;
76 }
77
78
79 static int
80 fft5d_limits(fft5d_plan                p,
81              ivec                      local_ndata,
82              ivec                      local_offset,
83              ivec                      local_size)
84 {
85     local_offset[2] = 0;
86     local_offset[1] = p->oM[0];  /*=p->coor[0]*p->MG/p->P[0]; */
87     local_offset[0] = p->oK[0];  /*=p->coor[1]*p->KG/p->P[1]; */
88
89     local_ndata[2] = p->rC[0];
90     local_ndata[1] = p->pM[0];
91     local_ndata[0] = p->pK[0];
92
93     if ((!(p->flags&FFT5D_BACKWARD)) && (p->flags&FFT5D_REALCOMPLEX))
94     {
95         //C is length in multiples of complex local_size in multiples of real
96         local_size[2] = p->C[0]*2;
97     }
98     else
99     {
100         local_size[2] = p->C[0];
101     }
102     local_size[1] = p->pM[0];
103     local_size[0] = p->pK[0];
104     return 0;
105 }
106
107 int
108 gmx_parallel_3dfft_real_limits(gmx_parallel_3dfft_t      pfft_setup,
109                                ivec                      local_ndata,
110                                ivec                      local_offset,
111                                ivec                      local_size)
112 {
113     return fft5d_limits(pfft_setup->p1, local_ndata, local_offset, local_size);
114 }
115
116 static void reorder_ivec_yzx(ivec v)
117 {
118     real tmp;
119
120     tmp   = v[0];
121     v[XX] = v[2];
122     v[ZZ] = v[1];
123     v[YY] = tmp;
124 }
125
126 int
127 gmx_parallel_3dfft_complex_limits(gmx_parallel_3dfft_t      pfft_setup,
128                                   ivec                      complex_order,
129                                   ivec                      local_ndata,
130                                   ivec                      local_offset,
131                                   ivec                      local_size)
132 {
133     int ret;
134
135     /* For now everything is in-order, but prepare to save communication by avoiding transposes */
136     complex_order[0] = 0;
137     complex_order[1] = 1;
138     complex_order[2] = 2;
139
140     ret = fft5d_limits(pfft_setup->p2, local_ndata, local_offset, local_size);
141
142     reorder_ivec_yzx(local_ndata);
143     reorder_ivec_yzx(local_offset);
144     reorder_ivec_yzx(local_size);
145
146     return ret;
147 }
148
149
150 int
151 gmx_parallel_3dfft_execute(gmx_parallel_3dfft_t    pfft_setup,
152                            enum gmx_fft_direction  dir,
153                            int                     thread,
154                            gmx_wallcycle_t         wcycle)
155 {
156     if ((!(pfft_setup->p1->flags&FFT5D_REALCOMPLEX)) ^ (dir == GMX_FFT_FORWARD || dir == GMX_FFT_BACKWARD))
157     {
158         gmx_fatal(FARGS, "Invalid transform. Plan and execution don't match regarding reel/complex");
159     }
160     if (dir == GMX_FFT_FORWARD || dir == GMX_FFT_REAL_TO_COMPLEX)
161     {
162         fft5d_execute(pfft_setup->p1, thread, wcycle);
163     }
164     else
165     {
166         fft5d_execute(pfft_setup->p2, thread, wcycle);
167     }
168     return 0;
169 }
170
171 int
172 gmx_parallel_3dfft_destroy(gmx_parallel_3dfft_t    pfft_setup)
173 {
174     if (pfft_setup)
175     {
176         fft5d_destroy(pfft_setup->p2);
177         fft5d_destroy(pfft_setup->p1);
178         sfree(pfft_setup);
179     }
180     return 0;
181 }