Fixing copyright issues and code contributors
[alexxy/gromacs.git] / src / mdlib / gmx_parallel_3dfft.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2005
5 * David van der Spoel, Erik Lindahl, University of Groningen.
6  * Copyright (c) 2012,2013, by the GROMACS development team, led by
7  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
8  * others, as listed in the AUTHORS file in the top-level source
9  * directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <stdlib.h>
42 #include <string.h>
43 #include <errno.h>
44
45 #ifdef GMX_LIB_MPI 
46 #include <mpi.h>
47 #endif
48 #ifdef GMX_THREAD_MPI
49 #include "tmpi.h"
50 #endif
51
52 #include "smalloc.h"
53 #include "gmx_parallel_3dfft.h"
54 #include "gmx_fft.h"
55 #include "gmxcomplex.h"
56 #include "gmx_fatal.h"
57
58 #include "fft5d.h"
59
60 struct gmx_parallel_3dfft  { 
61     fft5d_plan p1,p2;    
62 };
63
64 int
65 gmx_parallel_3dfft_init   (gmx_parallel_3dfft_t *    pfft_setup,
66                            ivec                      ndata,
67                                                    real **                   real_data,
68                                                    t_complex **              complex_data,
69                            MPI_Comm                  comm[2],
70                            int *                     slab2index_major,
71                            int *                     slab2index_minor,
72                            gmx_bool                      bReproducible,
73                            int nthreads)
74 {
75     int rN=ndata[2],M=ndata[1],K=ndata[0];
76     int flags = FFT5D_REALCOMPLEX | FFT5D_ORDER_YZ; /* FFT5D_DEBUG */
77     MPI_Comm rcomm[]={comm[1],comm[0]};
78     int Nb,Mb,Kb; /* dimension for backtransform (in starting order) */
79     t_complex *buf1, *buf2; /*intermediate buffers - used internally.*/
80     
81     snew(*pfft_setup,1);
82     if (bReproducible) flags |= FFT5D_NOMEASURE; 
83     
84     if (!(flags&FFT5D_ORDER_YZ)) { 
85         Nb=M;Mb=K;Kb=rN;                
86     } else {
87         Nb=K;Mb=rN;Kb=M;  /* currently always true because ORDER_YZ always set */
88     }
89     
90     (*pfft_setup)->p1 = fft5d_plan_3d(rN,M,K,rcomm, flags, (t_complex**)real_data, complex_data, &buf1, &buf2, nthreads);
91     
92     (*pfft_setup)->p2 = fft5d_plan_3d(Nb,Mb,Kb,rcomm,
93                                       (flags|FFT5D_BACKWARD|FFT5D_NOMALLOC)^FFT5D_ORDER_YZ, complex_data, (t_complex**)real_data, &buf1, &buf2, nthreads);
94     
95     return (*pfft_setup)->p1 != 0 && (*pfft_setup)->p2 !=0;
96 }
97
98
99 static int
100 fft5d_limits(fft5d_plan p, 
101              ivec                      local_ndata,
102              ivec                      local_offset,
103              ivec                      local_size) 
104 {
105     int N1,M0,K0,K1,*coor;
106     fft5d_local_size(p,&N1,&M0,&K0,&K1,&coor);  /* M0=MG/P[0], K1=KG/P[1], NG,MG,KG global sizes */
107     
108     local_offset[2]=0;
109     local_offset[1]=p->oM[0];  /*=p->coor[0]*p->MG/p->P[0]; */
110     local_offset[0]=p->oK[0];  /*=p->coor[1]*p->KG/p->P[1]; */
111     
112     local_ndata[2]=p->rC[0];
113     local_ndata[1]=p->pM[0]; 
114     local_ndata[0]=p->pK[0]; 
115     
116     if ((!(p->flags&FFT5D_BACKWARD)) && (p->flags&FFT5D_REALCOMPLEX)) {
117         local_size[2]=p->C[0]*2;
118     } else {
119         local_size[2]=p->C[0];
120     }
121     local_size[1]=p->pM[0]; 
122     local_size[0]=p->pK[0]; 
123     return 0;
124 }
125
126 int
127 gmx_parallel_3dfft_real_limits(gmx_parallel_3dfft_t      pfft_setup,
128                                ivec                      local_ndata,
129                                ivec                      local_offset,
130                                ivec                      local_size) {
131     return fft5d_limits(pfft_setup->p1,local_ndata,local_offset,local_size);
132 }
133
134 static void reorder_ivec_yzx(ivec v)
135 {
136     real tmp;
137
138     tmp   = v[0];
139     v[XX] = v[2];
140     v[ZZ] = v[1];
141     v[YY] = tmp;
142 }
143
144 int
145 gmx_parallel_3dfft_complex_limits(gmx_parallel_3dfft_t      pfft_setup,
146                                   ivec                      complex_order,
147                                   ivec                      local_ndata,
148                                   ivec                      local_offset,
149                                   ivec                      local_size) 
150 {
151     int ret;
152
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;
157
158     ret = fft5d_limits(pfft_setup->p2,local_ndata,local_offset,local_size);
159
160     reorder_ivec_yzx(local_ndata);
161     reorder_ivec_yzx(local_offset);
162     reorder_ivec_yzx(local_size);
163
164     return ret;
165 }
166
167
168 int
169 gmx_parallel_3dfft_execute(gmx_parallel_3dfft_t    pfft_setup,
170                            enum gmx_fft_direction  dir,
171                            void *                  in_data,
172                            void *                  out_data,
173                            int                     thread,
174                            gmx_wallcycle_t         wcycle) {
175     if ((!(pfft_setup->p1->flags&FFT5D_REALCOMPLEX)) ^ (dir==GMX_FFT_FORWARD ||dir==GMX_FFT_BACKWARD))
176     {
177         gmx_fatal(FARGS,"Invalid transform. Plan and execution don't match regarding reel/complex");
178     }
179     if (dir==GMX_FFT_FORWARD || dir==GMX_FFT_REAL_TO_COMPLEX) {
180         fft5d_execute(pfft_setup->p1,thread,wcycle);
181     } else {
182         fft5d_execute(pfft_setup->p2,thread,wcycle);
183     }
184     return 0;
185 }
186
187 int
188 gmx_parallel_3dfft_destroy(gmx_parallel_3dfft_t    pfft_setup) {
189     fft5d_destroy(pfft_setup->p2);
190     fft5d_destroy(pfft_setup->p1);
191     sfree(pfft_setup);
192     return 0;
193 }
194
195
196
197
198