fef1ac6083227088a1adfe470d5f8310d3f930f6
[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_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) flags |= FFT5D_NOMEASURE; 
64     
65     if (!(flags&FFT5D_ORDER_YZ)) { 
66         Nb=M;Mb=K;Kb=rN;                
67     } else {
68         Nb=K;Mb=rN;Kb=M;  /* currently always true because ORDER_YZ always set */
69     }
70     
71     (*pfft_setup)->p1 = fft5d_plan_3d(rN,M,K,rcomm, flags, (t_complex**)real_data, complex_data, &buf1, &buf2, nthreads);
72     
73     (*pfft_setup)->p2 = fft5d_plan_3d(Nb,Mb,Kb,rcomm,
74                                       (flags|FFT5D_BACKWARD|FFT5D_NOMALLOC)^FFT5D_ORDER_YZ, complex_data, (t_complex**)real_data, &buf1, &buf2, nthreads);
75     
76     return (*pfft_setup)->p1 != 0 && (*pfft_setup)->p2 !=0;
77 }
78
79
80 static int
81 fft5d_limits(fft5d_plan p, 
82              ivec                      local_ndata,
83              ivec                      local_offset,
84              ivec                      local_size) 
85 {
86     int N1,M0,K0,K1,*coor;
87     fft5d_local_size(p,&N1,&M0,&K0,&K1,&coor);  /* M0=MG/P[0], K1=KG/P[1], NG,MG,KG global sizes */
88     
89     local_offset[2]=0;
90     local_offset[1]=p->oM[0];  /*=p->coor[0]*p->MG/p->P[0]; */
91     local_offset[0]=p->oK[0];  /*=p->coor[1]*p->KG/p->P[1]; */
92     
93     local_ndata[2]=p->rC[0];
94     local_ndata[1]=p->pM[0]; 
95     local_ndata[0]=p->pK[0]; 
96     
97     if ((!(p->flags&FFT5D_BACKWARD)) && (p->flags&FFT5D_REALCOMPLEX)) {
98         local_size[2]=p->C[0]*2;
99     } else {
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     return fft5d_limits(pfft_setup->p1,local_ndata,local_offset,local_size);
113 }
114
115 static void reorder_ivec_yzx(ivec v)
116 {
117     real tmp;
118
119     tmp   = v[0];
120     v[XX] = v[2];
121     v[ZZ] = v[1];
122     v[YY] = tmp;
123 }
124
125 int
126 gmx_parallel_3dfft_complex_limits(gmx_parallel_3dfft_t      pfft_setup,
127                                   ivec                      complex_order,
128                                   ivec                      local_ndata,
129                                   ivec                      local_offset,
130                                   ivec                      local_size) 
131 {
132     int ret;
133
134     /* For now everything is in-order, but prepare to save communication by avoiding transposes */
135     complex_order[0] = 0;
136     complex_order[1] = 1;
137     complex_order[2] = 2;
138
139     ret = fft5d_limits(pfft_setup->p2,local_ndata,local_offset,local_size);
140
141     reorder_ivec_yzx(local_ndata);
142     reorder_ivec_yzx(local_offset);
143     reorder_ivec_yzx(local_size);
144
145     return ret;
146 }
147
148
149 int
150 gmx_parallel_3dfft_execute(gmx_parallel_3dfft_t    pfft_setup,
151                            enum gmx_fft_direction  dir,
152                            void *                  in_data,
153                            void *                  out_data,
154                            int                     thread,
155                            gmx_wallcycle_t         wcycle) {
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         fft5d_execute(pfft_setup->p1,thread,wcycle);
162     } else {
163         fft5d_execute(pfft_setup->p2,thread,wcycle);
164     }
165     return 0;
166 }
167
168 int
169 gmx_parallel_3dfft_destroy(gmx_parallel_3dfft_t    pfft_setup) {
170     fft5d_destroy(pfft_setup->p2);
171     fft5d_destroy(pfft_setup->p1);
172     sfree(pfft_setup);
173     return 0;
174 }
175
176
177
178
179