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