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