Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_internal.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 #ifndef _nbnxn_internal_h
37 #define _nbnxn_internal_h
38
39 #include "gromacs/domdec/domdec.h"
40 #include "gromacs/legacyheaders/typedefs.h"
41 #include "gromacs/mdlib/nbnxn_pairlist.h"
42 #include "gromacs/mdlib/nbnxn_simd.h"
43 #include "gromacs/timing/cyclecounter.h"
44
45
46 /* Bounding box calculations are (currently) always in single precision, so
47  * we only need to check for single precision support here.
48  * This uses less (cache-)memory and SIMD is faster, at least on x86.
49  */
50 #ifdef GMX_SIMD4_HAVE_FLOAT
51 #define NBNXN_SEARCH_BB_SIMD4
52 #endif
53
54
55 #ifdef __cplusplus
56 extern "C" {
57 #endif
58
59
60 #ifdef GMX_NBNXN_SIMD
61 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
62 #define NBNXN_MEM_ALIGN  (GMX_SIMD_REAL_WIDTH*sizeof(real))
63 #else
64 /* No alignment required, but set it so we can call the same routines */
65 #define NBNXN_MEM_ALIGN  32
66 #endif
67
68
69 #ifdef NBNXN_SEARCH_BB_SIMD4
70 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
71 #define NBNXN_SEARCH_BB_MEM_ALIGN  (GMX_SIMD4_WIDTH*sizeof(float))
72 #else
73 /* No alignment required, but set it so we can call the same routines */
74 #define NBNXN_SEARCH_BB_MEM_ALIGN  32
75 #endif
76
77
78 /* Pair search box lower and upper corner in x,y,z.
79  * Store this in 4 iso 3 reals, which is useful with 4-wide SIMD.
80  * To avoid complicating the code we also use 4 without 4-wide SIMD.
81  */
82 #define NNBSBB_C         4
83 /* Pair search box lower and upper bound in z only. */
84 #define NNBSBB_D         2
85 /* Pair search box lower and upper corner x,y,z indices, entry 3 is unused */
86 #define BB_X  0
87 #define BB_Y  1
88 #define BB_Z  2
89
90 /* Bounding box for a nbnxn atom cluster */
91 typedef struct {
92     float lower[NNBSBB_C];
93     float upper[NNBSBB_C];
94 } nbnxn_bb_t;
95
96
97 /* A pair-search grid struct for one domain decomposition zone */
98 typedef struct {
99     rvec          c0;               /* The lower corner of the (local) grid        */
100     rvec          c1;               /* The upper corner of the (local) grid        */
101     rvec          size;             /* c1 - c0                                     */
102     real          atom_density;     /* The atom number density for the local grid  */
103
104     gmx_bool      bSimple;          /* Is this grid simple or super/sub            */
105     int           na_c;             /* Number of atoms per cluster                 */
106     int           na_cj;            /* Number of atoms for list j-clusters         */
107     int           na_sc;            /* Number of atoms per super-cluster           */
108     int           na_c_2log;        /* 2log of na_c                                */
109
110     int           ncx;              /* Number of (super-)cells along x             */
111     int           ncy;              /* Number of (super-)cells along y             */
112     int           nc;               /* Total number of (super-)cells               */
113
114     real          sx;               /* x-size of a (super-)cell                    */
115     real          sy;               /* y-size of a (super-)cell                    */
116     real          inv_sx;           /* 1/sx                                        */
117     real          inv_sy;           /* 1/sy                                        */
118
119     int           cell0;            /* Index in nbs->cell corresponding to cell 0  */
120
121     int          *cxy_na;           /* The number of atoms for each column in x,y  */
122     int          *cxy_ind;          /* Grid (super)cell index, offset from cell0   */
123     int           cxy_nalloc;       /* Allocation size for cxy_na and cxy_ind      */
124
125     int          *nsubc;            /* The number of sub cells for each super cell */
126     float        *bbcz;             /* Bounding boxes in z for the super cells     */
127     nbnxn_bb_t   *bb;               /* 3D bounding boxes for the sub cells         */
128     nbnxn_bb_t   *bbj;              /* 3D j-bounding boxes for the case where      *
129                                      * the i- and j-cluster sizes are different    */
130     float        *pbb;              /* 3D b. boxes in xxxx format per super cell   */
131     int          *flags;            /* Flag for the super cells                    */
132     unsigned int *fep;              /* FEP signal bits for sub cells               */
133     int           nc_nalloc;        /* Allocation size for the pointers above      */
134
135     float        *bbcz_simple;      /* bbcz for simple grid converted from super   */
136     nbnxn_bb_t   *bb_simple;        /* bb for simple grid converted from super     */
137     int          *flags_simple;     /* flags for simple grid converted from super  */
138     int           nc_nalloc_simple; /* Allocation size for the pointers above   */
139
140     int           nsubc_tot;        /* Total number of subcell, used for printing  */
141 } nbnxn_grid_t;
142
143 #ifdef GMX_NBNXN_SIMD
144
145 typedef struct nbnxn_x_ci_simd_4xn {
146     /* The i-cluster coordinates for simple search */
147     gmx_simd_real_t ix_S0, iy_S0, iz_S0;
148     gmx_simd_real_t ix_S1, iy_S1, iz_S1;
149     gmx_simd_real_t ix_S2, iy_S2, iz_S2;
150     gmx_simd_real_t ix_S3, iy_S3, iz_S3;
151 } nbnxn_x_ci_simd_4xn_t;
152
153 typedef struct nbnxn_x_ci_simd_2xnn {
154     /* The i-cluster coordinates for simple search */
155     gmx_simd_real_t ix_S0, iy_S0, iz_S0;
156     gmx_simd_real_t ix_S2, iy_S2, iz_S2;
157 } nbnxn_x_ci_simd_2xnn_t;
158
159 #endif
160
161 /* Working data for the actual i-supercell during pair search */
162 typedef struct nbnxn_list_work {
163     gmx_cache_protect_t     cp0;    /* Protect cache between threads               */
164
165     nbnxn_bb_t             *bb_ci;  /* The bounding boxes, pbc shifted, for each cluster */
166     float                  *pbb_ci; /* As bb_ci, but in xxxx packed format               */
167     real                   *x_ci;   /* The coordinates, pbc shifted, for each atom       */
168 #ifdef GMX_NBNXN_SIMD
169     nbnxn_x_ci_simd_4xn_t  *x_ci_simd_4xn;
170     nbnxn_x_ci_simd_2xnn_t *x_ci_simd_2xnn;
171 #endif
172     int                     cj_ind;          /* The current cj_ind index for the current list     */
173     int                     cj4_init;        /* The first unitialized cj4 block                   */
174
175     float                  *d2;              /* Bounding box distance work array                  */
176
177     nbnxn_cj_t             *cj;              /* The j-cell list                                   */
178     int                     cj_nalloc;       /* Allocation size of cj                             */
179
180     int                     ncj_noq;         /* Nr. of cluster pairs without Coul for flop count  */
181     int                     ncj_hlj;         /* Nr. of cluster pairs with 1/2 LJ for flop count   */
182
183     int                    *sort;            /* Sort index                    */
184     int                     sort_nalloc;     /* Allocation size of sort       */
185
186     nbnxn_sci_t            *sci_sort;        /* Second sci array, for sorting */
187     int                     sci_sort_nalloc; /* Allocation size of sci_sort   */
188
189     gmx_cache_protect_t     cp1;             /* Protect cache between threads               */
190 } nbnxn_list_work_t;
191
192 /* Function type for setting the i-atom coordinate working data */
193 typedef void
194     gmx_icell_set_x_t (int ci,
195                        real shx, real shy, real shz,
196                        int na_c,
197                        int stride, const real *x,
198                        nbnxn_list_work_t *work);
199
200 static gmx_icell_set_x_t icell_set_x_simple;
201 #ifdef GMX_NBNXN_SIMD
202 static gmx_icell_set_x_t icell_set_x_simple_simd_4xn;
203 static gmx_icell_set_x_t icell_set_x_simple_simd_2xnn;
204 #endif
205 static gmx_icell_set_x_t icell_set_x_supersub;
206 #ifdef NBNXN_SEARCH_SSE
207 static gmx_icell_set_x_t icell_set_x_supersub_sse8;
208 #endif
209
210 /* Local cycle count struct for profiling */
211 typedef struct {
212     int          count;
213     gmx_cycles_t c;
214     gmx_cycles_t start;
215 } nbnxn_cycle_t;
216
217 /* Local cycle count enum for profiling */
218 enum {
219     enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
220 };
221
222 /* Thread-local work struct, contains part of nbnxn_grid_t */
223 typedef struct {
224     gmx_cache_protect_t  cp0;
225
226     int                 *cxy_na;
227     int                  cxy_na_nalloc;
228
229     int                 *sort_work;
230     int                  sort_work_nalloc;
231
232     nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
233
234     int                  ndistc;       /* Number of distance checks for flop counting */
235
236     t_nblist            *nbl_fep;      /* Temporary FEP list for load balancing */
237
238     nbnxn_cycle_t        cc[enbsCCnr];
239
240     gmx_cache_protect_t  cp1;
241 } nbnxn_search_work_t;
242
243 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
244 typedef struct nbnxn_search {
245     gmx_bool            bFEP;            /* Do we have perturbed atoms? */
246     int                 ePBC;            /* PBC type enum                              */
247     matrix              box;             /* The periodic unit-cell                     */
248
249     gmx_bool            DomDec;          /* Are we doing domain decomposition?         */
250     ivec                dd_dim;          /* Are we doing DD in x,y,z?                  */
251     gmx_domdec_zones_t *zones;           /* The domain decomposition zones        */
252
253     int                 ngrid;           /* The number of grids, equal to #DD-zones    */
254     nbnxn_grid_t       *grid;            /* Array of grids, size ngrid                 */
255     int                *cell;            /* Actual allocated cell array for all grids  */
256     int                 cell_nalloc;     /* Allocation size of cell                    */
257     int                *a;               /* Atom index for grid, the inverse of cell   */
258     int                 a_nalloc;        /* Allocation size of a                       */
259
260     int                 natoms_local;    /* The local atoms run from 0 to natoms_local */
261     int                 natoms_nonlocal; /* The non-local atoms run from natoms_local
262                                           * to natoms_nonlocal */
263
264     gmx_bool             print_cycles;
265     int                  search_count;
266     nbnxn_cycle_t        cc[enbsCCnr];
267
268     gmx_icell_set_x_t   *icell_set_x; /* Function for setting i-coords    */
269
270     int                  nthread_max; /* Maximum number of threads for pair-search  */
271     nbnxn_search_work_t *work;        /* Work array, size nthread_max          */
272 } nbnxn_search_t_t;
273
274
275 static void nbs_cycle_start(nbnxn_cycle_t *cc)
276 {
277     cc->start = gmx_cycles_read();
278 }
279
280 static void nbs_cycle_stop(nbnxn_cycle_t *cc)
281 {
282     cc->c += gmx_cycles_read() - cc->start;
283     cc->count++;
284 }
285
286
287 #ifdef __cplusplus
288 }
289 #endif
290
291 #endif